61#ifdef MFEM_USE_SIMMETRIX
68#include <apfConvert.h>
73#error This example requires that MFEM is built with MFEM_USE_PUMI=YES
79int main(
int argc,
char *argv[])
86 const char *mesh_file =
"../../data/pumi/serial/pillbox.smb";
87 const char *boundary_file =
"../../data/pumi/serial/boundary.mesh";
88#ifdef MFEM_USE_SIMMETRIX
89 const char *model_file =
"../../data/pumi/geom/pillbox.smd";
91 const char *model_file =
"../../data/pumi/geom/pillbox.dmg";
94 bool static_cond =
false;
95 bool visualization = 1;
99 args.
AddOption(&mesh_file,
"-m",
"--mesh",
100 "Mesh file to use.");
102 "Finite element order (polynomial degree).");
103 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
104 "--no-static-condensation",
"Enable static condensation.");
105 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
106 "--no-visualization",
107 "Enable or disable GLVis visualization.");
108 args.
AddOption(&model_file,
"-p",
"--parasolid",
109 "Parasolid model to use.");
110 args.
AddOption(&geom_order,
"-go",
"--geometry_order",
111 "Geometric order of the model");
112 args.
AddOption(&boundary_file,
"-bf",
"--txt",
113 "txt file containing boundary tags");
124#ifdef MFEM_USE_SIMMETRIX
125 Sim_readLicenseFile(0);
131 apf::Mesh2* pumi_mesh;
132 pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
137 crv::BezierCurver bc(pumi_mesh, geom_order, 0);
146 getline(input_bdr, bdr_tags);
148 cout <<
" the boundary tag is : " << bdr_tags << endl;
151 if (bdr_tags ==
"Dirichlet")
153 input_bdr >> numOfent;
154 cout <<
" num of Dirichlet bdr conditions : " << numOfent << endl;
156 for (
int kk = 0; kk < numOfent; kk++)
158 input_bdr >> Dirichlet[kk];
165 input_bdr >> bdr_tags;
167 cout <<
" the boundary tag is : " << bdr_tags << endl;
168 if (bdr_tags ==
"Load")
170 input_bdr >> numOfent;
172 cout <<
" num of load bdr conditions : " << numOfent << endl;
173 for (
int kk = 0; kk < numOfent; kk++)
175 input_bdr >> load_bdr[kk];
187 apf::MeshIterator* itr = pumi_mesh->begin(
dim-1);
188 apf::MeshEntity* ent ;
190 while ((ent = pumi_mesh->iterate(itr)))
192 apf::ModelEntity *me = pumi_mesh->toModel(ent);
193 if (pumi_mesh->getModelType(me) == (
dim-1))
197 int tag = pumi_mesh->getModelTag(me);
198 if (Dirichlet.
Find(tag) != -1)
203 else if (load_bdr.
Find(tag) != -1)
216 for (
int el = 0; el < mesh->
GetNE(); el++)
220 if (cent(0) <= -0.05)
224 else if (cent(0) >= 0.05)
236 cerr <<
"\nInput mesh should have at least two materials and "
237 <<
"two boundary attributes! (See schematic in ex2.cpp)\n"
248 (int)floor(log(5000./mesh->
GetNE())/log(2.)/
dim);
249 for (
int l = 0; l < ref_levels; l++)
265 fespace = mesh->
GetNodes()->FESpace();
272 cout <<
"Number of finite element unknowns: " << fespace->
GetTrueVSize()
273 << endl <<
"Assembling: " << flush;
293 for (
int i = 0; i <
dim-1; i++)
300 pull_force(1) = -3.0e-2;
307 cout <<
"r.h.s. ... " << flush;
321 lambda(0) = lambda(1)*10;
322 lambda(1) = lambda(1)*100;
337 cout <<
"matrix ... " << flush;
338 if (static_cond) {
a->EnableStaticCondensation(); }
343 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
344 cout <<
"done." << endl;
346 cout <<
"Size of linear system: " << A.
Height() << endl;
348#ifndef MFEM_USE_SUITESPARSE
352 PCG(A, M, B, X, 1, 500, 1e-8, 0.0);
356 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
358 umf_solver.
Mult(B, X);
362 a->RecoverFEMSolution(X, *
b, x);
383 ofstream mesh_ofs(
"displaced.mesh");
384 mesh_ofs.precision(8);
385 mesh->
Print(mesh_ofs);
386 ofstream sol_ofs(
"sol.gf");
387 sol_ofs.precision(8);
398 sol_sock.precision(8);
399 sol_sock <<
"solution\n" << *mesh << x << flush;
412 pumi_mesh->destroyNative();
413 apf::destroyMesh(pumi_mesh);
415#ifdef MFEM_USE_SIMMETRIX
417 Sim_unregisterAllKeys();
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Find(const T &el) const
Return the first index where 'el' is found; return -1 if not found.
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
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.
const IntegrationPoint & GetCenter(int GeomType) const
Return the center of the given Geometry::Type, GeomType.
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.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void SetAttribute(int i, int attr)
Set the attribute of element i.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Print the mesh to the given stream using the default MFEM mesh format.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
virtual void SetAttributes(bool elem_attrs_changed=true, bool bdr_face_attrs_changed=true)
Determine the sets of unique attribute values in domain if elem_attrs_changed and boundary elements i...
void GetNodes(Vector &node_coord) const
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Geometry::Type GetElementBaseGeometry(int i) const
Array< int > attributes
A list of all unique element attributes used by the Mesh.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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.
Base class for PUMI meshes.
Direct sparse solver using UMFPACK.
real_t Control[UMFPACK_CONTROL]
void SetOperator(const Operator &op) override
Factorize the given Operator op which must be a SparseMatrix.
void Mult(const Vector &b, Vector &x) const override
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)
void filter_dos(std::string &line)
Check for, and remove, a trailing '\r' from and std::string.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void skip_comment_lines(std::istream &is, const char comment_char)
Check if the stream starts with comment_char. If so skip it.
std::array< int, NCMesh::MaxFaceNodes > nodes