47int main(
int argc,
char *argv[])
50 const char *mesh_file =
"../data/beam-tri.mesh";
52 bool static_cond =
false;
53 bool visualization = 1;
56 args.
AddOption(&mesh_file,
"-m",
"--mesh",
59 "Finite element order (polynomial degree).");
60 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
61 "--no-static-condensation",
"Enable static condensation.");
62 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
64 "Enable or disable GLVis visualization.");
75 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
80 cerr <<
"\nInput mesh should have at least two materials and "
81 <<
"two boundary attributes! (See schematic in ex2.cpp)\n"
99 (int)floor(log(5000./mesh->
GetNE())/log(2.)/
dim);
100 for (
int l = 0; l < ref_levels; l++)
116 fespace = mesh->
GetNodes()->FESpace();
123 cout <<
"Number of finite element unknowns: " << fespace->
GetTrueVSize()
124 << endl <<
"Assembling: " << flush;
144 for (
int i = 0; i <
dim-1; i++)
151 pull_force(1) = -1.0e-2;
157 cout <<
"r.h.s. ... " << flush;
171 lambda(0) = lambda(1)*50;
185 cout <<
"matrix ... " << flush;
186 if (static_cond) {
a->EnableStaticCondensation(); }
191 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
192 cout <<
"done." << endl;
194 cout <<
"Size of linear system: " << A.
Height() << endl;
196#ifndef MFEM_USE_SUITESPARSE
200 PCG(A, M, B, X, 1, 500, 1e-8, 0.0);
204 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
206 umf_solver.
Mult(B, X);
210 a->RecoverFEMSolution(X, *
b, x);
231 ofstream mesh_ofs(
"displaced.mesh");
232 mesh_ofs.precision(8);
233 mesh->
Print(mesh_ofs);
234 ofstream sol_ofs(
"sol.gf");
235 sol_ofs.precision(8);
246 sol_sock.precision(8);
247 sol_sock <<
"solution\n" << *mesh << x << flush;
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
A coefficient that is constant across space and time.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Arbitrary order H1-conforming (continuous) finite elements.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
void DegreeElevate(int rel_degree, int degree=16)
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Array< int > attributes
A list of all unique element attributes used by the Mesh.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
bool Good() const
Return true if the command line options were parsed successfully.
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
Direct sparse solver using UMFPACK.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
real_t Control[UMFPACK_CONTROL]
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using UMFPACK.
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
std::function< real_t(const Vector &)> f(real_t mass_coeff)