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