1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261
|
// MFEM Example 2
//
// Compile with: make ex2
//
// Sample runs: ex2 -m ../data/beam-tri.mesh
// ex2 -m ../data/beam-quad.mesh
// ex2 -m ../data/beam-tet.mesh
// ex2 -m ../data/beam-hex.mesh
// ex2 -m ../data/beam-wedge.mesh
// ex2 -m ../data/beam-quad.mesh -o 3 -sc
// ex2 -m ../data/beam-quad-nurbs.mesh
// ex2 -m ../data/beam-hex-nurbs.mesh
//
// Description: This example code solves a simple linear elasticity problem
// describing a multi-material cantilever beam.
//
// Specifically, we approximate the weak form of -div(sigma(u))=0
// where sigma(u)=lambda*div(u)*I+mu*(grad*u+u*grad) is the stress
// tensor corresponding to displacement field u, and lambda and mu
// are the material Lame constants. The boundary conditions are
// u=0 on the fixed part of the boundary with attribute 1, and
// sigma(u).n=f on the remainder with f being a constant pull down
// vector on boundary elements with attribute 2, and zero
// otherwise. The geometry of the domain is assumed to be as
// follows:
//
// +----------+----------+
// boundary --->| material | material |<--- boundary
// attribute 1 | 1 | 2 | attribute 2
// (fixed) +----------+----------+ (pull down)
//
// The example demonstrates the use of high-order and NURBS vector
// finite element spaces with the linear elasticity bilinear form,
// meshes with curved elements, and the definition of piece-wise
// constant and vector coefficient objects. Static condensation is
// also illustrated.
//
// We recommend viewing Example 1 before viewing this example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
const char *mesh_file = "../data/beam-tri.mesh";
int order = 1;
bool static_cond = false;
bool visualization = 1;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree).");
args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
"--no-static-condensation", "Enable static condensation.");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
args.PrintOptions(cout);
// 2. Read the mesh from the given mesh file. We can handle triangular,
// quadrilateral, tetrahedral or hexahedral elements with the same code.
Mesh *mesh = new Mesh(mesh_file, 1, 1);
int dim = mesh->Dimension();
if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
{
cerr << "\nInput mesh should have at least two materials and "
<< "two boundary attributes! (See schematic in ex2.cpp)\n"
<< endl;
return 3;
}
// 3. Select the order of the finite element discretization space. For NURBS
// meshes, we increase the order by degree elevation.
if (mesh->NURBSext)
{
mesh->DegreeElevate(order, order);
}
// 4. Refine the mesh to increase the resolution. In this example we do
// 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
// largest number that gives a final mesh with no more than 5,000
// elements.
{
int ref_levels =
(int)floor(log(5000./mesh->GetNE())/log(2.)/dim);
for (int l = 0; l < ref_levels; l++)
{
mesh->UniformRefinement();
}
}
// 5. Define a finite element space on the mesh. Here we use vector finite
// elements, i.e. dim copies of a scalar finite element space. The vector
// dimension is specified by the last argument of the FiniteElementSpace
// constructor. For NURBS meshes, we use the (degree elevated) NURBS space
// associated with the mesh nodes.
FiniteElementCollection *fec;
FiniteElementSpace *fespace;
if (mesh->NURBSext)
{
fec = NULL;
fespace = mesh->GetNodes()->FESpace();
}
else
{
fec = new H1_FECollection(order, dim);
fespace = new FiniteElementSpace(mesh, fec, dim);
}
cout << "Number of finite element unknowns: " << fespace->GetTrueVSize()
<< endl << "Assembling: " << flush;
// 6. Determine the list of true (i.e. conforming) essential boundary dofs.
// In this example, the boundary conditions are defined by marking only
// boundary attribute 1 from the mesh as essential and converting it to a
// list of true dofs.
Array<int> ess_tdof_list, ess_bdr(mesh->bdr_attributes.Max());
ess_bdr = 0;
ess_bdr[0] = 1;
fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
// 7. Set up the linear form b(.) which corresponds to the right-hand side of
// the FEM linear system. In this case, b_i equals the boundary integral
// of f*phi_i where f represents a "pull down" force on the Neumann part
// of the boundary and phi_i are the basis functions in the finite element
// fespace. The force is defined by the VectorArrayCoefficient object f,
// which is a vector of Coefficient objects. The fact that f is non-zero
// on boundary attribute 2 is indicated by the use of piece-wise constants
// coefficient for its last component.
VectorArrayCoefficient f(dim);
for (int i = 0; i < dim-1; i++)
{
f.Set(i, new ConstantCoefficient(0.0));
}
{
Vector pull_force(mesh->bdr_attributes.Max());
pull_force = 0.0;
pull_force(1) = -1.0e-2;
f.Set(dim-1, new PWConstCoefficient(pull_force));
}
LinearForm *b = new LinearForm(fespace);
b->AddBoundaryIntegrator(new VectorBoundaryLFIntegrator(f));
cout << "r.h.s. ... " << flush;
b->Assemble();
// 8. Define the solution vector x as a finite element grid function
// corresponding to fespace. Initialize x with initial guess of zero,
// which satisfies the boundary conditions.
GridFunction x(fespace);
x = 0.0;
// 9. Set up the bilinear form a(.,.) on the finite element space
// corresponding to the linear elasticity integrator with piece-wise
// constants coefficient lambda and mu.
Vector lambda(mesh->attributes.Max());
lambda = 1.0;
lambda(0) = lambda(1)*50;
PWConstCoefficient lambda_func(lambda);
Vector mu(mesh->attributes.Max());
mu = 1.0;
mu(0) = mu(1)*50;
PWConstCoefficient mu_func(mu);
BilinearForm *a = new BilinearForm(fespace);
a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func,mu_func));
// 10. Assemble the bilinear form and the corresponding linear system,
// applying any necessary transformations such as: eliminating boundary
// conditions, applying conforming constraints for non-conforming AMR,
// static condensation, etc.
cout << "matrix ... " << flush;
if (static_cond) { a->EnableStaticCondensation(); }
a->Assemble();
SparseMatrix A;
Vector B, X;
a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
cout << "done." << endl;
cout << "Size of linear system: " << A.Height() << endl;
#ifndef MFEM_USE_SUITESPARSE
// 11. Define a simple symmetric Gauss-Seidel preconditioner and use it to
// solve the system Ax=b with PCG.
GSSmoother M(A);
PCG(A, M, B, X, 1, 500, 1e-8, 0.0);
#else
// 11. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
UMFPackSolver umf_solver;
umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
umf_solver.SetOperator(A);
umf_solver.Mult(B, X);
#endif
// 12. Recover the solution as a finite element grid function.
a->RecoverFEMSolution(X, *b, x);
// 13. For non-NURBS meshes, make the mesh curved based on the finite element
// space. This means that we define the mesh elements through a fespace
// based transformation of the reference element. This allows us to save
// the displaced mesh as a curved mesh when using high-order finite
// element displacement field. We assume that the initial mesh (read from
// the file) is not higher order curved mesh compared to the chosen FE
// space.
if (!mesh->NURBSext)
{
mesh->SetNodalFESpace(fespace);
}
// 14. Save the displaced mesh and the inverted solution (which gives the
// backward displacements to the original grid). This output can be
// viewed later using GLVis: "glvis -m displaced.mesh -g sol.gf".
{
GridFunction *nodes = mesh->GetNodes();
*nodes += x;
x *= -1;
ofstream mesh_ofs("displaced.mesh");
mesh_ofs.precision(8);
mesh->Print(mesh_ofs);
ofstream sol_ofs("sol.gf");
sol_ofs.precision(8);
x.Save(sol_ofs);
}
// 15. Send the above data by socket to a GLVis server. Use the "n" and "b"
// keys in GLVis to visualize the displacements.
if (visualization)
{
char vishost[] = "localhost";
int visport = 19916;
socketstream sol_sock(vishost, visport);
sol_sock.precision(8);
sol_sock << "solution\n" << *mesh << x << flush;
}
// 16. Free the used memory.
delete a;
delete b;
if (fec)
{
delete fespace;
delete fec;
}
delete mesh;
return 0;
}
|