// MFEM Example 1 // PUMI Modification // // Compile with: make ex1 // // Sample runs: // ex1 -m ../../data/pumi/serial/Kova.smb -p ../../data/pumi/geom/Kova.dmg // // Note: Example models + meshes for the PUMI examples can be downloaded // from github.com/mfem/data/pumi. After downloading we recommend // creating a symbolic link to the above directory in ../../data. // // Description: This example code demonstrates the use of MFEM to define a // simple finite element discretization of the Laplace problem // -Delta u = 1 with homogeneous Dirichlet boundary conditions. // Specifically, we discretize using a FE space of the specified // order, or if order < 1 using an isoparametric/isogeometric // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for // NURBS mesh, etc.) // // The example highlights the use of mesh refinement, finite // element grid functions, as well as linear and bilinear forms // corresponding to the left-hand side and right-hand side of the // discrete linear system. We also cover the explicit elimination // of essential boundary conditions, static condensation, and the // optional connection to the GLVis tool for visualization. // // This PUMI modification demonstrates how PUMI's API can be used // to load a PUMI mesh classified on a geometric model and then // convert it to the MFEM mesh format. The inputs are a Parasolid // model, "*.xmt_txt" and a SCOREC mesh "*.smb". The option "-o" // is used for the Finite Element order and "-go" is used for the // geometry order. Note that they can be used independently, i.e. // "-o 8 -go 3" solves for 8th order FE on a third order geometry. // // NOTE: Model/Mesh files for this example are in the (large) data file // repository of MFEM here https://github.com/mfem/data under the // folder named "pumi", which consists of the following sub-folders: // a) geom --> model files // b) parallel --> parallel pumi mesh files // c) serial --> serial pumi mesh files #include "mfem.hpp" #include #include #ifdef MFEM_USE_SIMMETRIX #include #include #endif #include #include #include #include #include #include #ifndef MFEM_USE_PUMI #error This example requires that MFEM is built with MFEM_USE_PUMI=YES #endif using namespace std; using namespace mfem; int main(int argc, char *argv[]) { // 1. Initialize MPI (required by PUMI) and HYPRE. Mpi::Init(argc, argv); int num_procs = Mpi::WorldSize(); int myid = Mpi::WorldRank(); Hypre::Init(); // 2. Parse command-line options. const char *mesh_file = "../../data/pumi/serial/Kova.smb"; #ifdef MFEM_USE_SIMMETRIX const char *model_file = "../../data/pumi/geom/Kova.x_t"; #else const char *model_file = "../../data/pumi/geom/Kova.dmg"; #endif int order = 1; bool static_cond = false; bool visualization = 1; int geom_order = 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) or -1 for" " isoparametric space."); 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.AddOption(&model_file, "-p", "--parasolid", "Parasolid model to use."); args.AddOption(&geom_order, "-go", "--geometry_order", "Geometric order of the model"); args.Parse(); if (!args.Good()) { if (myid == 0) { args.PrintUsage(cout); } return 1; } if (myid == 0) { args.PrintOptions(cout); } // 3. Read the SCOREC Mesh. PCU_Comm_Init(); #ifdef MFEM_USE_SIMMETRIX Sim_readLicenseFile(0); gmi_sim_start(); gmi_register_sim(); #endif gmi_register_mesh(); apf::Mesh2* pumi_mesh; pumi_mesh = apf::loadMdsMesh(model_file, mesh_file); // 4. Increase the geometry order if necessary. if (geom_order > 1) { crv::BezierCurver bc(pumi_mesh, geom_order, 2); bc.run(); } pumi_mesh->verify(); // 5. Create the MFEM mesh object from the PUMI mesh. We can handle // triangular and tetrahedral meshes. Other inputs are the same as the // MFEM default constructor. Mesh *mesh = new PumiMesh(pumi_mesh, 1, 1); int dim = mesh->Dimension(); // 6. 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 50,000 // elements. { int ref_levels = (int)floor(log(50000./mesh->GetNE())/log(2.)/dim); for (int l = 0; l < ref_levels; l++) { mesh->UniformRefinement(); } } // 7. Define a finite element space on the mesh. Here we use continuous // Lagrange finite elements of the specified order. If order < 1, we // instead use an isoparametric/isogeometric space. FiniteElementCollection *fec; if (order > 0) { fec = new H1_FECollection(order, dim); } else if (mesh->GetNodes()) { fec = mesh->GetNodes()->OwnFEC(); cout << "Using isoparametric FEs: " << fec->Name() << endl; } else { fec = new H1_FECollection(order = 1, dim); } FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec); cout << "Number of finite element unknowns: " << fespace->GetTrueVSize() << endl; // 8. Determine the list of true (i.e. conforming) essential boundary dofs. // In this example, the boundary conditions are defined by marking all // the boundary attributes from the mesh as essential (Dirichlet) and // converting them to a list of true dofs. Array ess_tdof_list; if (mesh->bdr_attributes.Size()) { Array ess_bdr(mesh->bdr_attributes.Max()); ess_bdr = 1; fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list); } // 9. Set up the linear form b(.) which corresponds to the right-hand side of // the FEM linear system, which in this case is (1,phi_i) where phi_i are // the basis functions in the finite element fespace. LinearForm *b = new LinearForm(fespace); ConstantCoefficient one(1.0); b->AddDomainIntegrator(new DomainLFIntegrator(one)); b->Assemble(); // 10. 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; // 11. Set up the bilinear form a(.,.) on the finite element space // corresponding to the Laplacian operator -Delta, by adding the // Diffusion domain integrator. BilinearForm *a = new BilinearForm(fespace); a->AddDomainIntegrator(new DiffusionIntegrator(one)); // 12. 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. if (static_cond) { a->EnableStaticCondensation(); } a->Assemble(); SparseMatrix A; Vector B, X; a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B); cout << "Size of linear system: " << A.Height() << endl; #ifndef MFEM_USE_SUITESPARSE // 13. Define a simple symmetric Gauss-Seidel preconditioner and use it to // solve the system A X = B with PCG. GSSmoother M(A); PCG(A, M, B, X, 1, 200, 1e-12, 0.0); #else // 13. 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 // 14. Recover the solution as a finite element grid function. a->RecoverFEMSolution(X, *b, x); // 15. Save the refined mesh and the solution. This output can be viewed later // using GLVis: "glvis -m refined.mesh -g sol.gf". ofstream mesh_ofs("refined.mesh"); mesh_ofs.precision(8); mesh->Print(mesh_ofs); ofstream sol_ofs("sol.gf"); sol_ofs.precision(8); x.Save(sol_ofs); // 16. Send the solution by socket to a GLVis server. if (visualization) { char vishost[] = "localhost"; int visport = 19916; socketstream sol_sock(vishost, visport); sol_sock.precision(8); sol_sock << "solution\n" << *mesh << x << flush; } // 17. Free the used memory. delete a; delete b; delete fespace; if (order > 0) { delete fec; } delete mesh; pumi_mesh->destroyNative(); apf::destroyMesh(pumi_mesh); PCU_Comm_Free(); #ifdef MFEM_USE_SIMMETRIX gmi_sim_stop(); Sim_unregisterAllKeys(); #endif return 0; }