47 #ifdef MFEM_USE_SIMMETRIX 54 #include <apfConvert.h> 59 #error This example requires that MFEM is built with MFEM_USE_PUMI=YES 65 int main(
int argc,
char *argv[])
68 Mpi::Init(argc, argv);
69 int num_procs = Mpi::WorldSize();
70 int myid = Mpi::WorldRank();
74 const char *mesh_file =
"../../data/pumi/serial/Kova.smb";
75 #ifdef MFEM_USE_SIMMETRIX 76 const char *model_file =
"../../data/pumi/geom/Kova.x_t";
78 const char *model_file =
"../../data/pumi/geom/Kova.dmg";
81 bool static_cond =
false;
82 bool visualization = 1;
86 args.
AddOption(&mesh_file,
"-m",
"--mesh",
89 "Finite element order (polynomial degree) or -1 for" 90 " isoparametric space.");
91 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
92 "--no-static-condensation",
"Enable static condensation.");
93 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
95 "Enable or disable GLVis visualization.");
96 args.
AddOption(&model_file,
"-p",
"--parasolid",
97 "Parasolid model to use.");
98 args.
AddOption(&geom_order,
"-go",
"--geometry_order",
99 "Geometric order of the model");
116 #ifdef MFEM_USE_SIMMETRIX 117 Sim_readLicenseFile(0);
123 apf::Mesh2* pumi_mesh;
124 pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
129 crv::BezierCurver bc(pumi_mesh, geom_order, 2);
147 (int)floor(log(50000./mesh->
GetNE())/log(2.)/
dim);
148 for (
int l = 0; l < ref_levels; l++)
165 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
172 cout <<
"Number of finite element unknowns: " 211 if (static_cond) {
a->EnableStaticCondensation(); }
216 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
218 cout <<
"Size of linear system: " << A.
Height() << endl;
220 #ifndef MFEM_USE_SUITESPARSE 224 PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
228 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
230 umf_solver.
Mult(B, X);
234 a->RecoverFEMSolution(X, *
b, x);
238 ofstream mesh_ofs(
"refined.mesh");
239 mesh_ofs.precision(8);
240 mesh->
Print(mesh_ofs);
241 ofstream sol_ofs(
"sol.gf");
242 sol_ofs.precision(8);
251 sol_sock.precision(8);
252 sol_sock <<
"solution\n" << *mesh << x << flush;
259 if (order > 0) {
delete fec; }
262 pumi_mesh->destroyNative();
263 apf::destroyMesh(pumi_mesh);
265 #ifdef MFEM_USE_SIMMETRIX 267 Sim_unregisterAllKeys();
Class for domain integration L(v) := (f, v)
Class for grid function - Vector with associated FE space.
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
int Dimension() const
Dimension of the reference space used within the elements.
void PrintUsage(std::ostream &out) const
Print the usage message.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Base class for PUMI meshes.
bool Good() const
Return true if the command line options were parsed successfully.
int main(int argc, char *argv[])
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
virtual const char * Name() const
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double Control[UMFPACK_CONTROL]
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int GetNE() const
Returns number of elements.
int Size() const
Return the logical size of the array.
virtual void Print(std::ostream &os=mfem::out) const
Arbitrary order H1-conforming (continuous) finite elements.
void GetNodes(Vector &node_coord) const
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.