// MFEM Example 1 - Parallel Version // PUMI Modification // // Compile with: make ex1p // // Sample runs: // mpirun -np 8 ex1p -m ../../data/pumi/parallel/Kova/Kova100k_8.smb // -p ../../data/pumi/geom/Kova.dmg -o 1 -go 2 // // 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 parallel PUMI mesh classified on a geometric model // and then generate the corresponding parallel MFEM mesh. The // example also performs a "uniform" refinement, similar to the // MFEM examples, for coarse meshes. However, the refinement is // performed using the PUMI API. The inputs are a Parasolid // model, "*.xmt_txt" and SCOREC parallel meshes "*.smb". The // option "-o" is used for the Finite Element order and "-go" for // the geometry order. Note that they can be used independently: // "-o 8 -go 3" solves for 8th order FE on 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 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/parallel/Kova/Kova100k_8.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 and refine the mesh if necessary. Parallel // uniform refinement is performed if the total number of elements is less // than 10,000. int dim = pumi_mesh->getDimension(); int nEle = pumi_mesh->count(dim); int ref_levels = (int)floor(log(10000./nEle)/log(2.)/dim); if (geom_order > 1) { crv::BezierCurver bc(pumi_mesh, geom_order, 2); bc.run(); } // Perform Uniform refinement if (ref_levels > 1) { auto uniInput = ma::configureUniformRefine(pumi_mesh, ref_levels); if (geom_order > 1) { crv::adapt(uniInput); } else { ma::adapt(uniInput); } } pumi_mesh->verify(); // 5. Create the parallel MFEM mesh object from the parallel PUMI mesh. // We can handle triangular and tetrahedral meshes. Note that the // mesh resolution is performed on the PUMI mesh. ParMesh *pmesh = new ParPumiMesh(MPI_COMM_WORLD, pumi_mesh); // 6. Define a parallel finite element space on the parallel 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 (pmesh->GetNodes()) { fec = pmesh->GetNodes()->OwnFEC(); if (myid == 0) { cout << "Using isoparametric FEs: " << fec->Name() << endl; } } else { fec = new H1_FECollection(order = 1, dim); } ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec); HYPRE_BigInt size = fespace->GlobalTrueVSize(); if (myid == 0) { cout << "Number of finite element unknowns: " << size << endl; } // 7. Determine the list of true (i.e. parallel 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 (pmesh->bdr_attributes.Size()) { Array ess_bdr(pmesh->bdr_attributes.Max()); ess_bdr = 1; fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list); } // 8. Set up the parallel 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 fespace. ParLinearForm *b = new ParLinearForm(fespace); ConstantCoefficient one(1.0); b->AddDomainIntegrator(new DomainLFIntegrator(one)); b->Assemble(); // 9. Define the solution vector x as a parallel finite element grid function // corresponding to fespace. Initialize x with initial guess of zero, // which satisfies the boundary conditions. ParGridFunction x(fespace); x = 0.0; // 10. Set up the parallel bilinear form a(.,.) on the finite element space // corresponding to the Laplacian operator -Delta, by adding the Diffusion // domain integrator. ParBilinearForm *a = new ParBilinearForm(fespace); a->AddDomainIntegrator(new DiffusionIntegrator(one)); // 11. Assemble the parallel bilinear form and the corresponding linear // system, applying any necessary transformations such as: parallel // assembly, eliminating boundary conditions, applying conforming // constraints for non-conforming AMR, static condensation, etc. if (static_cond) { a->EnableStaticCondensation(); } a->Assemble(); HypreParMatrix A; Vector B, X; a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B); if (myid == 0) { cout << "Size of linear system: " << A.GetGlobalNumRows() << endl; } // 12. Define and apply a parallel PCG solver for AX=B with the BoomerAMG // preconditioner from hypre. HypreSolver *amg = new HypreBoomerAMG(A); HyprePCG *pcg = new HyprePCG(A); pcg->SetTol(1e-12); pcg->SetMaxIter(200); pcg->SetPrintLevel(2); pcg->SetPreconditioner(*amg); pcg->Mult(B, X); // 13. Recover the parallel grid function corresponding to X. This is the // local finite element solution on each processor. a->RecoverFEMSolution(X, *b, x); // 14. Save the refined mesh and the solution in parallel. This output can // be viewed later using GLVis: "glvis -np -m mesh -g sol". { ostringstream mesh_name, sol_name; mesh_name << "mesh." << setfill('0') << setw(6) << myid; sol_name << "sol." << setfill('0') << setw(6) << myid; ofstream mesh_ofs(mesh_name.str().c_str()); mesh_ofs.precision(8); pmesh->Print(mesh_ofs); ofstream sol_ofs(sol_name.str().c_str()); sol_ofs.precision(8); x.Save(sol_ofs); } // 15. Send the solution by socket to a GLVis server. if (visualization) { char vishost[] = "localhost"; int visport = 19916; socketstream sol_sock(vishost, visport); sol_sock << "parallel " << num_procs << " " << myid << "\n"; sol_sock.precision(8); sol_sock << "solution\n" << *pmesh << x << flush; } // 16. Free the used memory. delete pcg; delete amg; delete a; delete b; delete fespace; if (order > 0) { delete fec; } delete pmesh; pumi_mesh->destroyNative(); apf::destroyMesh(pumi_mesh); PCU_Comm_Free(); #ifdef MFEM_USE_SIMMETRIX gmi_sim_stop(); Sim_unregisterAllKeys(); #endif return 0; }