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