1 // libCEED + MFEM Example: BP1 2 // 3 // This example illustrates a simple usage of libCEED with the MFEM (mfem.org) 4 // finite element library. 5 // 6 // The example reads a mesh from a file and solves a simple linear system with a 7 // mass matrix (L2-projection of a given analytic function provided by 8 // 'solution'). The mass matrix required for performing the projection is 9 // expressed as a new class, CeedMassOperator, derived from mfem::Operator. 10 // Internally, CeedMassOperator uses a CeedOperator object constructed based on 11 // an mfem::FiniteElementSpace. All libCEED objects use a Ceed device object 12 // constructed based on a command line argument (-ceed). 13 // 14 // The mass matrix is inverted using a simple conjugate gradient algorithm 15 // corresponding to CEED BP1, see http://ceed.exascaleproject.org/bps. Arbitrary 16 // mesh and solution orders in 1D, 2D and 3D are supported from the same code. 17 // 18 // Build with: 19 // 20 // make bp1 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>] 21 // 22 // Sample runs: 23 // 24 // ./bp1 25 // ./bp1 -ceed /cpu/self 26 // ./bp1 -ceed /gpu/cuda 27 // ./bp1 -m ../../../mfem/data/fichera.mesh 28 // ./bp1 -m ../../../mfem/data/star.vtk -o 3 29 // ./bp1 -m ../../../mfem/data/inline-segment.mesh -o 8 30 31 /// @file 32 /// MFEM mass operator based on libCEED 33 34 #include "bp1.hpp" 35 36 #include <ceed.h> 37 38 #include <mfem.hpp> 39 40 /// Continuous function to project on the discrete FE space 41 double solution(const mfem::Vector &pt) { 42 return pt.Norml2(); // distance to the origin 43 } 44 45 //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000 --order 4 46 int main(int argc, char *argv[]) { 47 // 1. Parse command-line options. 48 const char *ceed_spec = "/cpu/self"; 49 #ifndef MFEM_DIR 50 const char *mesh_file = "../../../mfem/data/star.mesh"; 51 #else 52 const char *mesh_file = MFEM_DIR "/data/star.mesh"; 53 #endif 54 int order = 1; 55 bool visualization = true; 56 bool test = false; 57 double max_nnodes = 50000; 58 59 mfem::OptionsParser args(argc, argv); 60 args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 61 args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 62 args.AddOption(&order, "-o", "--order", "Finite element order (polynomial degree)."); 63 args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)"); 64 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", "--no-visualization", "Enable or disable GLVis visualization."); 65 args.AddOption(&test, "-t", "--test", "-no-test", "--no-test", "Enable or disable test mode."); 66 args.Parse(); 67 if (!args.Good()) { 68 args.PrintUsage(std::cout); 69 return 1; 70 } 71 if (!test) { 72 args.PrintOptions(std::cout); 73 } 74 75 // 2. Initialize a Ceed device object using the given Ceed specification. 76 Ceed ceed; 77 CeedInit(ceed_spec, &ceed); 78 79 // 3. Read the mesh from the given mesh file. 80 mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1); 81 int dim = mesh->Dimension(); 82 83 // 4. Refine the mesh to increase the resolution. In this example we do 84 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the 85 // largest number that gives a final system with no more than 50,000 86 // unknowns, approximately. 87 { 88 int ref_levels = (int)floor((log(max_nnodes / mesh->GetNE()) - dim * log(order)) / log(2.) / dim); 89 for (int l = 0; l < ref_levels; l++) { 90 mesh->UniformRefinement(); 91 } 92 } 93 if (mesh->GetNodalFESpace() == NULL) { 94 mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 95 } 96 if (mesh->NURBSext) { 97 mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 98 } 99 100 // 5. Define a finite element space on the mesh. Here we use continuous 101 // Lagrange finite elements of the specified order. 102 MFEM_VERIFY(order > 0, "invalid order"); 103 mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 104 mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 105 if (!test) { 106 std::cout << "Number of finite element unknowns: " << fespace->GetTrueVSize() << std::endl; 107 } 108 109 // 6. Construct a rhs vector using the linear form f(v) = (solution, v), where 110 // v is a test function. 111 mfem::LinearForm b(fespace); 112 mfem::FunctionCoefficient sol_coeff(solution); 113 b.AddDomainIntegrator(new mfem::DomainLFIntegrator(sol_coeff)); 114 b.Assemble(); 115 116 // 7. Construct a CeedMassOperator utilizing the 'ceed' device and using the 117 // 'fespace' object to extract data needed by the Ceed objects. 118 CeedMassOperator mass(ceed, fespace); 119 120 // 8. Solve the discrete system using the conjugate gradients (CG) method. 121 mfem::CGSolver cg; 122 cg.SetRelTol(1e-6); 123 cg.SetMaxIter(100); 124 if (test) { 125 cg.SetPrintLevel(0); 126 } else { 127 cg.SetPrintLevel(3); 128 } 129 cg.SetOperator(mass); 130 131 mfem::GridFunction sol(fespace); 132 sol = 0.0; 133 cg.Mult(b, sol); 134 135 // 9. Compute and print the L2 projection error. 136 double err_l2 = sol.ComputeL2Error(sol_coeff); 137 if (!test) { 138 std::cout << "L2 projection error: " << err_l2 << std::endl; 139 } else { 140 if (fabs(sol.ComputeL2Error(sol_coeff)) > 2e-4) { 141 std::cout << "Error too large: " << err_l2 << std::endl; 142 } 143 } 144 145 // 10. Open a socket connection to GLVis and send the mesh and solution for 146 // visualization. 147 if (visualization) { 148 char vishost[] = "localhost"; 149 int visport = 19916; 150 mfem::socketstream sol_sock(vishost, visport); 151 sol_sock.precision(8); 152 sol_sock << "solution\n" << *mesh << sol << std::flush; 153 } 154 155 // 11. Free memory and exit. 156 delete fespace; 157 delete fec; 158 delete mesh; 159 CeedDestroy(&ceed); 160 return 0; 161 } 162