1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors. 2ea61e9acSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ea61e9acSJeremy L Thompson // 4ea61e9acSJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5ea61e9acSJeremy L Thompson // 6ea61e9acSJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7ea61e9acSJeremy L Thompson 8182fbe45STzanio // libCEED + MFEM Example: BP1 9182fbe45STzanio // 10ea61e9acSJeremy L Thompson // This example illustrates a simple usage of libCEED with the MFEM (mfem.org) finite element library. 11182fbe45STzanio // 12ea61e9acSJeremy L Thompson // The example reads a mesh from a file and solves a simple linear system with a mass matrix (L2-projection of a given analytic function provided by 13ea61e9acSJeremy L Thompson // 'solution'). The mass matrix required for performing the projection is expressed as a new class, CeedMassOperator, derived from mfem::Operator. 14ea61e9acSJeremy L Thompson // Internally, CeedMassOperator uses a CeedOperator object constructed based on an mfem::FiniteElementSpace. 15ea61e9acSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command line argument (-ceed). 16182fbe45STzanio // 17ea61e9acSJeremy L Thompson // The mass matrix is inverted using a simple conjugate gradient algorithm corresponding to CEED BP1, see http://ceed.exascaleproject.org/bps. 18ea61e9acSJeremy L Thompson // Arbitrary mesh and solution orders in 1D, 2D and 3D are supported from the same code. 19182fbe45STzanio // 20182fbe45STzanio // Build with: 21182fbe45STzanio // 22182fbe45STzanio // make bp1 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>] 23182fbe45STzanio // 24182fbe45STzanio // Sample runs: 25182fbe45STzanio // 2666087c08SValeria Barra // ./bp1 2766087c08SValeria Barra // ./bp1 -ceed /cpu/self 2828688798Sjeremylt // ./bp1 -ceed /gpu/cuda 2966087c08SValeria Barra // ./bp1 -m ../../../mfem/data/fichera.mesh 3066087c08SValeria Barra // ./bp1 -m ../../../mfem/data/star.vtk -o 3 3166087c08SValeria Barra // ./bp1 -m ../../../mfem/data/inline-segment.mesh -o 8 32182fbe45STzanio 335d6bafb2Sjeremylt /// @file 345d6bafb2Sjeremylt /// MFEM mass operator based on libCEED 355d6bafb2Sjeremylt 36c0c38e35SVeselin Dobrev #include "bp1.hpp" 37182fbe45STzanio 382b730f8bSJeremy L Thompson #include <ceed.h> 392b730f8bSJeremy L Thompson 402b730f8bSJeremy L Thompson #include <mfem.hpp> 412b730f8bSJeremy L Thompson 42182fbe45STzanio /// Continuous function to project on the discrete FE space 43182fbe45STzanio double solution(const mfem::Vector &pt) { 44182fbe45STzanio return pt.Norml2(); // distance to the origin 45182fbe45STzanio } 46182fbe45STzanio 4743dae957SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000 --order 4 48182fbe45STzanio int main(int argc, char *argv[]) { 49182fbe45STzanio // 1. Parse command-line options. 50182fbe45STzanio const char *ceed_spec = "/cpu/self"; 51c0c38e35SVeselin Dobrev #ifndef MFEM_DIR 52182fbe45STzanio const char *mesh_file = "../../../mfem/data/star.mesh"; 53c0c38e35SVeselin Dobrev #else 54c0c38e35SVeselin Dobrev const char *mesh_file = MFEM_DIR "/data/star.mesh"; 55c0c38e35SVeselin Dobrev #endif 56182fbe45STzanio int order = 1; 57182fbe45STzanio bool visualization = true; 58dc00e230Sjeremylt bool test = false; 59e2b2c771Svaleria double max_nnodes = 50000; 60182fbe45STzanio 61182fbe45STzanio mfem::OptionsParser args(argc, argv); 62182fbe45STzanio args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 63182fbe45STzanio args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 642b730f8bSJeremy L Thompson args.AddOption(&order, "-o", "--order", "Finite element order (polynomial degree)."); 65e2b2c771Svaleria args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)"); 662b730f8bSJeremy L Thompson args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", "--no-visualization", "Enable or disable GLVis visualization."); 672b730f8bSJeremy L Thompson args.AddOption(&test, "-t", "--test", "-no-test", "--no-test", "Enable or disable test mode."); 68182fbe45STzanio args.Parse(); 69182fbe45STzanio if (!args.Good()) { 70182fbe45STzanio args.PrintUsage(std::cout); 71182fbe45STzanio return 1; 72182fbe45STzanio } 73dc00e230Sjeremylt if (!test) { 74182fbe45STzanio args.PrintOptions(std::cout); 75dc00e230Sjeremylt } 76182fbe45STzanio 77182fbe45STzanio // 2. Initialize a Ceed device object using the given Ceed specification. 78182fbe45STzanio Ceed ceed; 79182fbe45STzanio CeedInit(ceed_spec, &ceed); 80182fbe45STzanio 81182fbe45STzanio // 3. Read the mesh from the given mesh file. 82182fbe45STzanio mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1); 83182fbe45STzanio int dim = mesh->Dimension(); 84182fbe45STzanio 85ea61e9acSJeremy L Thompson // 4. Refine the mesh to increase the resolution. 86ea61e9acSJeremy L Thompson // In this example we do 'ref_levels' of uniform refinement. 87ea61e9acSJeremy L Thompson // We choose 'ref_levels' to be the largest number that gives a final system with no more than 50,000 unknowns, approximately. 88182fbe45STzanio { 892b730f8bSJeremy L Thompson int ref_levels = (int)floor((log(max_nnodes / mesh->GetNE()) - dim * log(order)) / log(2.) / dim); 90182fbe45STzanio for (int l = 0; l < ref_levels; l++) { 91182fbe45STzanio mesh->UniformRefinement(); 92182fbe45STzanio } 93182fbe45STzanio } 94182fbe45STzanio if (mesh->GetNodalFESpace() == NULL) { 95182fbe45STzanio mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 96182fbe45STzanio } 97182fbe45STzanio if (mesh->NURBSext) { 98182fbe45STzanio mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 99182fbe45STzanio } 100182fbe45STzanio 101ea61e9acSJeremy L Thompson // 5. Define a finite element space on the mesh. 102ea61e9acSJeremy L Thompson // Here we use continuous Lagrange finite elements of the specified order. 103182fbe45STzanio MFEM_VERIFY(order > 0, "invalid order"); 104182fbe45STzanio mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 105182fbe45STzanio mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 106dc00e230Sjeremylt if (!test) { 1072b730f8bSJeremy L Thompson std::cout << "Number of finite element unknowns: " << fespace->GetTrueVSize() << std::endl; 108dc00e230Sjeremylt } 109182fbe45STzanio 110ea61e9acSJeremy L Thompson // 6. Construct a rhs vector using the linear form f(v) = (solution, v), where 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 116ea61e9acSJeremy L Thompson // 7. Construct a CeedMassOperator utilizing the 'ceed' device and using the 'fespace' object to extract data needed by the Ceed objects. 117182fbe45STzanio CeedMassOperator mass(ceed, fespace); 118182fbe45STzanio 119182fbe45STzanio // 8. Solve the discrete system using the conjugate gradients (CG) method. 120182fbe45STzanio mfem::CGSolver cg; 121182fbe45STzanio cg.SetRelTol(1e-6); 122182fbe45STzanio cg.SetMaxIter(100); 123dc00e230Sjeremylt if (test) { 124dc00e230Sjeremylt cg.SetPrintLevel(0); 125dc00e230Sjeremylt } else { 126182fbe45STzanio cg.SetPrintLevel(3); 127dc00e230Sjeremylt } 128182fbe45STzanio cg.SetOperator(mass); 129182fbe45STzanio 130182fbe45STzanio mfem::GridFunction sol(fespace); 131182fbe45STzanio sol = 0.0; 132182fbe45STzanio cg.Mult(b, sol); 133182fbe45STzanio 134182fbe45STzanio // 9. Compute and print the L2 projection error. 1359647a07eSDavid Medina double err_l2 = sol.ComputeL2Error(sol_coeff); 136dc00e230Sjeremylt if (!test) { 1372b730f8bSJeremy L Thompson std::cout << "L2 projection error: " << err_l2 << std::endl; 138dc00e230Sjeremylt } else { 139f063656dSJed Brown if (fabs(sol.ComputeL2Error(sol_coeff)) > 2e-4) { 1409647a07eSDavid Medina std::cout << "Error too large: " << err_l2 << std::endl; 141dc00e230Sjeremylt } 1429b872752Sjeremylt } 143182fbe45STzanio 144ea61e9acSJeremy L Thompson // 10. Open a socket connection to GLVis and send the mesh and solution for visualization. 145182fbe45STzanio if (visualization) { 146182fbe45STzanio char vishost[] = "localhost"; 147182fbe45STzanio int visport = 19916; 148182fbe45STzanio mfem::socketstream sol_sock(vishost, visport); 149182fbe45STzanio sol_sock.precision(8); 150182fbe45STzanio sol_sock << "solution\n" << *mesh << sol << std::flush; 151182fbe45STzanio } 152182fbe45STzanio 153182fbe45STzanio // 11. Free memory and exit. 154182fbe45STzanio delete fespace; 155182fbe45STzanio delete fec; 156182fbe45STzanio delete mesh; 157182fbe45STzanio CeedDestroy(&ceed); 158182fbe45STzanio return 0; 159182fbe45STzanio } 160