1182fbe45STzanio // libCEED + MFEM Example: BP1 2182fbe45STzanio // 3182fbe45STzanio // This example illustrates a simple usage of libCEED with the MFEM (mfem.org) 4182fbe45STzanio // finite element library. 5182fbe45STzanio // 6182fbe45STzanio // The example reads a mesh from a file and solves a simple linear system with a 7182fbe45STzanio // mass matrix (L2-projection of a given analytic function provided by 8182fbe45STzanio // 'solution'). The mass matrix required for performing the projection is 9182fbe45STzanio // expressed as a new class, CeedMassOperator, derived from mfem::Operator. 10182fbe45STzanio // Internally, CeedMassOperator uses a CeedOperator object constructed based on 11182fbe45STzanio // an mfem::FiniteElementSpace. All libCEED objects use a Ceed device object 12182fbe45STzanio // constructed based on a command line argument (-ceed). 13182fbe45STzanio // 14182fbe45STzanio // The mass matrix is inverted using a simple conjugate gradient algorithm 15182fbe45STzanio // corresponding to CEED BP1, see http://ceed.exascaleproject.org/bps. Arbitrary 16182fbe45STzanio // mesh and solution orders in 1D, 2D and 3D are supported from the same code. 17182fbe45STzanio // 18182fbe45STzanio // Build with: 19182fbe45STzanio // 20182fbe45STzanio // make bp1 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>] 21182fbe45STzanio // 22182fbe45STzanio // Sample runs: 23182fbe45STzanio // 2466087c08SValeria Barra // ./bp1 2566087c08SValeria Barra // ./bp1 -ceed /cpu/self 2628688798Sjeremylt // ./bp1 -ceed /gpu/cuda 2766087c08SValeria Barra // ./bp1 -m ../../../mfem/data/fichera.mesh 2866087c08SValeria Barra // ./bp1 -m ../../../mfem/data/star.vtk -o 3 2966087c08SValeria Barra // ./bp1 -m ../../../mfem/data/inline-segment.mesh -o 8 30182fbe45STzanio 315d6bafb2Sjeremylt /// @file 325d6bafb2Sjeremylt /// MFEM mass operator based on libCEED 335d6bafb2Sjeremylt 34c0c38e35SVeselin Dobrev #include "bp1.hpp" 35182fbe45STzanio 36*2b730f8bSJeremy L Thompson #include <ceed.h> 37*2b730f8bSJeremy L Thompson 38*2b730f8bSJeremy L Thompson #include <mfem.hpp> 39*2b730f8bSJeremy L Thompson 40182fbe45STzanio /// Continuous function to project on the discrete FE space 41182fbe45STzanio double solution(const mfem::Vector &pt) { 42182fbe45STzanio return pt.Norml2(); // distance to the origin 43182fbe45STzanio } 44182fbe45STzanio 4543dae957SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000 --order 4 46182fbe45STzanio int main(int argc, char *argv[]) { 47182fbe45STzanio // 1. Parse command-line options. 48182fbe45STzanio const char *ceed_spec = "/cpu/self"; 49c0c38e35SVeselin Dobrev #ifndef MFEM_DIR 50182fbe45STzanio const char *mesh_file = "../../../mfem/data/star.mesh"; 51c0c38e35SVeselin Dobrev #else 52c0c38e35SVeselin Dobrev const char *mesh_file = MFEM_DIR "/data/star.mesh"; 53c0c38e35SVeselin Dobrev #endif 54182fbe45STzanio int order = 1; 55182fbe45STzanio bool visualization = true; 56dc00e230Sjeremylt bool test = false; 57e2b2c771Svaleria double max_nnodes = 50000; 58182fbe45STzanio 59182fbe45STzanio mfem::OptionsParser args(argc, argv); 60182fbe45STzanio args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 61182fbe45STzanio args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 62*2b730f8bSJeremy L Thompson args.AddOption(&order, "-o", "--order", "Finite element order (polynomial degree)."); 63e2b2c771Svaleria args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)"); 64*2b730f8bSJeremy L Thompson args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", "--no-visualization", "Enable or disable GLVis visualization."); 65*2b730f8bSJeremy L Thompson args.AddOption(&test, "-t", "--test", "-no-test", "--no-test", "Enable or disable test mode."); 66182fbe45STzanio args.Parse(); 67182fbe45STzanio if (!args.Good()) { 68182fbe45STzanio args.PrintUsage(std::cout); 69182fbe45STzanio return 1; 70182fbe45STzanio } 71dc00e230Sjeremylt if (!test) { 72182fbe45STzanio args.PrintOptions(std::cout); 73dc00e230Sjeremylt } 74182fbe45STzanio 75182fbe45STzanio // 2. Initialize a Ceed device object using the given Ceed specification. 76182fbe45STzanio Ceed ceed; 77182fbe45STzanio CeedInit(ceed_spec, &ceed); 78182fbe45STzanio 79182fbe45STzanio // 3. Read the mesh from the given mesh file. 80182fbe45STzanio mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1); 81182fbe45STzanio int dim = mesh->Dimension(); 82182fbe45STzanio 83182fbe45STzanio // 4. Refine the mesh to increase the resolution. In this example we do 84182fbe45STzanio // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the 85182fbe45STzanio // largest number that gives a final system with no more than 50,000 86182fbe45STzanio // unknowns, approximately. 87182fbe45STzanio { 88*2b730f8bSJeremy L Thompson int ref_levels = (int)floor((log(max_nnodes / mesh->GetNE()) - dim * log(order)) / log(2.) / dim); 89182fbe45STzanio for (int l = 0; l < ref_levels; l++) { 90182fbe45STzanio mesh->UniformRefinement(); 91182fbe45STzanio } 92182fbe45STzanio } 93182fbe45STzanio if (mesh->GetNodalFESpace() == NULL) { 94182fbe45STzanio mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 95182fbe45STzanio } 96182fbe45STzanio if (mesh->NURBSext) { 97182fbe45STzanio mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 98182fbe45STzanio } 99182fbe45STzanio 100182fbe45STzanio // 5. Define a finite element space on the mesh. Here we use continuous 101182fbe45STzanio // Lagrange finite elements of the specified order. 102182fbe45STzanio MFEM_VERIFY(order > 0, "invalid order"); 103182fbe45STzanio mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 104182fbe45STzanio mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 105dc00e230Sjeremylt if (!test) { 106*2b730f8bSJeremy L Thompson std::cout << "Number of finite element unknowns: " << fespace->GetTrueVSize() << std::endl; 107dc00e230Sjeremylt } 108182fbe45STzanio 109182fbe45STzanio // 6. Construct a rhs vector using the linear form f(v) = (solution, v), where 110182fbe45STzanio // v is a test function. 111182fbe45STzanio mfem::LinearForm b(fespace); 112182fbe45STzanio mfem::FunctionCoefficient sol_coeff(solution); 113182fbe45STzanio b.AddDomainIntegrator(new mfem::DomainLFIntegrator(sol_coeff)); 114182fbe45STzanio b.Assemble(); 115182fbe45STzanio 116182fbe45STzanio // 7. Construct a CeedMassOperator utilizing the 'ceed' device and using the 117182fbe45STzanio // 'fespace' object to extract data needed by the Ceed objects. 118182fbe45STzanio CeedMassOperator mass(ceed, fespace); 119182fbe45STzanio 120182fbe45STzanio // 8. Solve the discrete system using the conjugate gradients (CG) method. 121182fbe45STzanio mfem::CGSolver cg; 122182fbe45STzanio cg.SetRelTol(1e-6); 123182fbe45STzanio cg.SetMaxIter(100); 124dc00e230Sjeremylt if (test) { 125dc00e230Sjeremylt cg.SetPrintLevel(0); 126dc00e230Sjeremylt } else { 127182fbe45STzanio cg.SetPrintLevel(3); 128dc00e230Sjeremylt } 129182fbe45STzanio cg.SetOperator(mass); 130182fbe45STzanio 131182fbe45STzanio mfem::GridFunction sol(fespace); 132182fbe45STzanio sol = 0.0; 133182fbe45STzanio cg.Mult(b, sol); 134182fbe45STzanio 135182fbe45STzanio // 9. Compute and print the L2 projection error. 1369647a07eSDavid Medina double err_l2 = sol.ComputeL2Error(sol_coeff); 137dc00e230Sjeremylt if (!test) { 138*2b730f8bSJeremy L Thompson std::cout << "L2 projection error: " << err_l2 << std::endl; 139dc00e230Sjeremylt } else { 140f063656dSJed Brown if (fabs(sol.ComputeL2Error(sol_coeff)) > 2e-4) { 1419647a07eSDavid Medina std::cout << "Error too large: " << err_l2 << std::endl; 142dc00e230Sjeremylt } 1439b872752Sjeremylt } 144182fbe45STzanio 145182fbe45STzanio // 10. Open a socket connection to GLVis and send the mesh and solution for 146182fbe45STzanio // visualization. 147182fbe45STzanio if (visualization) { 148182fbe45STzanio char vishost[] = "localhost"; 149182fbe45STzanio int visport = 19916; 150182fbe45STzanio mfem::socketstream sol_sock(vishost, visport); 151182fbe45STzanio sol_sock.precision(8); 152182fbe45STzanio sol_sock << "solution\n" << *mesh << sol << std::flush; 153182fbe45STzanio } 154182fbe45STzanio 155182fbe45STzanio // 11. Free memory and exit. 156182fbe45STzanio delete fespace; 157182fbe45STzanio delete fec; 158182fbe45STzanio delete mesh; 159182fbe45STzanio CeedDestroy(&ceed); 160182fbe45STzanio return 0; 161182fbe45STzanio } 162