// MFEM Example 1 // AmgX Modification // // Compile with: make ex1 // // AmgX sample runs: // ex1 // ex1 -d cuda // ex1 --amgx-file multi_gs.json --amgx-solver // ex1 --amgx-file precon.json --amgx-preconditioner // ex1 --amgx-file multi_gs.json --amgx-solver -d cuda // ex1 --amgx-file precon.json --amgx-preconditioner -d cuda // // 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. #include "mfem.hpp" #include #include #ifndef MFEM_USE_AMGX #error This example requires that MFEM is built with MFEM_USE_AMGX=YES #endif using namespace std; using namespace mfem; int main(int argc, char *argv[]) { // 1. Parse command-line options. const char *mesh_file = "../../data/star.mesh"; int order = 1; bool static_cond = false; bool pa = false; const char *device_config = "cpu"; bool visualization = true; bool amgx_lib = true; bool amgx_solver = true; const char* amgx_json_file = ""; // JSON file for AmgX 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(&pa, "-pa", "--partial-assembly", "-no-pa", "--no-partial-assembly", "Enable Partial Assembly."); args.AddOption(&amgx_lib, "-amgx", "--amgx-lib", "-no-amgx", "--no-amgx-lib", "Use AmgX in example."); args.AddOption(&amgx_json_file, "--amgx-file", "--amgx-file", "AMGX solver config file (overrides --amgx-solver, --amgx-verbose)"); args.AddOption(&amgx_solver, "--amgx-solver", "--amgx-solver", "--amgx-preconditioner", "--amgx-preconditioner", "Configure AMGX as solver or preconditioner."); args.AddOption(&device_config, "-d", "--device", "Device configuration string, see Device::Configure()."); args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", "--no-visualization", "Enable or disable GLVis visualization."); args.Parse(); if (!args.Good()) { args.PrintUsage(cout); return 1; } args.PrintOptions(cout); // 2. Enable hardware devices such as GPUs, and programming models such as // CUDA, OCCA, RAJA and OpenMP based on command line options. Device device(device_config); device.Print(); // 3. Read the mesh from the given mesh file. We can handle triangular, // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with // the same code. Mesh mesh(mesh_file, 1, 1); int dim = mesh.Dimension(); // 4. 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(); } } // 5. 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; bool delete_fec; if (order > 0) { fec = new H1_FECollection(order, dim); delete_fec = true; } else if (mesh.GetNodes()) { fec = mesh.GetNodes()->OwnFEC(); delete_fec = false; cout << "Using isoparametric FEs: " << fec->Name() << endl; } else { fec = new H1_FECollection(order = 1, dim); delete_fec = true; } FiniteElementSpace fespace(&mesh, fec); cout << "Number of finite element unknowns: " << fespace.GetTrueVSize() << endl; // 6. 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); } // 7. 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(&fespace); ConstantCoefficient one(1.0); b.AddDomainIntegrator(new DomainLFIntegrator(one)); b.Assemble(); // 8. 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; // 9. 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(&fespace); if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); } a.AddDomainIntegrator(new DiffusionIntegrator(one)); // 10. 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(); OperatorPtr A; Vector B, X; a.FormLinearSystem(ess_tdof_list, x, b, A, X, B); cout << "Size of linear system: " << A->Height() << endl; // 11. Solve the linear system A X = B. if (pa) { // Jacobi preconditioning in partial assembly mode if (UsesTensorBasis(fespace)) { OperatorJacobiSmoother M(a, ess_tdof_list); PCG(*A, M, B, X, 1, 400, 1e-12, 0.0); } else { CG(*A, B, X, 1, 400, 1e-12, 0.0); } } else if (amgx_lib && strcmp(amgx_json_file,"") == 0) { bool amgx_verbose = false; AmgXSolver amgx(AmgXSolver::PRECONDITIONER, amgx_verbose); amgx.SetOperator(*A.As()); PCG(*A, amgx, B, X, 1, 200, 1e-12, 0.0); } else if (amgx_lib && strcmp(amgx_json_file,"") != 0) { AmgXSolver amgx; amgx.ReadParameters(amgx_json_file, AmgXSolver::EXTERNAL); amgx.InitSerial(); amgx.SetOperator(*A.As()); if (amgx_solver) { amgx.SetConvergenceCheck(true); amgx.Mult(B,X); } else { // Omit convergence check at the AmgX level when using as a // preconditioner. amgx.SetConvergenceCheck(false); PCG(*A.As(), amgx, B, X, 3, 40, 1e-12, 0.0); } } else { #ifndef MFEM_USE_SUITESPARSE // Use a simple symmetric Gauss-Seidel preconditioner with PCG. GSSmoother M((SparseMatrix&)(*A)); PCG(*A, M, B, X, 1, 200, 1e-12, 0.0); #else // 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 } // 12. Recover the solution as a finite element grid function. a.RecoverFEMSolution(X, b, x); // 13. 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); // 14. 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; } // 15. Free the used memory. if (delete_fec) { delete fec; } return 0; }