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[])
88 const char *mesh_file =
"../../data/pumi/serial/pillbox.smb";
89 const char *boundary_file =
"../../data/pumi/serial/boundary.mesh";
90#ifdef MFEM_USE_SIMMETRIX
91 const char *model_file =
"../../data/pumi/geom/pillbox.smd";
93 const char *model_file =
"../../data/pumi/geom/pillbox.dmg";
96 bool static_cond =
false;
97 bool visualization = 1;
101 args.
AddOption(&mesh_file,
"-m",
"--mesh",
102 "Mesh file to use.");
104 "Finite element order (polynomial degree).");
105 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
106 "--no-static-condensation",
"Enable static condensation.");
107 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
108 "--no-visualization",
109 "Enable or disable GLVis visualization.");
110 args.
AddOption(&model_file,
"-p",
"--parasolid",
111 "Parasolid model to use.");
112 args.
AddOption(&geom_order,
"-go",
"--geometry_order",
113 "Geometric order of the model");
114 args.
AddOption(&boundary_file,
"-bf",
"--txt",
115 "txt file containing boundary tags");
126#ifdef MFEM_USE_SIMMETRIX
127 Sim_readLicenseFile(0);
133 apf::Mesh2* pumi_mesh;
134 pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
139 crv::BezierCurver bc(pumi_mesh, geom_order, 0);
148 getline(input_bdr, bdr_tags);
150 cout <<
" the boundary tag is : " << bdr_tags << endl;
153 if (bdr_tags ==
"Dirichlet")
155 input_bdr >> numOfent;
156 cout <<
" num of Dirichlet bdr conditions : " << numOfent << endl;
158 for (
int kk = 0; kk < numOfent; kk++)
160 input_bdr >> Dirichlet[kk];
167 input_bdr >> bdr_tags;
169 cout <<
" the boundary tag is : " << bdr_tags << endl;
170 if (bdr_tags ==
"Load")
172 input_bdr >> numOfent;
174 cout <<
" num of load bdr conditions : " << numOfent << endl;
175 for (
int kk = 0; kk < numOfent; kk++)
177 input_bdr >> load_bdr[kk];
189 apf::MeshIterator* itr = pumi_mesh->begin(
dim-1);
190 apf::MeshEntity* ent ;
192 while ((ent = pumi_mesh->iterate(itr)))
194 apf::ModelEntity *me = pumi_mesh->toModel(ent);
195 if (pumi_mesh->getModelType(me) == (
dim-1))
199 int tag = pumi_mesh->getModelTag(me);
200 if (Dirichlet.
Find(tag) != -1)
205 else if (load_bdr.
Find(tag) != -1)
218 for (
int el = 0; el < mesh->
GetNE(); el++)
222 if (cent(0) <= -0.05)
226 else if (cent(0) >= 0.05)
238 cerr <<
"\nInput mesh should have at least two materials and "
239 <<
"two boundary attributes! (See schematic in ex2.cpp)\n"
250 (int)floor(log(5000./mesh->
GetNE())/log(2.)/
dim);
251 for (
int l = 0; l < ref_levels; l++)
267 fespace = mesh->
GetNodes()->FESpace();
274 cout <<
"Number of finite element unknowns: " << fespace->
GetTrueVSize()
275 << endl <<
"Assembling: " << flush;
295 for (
int i = 0; i <
dim-1; i++)
302 pull_force(1) = -3.0e-2;
309 cout <<
"r.h.s. ... " << flush;
323 lambda(0) = lambda(1)*10;
324 lambda(1) = lambda(1)*100;
339 cout <<
"matrix ... " << flush;
340 if (static_cond) {
a->EnableStaticCondensation(); }
345 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
346 cout <<
"done." << endl;
348 cout <<
"Size of linear system: " << A.
Height() << endl;
350#ifndef MFEM_USE_SUITESPARSE
354 PCG(A, M, B, X, 1, 500, 1e-8, 0.0);
358 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
360 umf_solver.
Mult(B, X);
364 a->RecoverFEMSolution(X, *
b, x);
385 ofstream mesh_ofs(
"displaced.mesh");
386 mesh_ofs.precision(8);
387 mesh->
Print(mesh_ofs);
388 ofstream sol_ofs(
"sol.gf");
389 sol_ofs.precision(8);
400 sol_sock.precision(8);
401 sol_sock <<
"solution\n" << *mesh << x << flush;
414 pumi_mesh->destroyNative();
415 apf::destroyMesh(pumi_mesh);
417#ifdef MFEM_USE_SIMMETRIX
419 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
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...
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.
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
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.
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)
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.