1182fbe45STzanio // libCEED + MFEM Example: BP3 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 linear system with a 7182fbe45STzanio // diffusion stiffness matrix (with a prescribed analytic solution, provided by 8182fbe45STzanio // the function 'solution'). The diffusion matrix is expressed as a new class, 9182fbe45STzanio // CeedDiffusionOperator, derived from mfem::Operator. Internally, 10182fbe45STzanio // CeedDiffusionOperator uses a CeedOperator object constructed based on an 11182fbe45STzanio // mfem::FiniteElementSpace. All libCEED objects use a Ceed logical device 12182fbe45STzanio // object constructed based on a command line argument. (-ceed). 13182fbe45STzanio // 14182fbe45STzanio // The linear system is inverted using the conjugate gradients algorithm 15182fbe45STzanio // corresponding to CEED BP3, 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 bp3 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>] 21182fbe45STzanio // 22182fbe45STzanio // Sample runs: 23182fbe45STzanio // 24*66087c08SValeria Barra // ./bp3 25*66087c08SValeria Barra // ./bp3 -ceed /cpu/self 26*66087c08SValeria Barra // ./bp3 -m ../../../mfem/data/fichera.mesh -o 4 27*66087c08SValeria Barra // ./bp3 -m ../../../mfem/data/square-disc-nurbs.mesh -o 6 28*66087c08SValeria Barra // ./bp3 -m ../../../mfem/data/inline-segment.mesh -o 8 29182fbe45STzanio 305d6bafb2Sjeremylt /// @file 315d6bafb2Sjeremylt /// MFEM diffusion operator based on libCEED 325d6bafb2Sjeremylt 33182fbe45STzanio #include <ceed.h> 34182fbe45STzanio #include <mfem.hpp> 35c0c38e35SVeselin Dobrev #include "bp3.hpp" 36182fbe45STzanio 37182fbe45STzanio /// Exact solution 38182fbe45STzanio double solution(const mfem::Vector &pt) { 39182fbe45STzanio static const double x[3] = { -0.32, 0.15, 0.24 }; 40182fbe45STzanio static const double k[3] = { 1.21, 1.45, 1.37 }; 41182fbe45STzanio double val = sin(M_PI*(x[0]+k[0]*pt(0))); 42182fbe45STzanio for (int d = 1; d < pt.Size(); d++) 43182fbe45STzanio val *= sin(M_PI*(x[d]+k[d]*pt(d))); 44182fbe45STzanio return val; 45182fbe45STzanio } 46182fbe45STzanio 47182fbe45STzanio /// Right-hand side 48182fbe45STzanio double rhs(const mfem::Vector &pt) { 49182fbe45STzanio static const double x[3] = { -0.32, 0.15, 0.24 }; 50182fbe45STzanio static const double k[3] = { 1.21, 1.45, 1.37 }; 51182fbe45STzanio double f[3], l[3], val, lap; 52182fbe45STzanio f[0] = sin(M_PI*(x[0]+k[0]*pt(0))); 53182fbe45STzanio l[0] = M_PI*M_PI*k[0]*k[0]*f[0]; 54182fbe45STzanio val = f[0]; 55182fbe45STzanio lap = l[0]; 56182fbe45STzanio for (int d = 1; d < pt.Size(); d++) { 57182fbe45STzanio f[d] = sin(M_PI*(x[d]+k[d]*pt(d))); 58182fbe45STzanio l[d] = M_PI*M_PI*k[d]*k[d]*f[d]; 59182fbe45STzanio lap = lap*f[d] + val*l[d]; 60182fbe45STzanio val = val*f[d]; 61182fbe45STzanio } 62182fbe45STzanio return lap; 63182fbe45STzanio } 64182fbe45STzanio 65f063656dSJed Brown //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000 66182fbe45STzanio int main(int argc, char *argv[]) { 67182fbe45STzanio // 1. Parse command-line options. 68182fbe45STzanio const char *ceed_spec = "/cpu/self"; 69c0c38e35SVeselin Dobrev #ifndef MFEM_DIR 70182fbe45STzanio const char *mesh_file = "../../../mfem/data/star.mesh"; 71c0c38e35SVeselin Dobrev #else 72c0c38e35SVeselin Dobrev const char *mesh_file = MFEM_DIR "/data/star.mesh"; 73c0c38e35SVeselin Dobrev #endif 74182fbe45STzanio int order = 2; 75182fbe45STzanio bool visualization = true; 76dc00e230Sjeremylt bool test = false; 77e2b2c771Svaleria double max_nnodes = 50000; 78182fbe45STzanio 79182fbe45STzanio mfem::OptionsParser args(argc, argv); 80182fbe45STzanio args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 81182fbe45STzanio args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 82182fbe45STzanio args.AddOption(&order, "-o", "--order", 83182fbe45STzanio "Finite element order (polynomial degree)."); 84e2b2c771Svaleria args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)"); 85182fbe45STzanio args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", 86182fbe45STzanio "--no-visualization", 87182fbe45STzanio "Enable or disable GLVis visualization."); 88dc00e230Sjeremylt args.AddOption(&test, "-t", "--test", "-no-test", 89dc00e230Sjeremylt "--no-test", 90dc00e230Sjeremylt "Enable or disable test mode."); 91182fbe45STzanio args.Parse(); 92182fbe45STzanio if (!args.Good()) { 93182fbe45STzanio args.PrintUsage(std::cout); 94182fbe45STzanio return 1; 95182fbe45STzanio } 96dc00e230Sjeremylt if (!test) { 97182fbe45STzanio args.PrintOptions(std::cout); 98dc00e230Sjeremylt } 99182fbe45STzanio 100182fbe45STzanio // 2. Initialize a Ceed device object using the given Ceed specification. 101182fbe45STzanio Ceed ceed; 102182fbe45STzanio CeedInit(ceed_spec, &ceed); 103182fbe45STzanio 104182fbe45STzanio // 3. Read the mesh from the given mesh file. 105182fbe45STzanio mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1); 106182fbe45STzanio int dim = mesh->Dimension(); 107182fbe45STzanio 108182fbe45STzanio // 4. Refine the mesh to increase the resolution. In this example we do 109182fbe45STzanio // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the 11031d4d2baSJed Brown // largest number that gives a final system with no more than 50,000 11131d4d2baSJed Brown // unknowns, approximately. 112182fbe45STzanio { 113182fbe45STzanio int ref_levels = 114e2b2c771Svaleria (int)floor((log(max_nnodes/mesh->GetNE())-dim*log(order))/log(2.)/dim); 115182fbe45STzanio for (int l = 0; l < ref_levels; l++) { 116182fbe45STzanio mesh->UniformRefinement(); 117182fbe45STzanio } 118182fbe45STzanio } 119182fbe45STzanio if (mesh->GetNodalFESpace() == NULL) { 120182fbe45STzanio mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 121182fbe45STzanio } 122182fbe45STzanio if (mesh->NURBSext) { 123182fbe45STzanio mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 124182fbe45STzanio } 125182fbe45STzanio 126182fbe45STzanio // 5. Define a finite element space on the mesh. Here we use continuous 127182fbe45STzanio // Lagrange finite elements of the specified order. 128182fbe45STzanio MFEM_VERIFY(order > 0, "invalid order"); 129182fbe45STzanio mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 130182fbe45STzanio mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 131dc00e230Sjeremylt if (!test) { 132182fbe45STzanio std::cout << "Number of finite element unknowns: " 133182fbe45STzanio << fespace->GetTrueVSize() << std::endl; 134dc00e230Sjeremylt } 135182fbe45STzanio 136182fbe45STzanio mfem::FunctionCoefficient sol_coeff(solution); 137182fbe45STzanio mfem::Array<int> ess_tdof_list; 138182fbe45STzanio mfem::GridFunction sol(fespace); 139182fbe45STzanio if (mesh->bdr_attributes.Size()) { 140182fbe45STzanio mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max()); 141182fbe45STzanio ess_bdr = 1; 142182fbe45STzanio fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list); 143182fbe45STzanio sol.ProjectBdrCoefficient(sol_coeff, ess_bdr); 144182fbe45STzanio } 145182fbe45STzanio 146182fbe45STzanio // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where 147182fbe45STzanio // v is a test function. 148182fbe45STzanio mfem::LinearForm b(fespace); 149182fbe45STzanio mfem::FunctionCoefficient rhs_coeff(rhs); 150182fbe45STzanio b.AddDomainIntegrator(new mfem::DomainLFIntegrator(rhs_coeff)); 151182fbe45STzanio b.Assemble(); 152182fbe45STzanio 153182fbe45STzanio // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using 154182fbe45STzanio // the 'fespace' object to extract data needed by the Ceed objects. 155182fbe45STzanio CeedDiffusionOperator diff(ceed, fespace); 156182fbe45STzanio 157182fbe45STzanio mfem::Operator *D; 158182fbe45STzanio mfem::Vector X, B; 159182fbe45STzanio diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B); 160182fbe45STzanio 161182fbe45STzanio // 8. Solve the discrete system using the conjugate gradients (CG) method. 162182fbe45STzanio mfem::CGSolver cg; 163182fbe45STzanio cg.SetRelTol(1e-6); 164182fbe45STzanio cg.SetMaxIter(1000); 165dc00e230Sjeremylt if (test) { 166dc00e230Sjeremylt cg.SetPrintLevel(0); 167dc00e230Sjeremylt } else { 168182fbe45STzanio cg.SetPrintLevel(3); 169dc00e230Sjeremylt } 170182fbe45STzanio cg.SetOperator(*D); 171182fbe45STzanio 172182fbe45STzanio cg.Mult(B, X); 173182fbe45STzanio 174182fbe45STzanio // 9. Compute and print the L2 norm of the error. 175dc00e230Sjeremylt if (!test) { 176dc00e230Sjeremylt std::cout << "L2 projection error: " << sol.ComputeL2Error(sol_coeff) 177182fbe45STzanio << std::endl; 178dc00e230Sjeremylt } else { 179f063656dSJed Brown if (fabs(sol.ComputeL2Error(sol_coeff))>2e-3) { 180dc00e230Sjeremylt std::cout << "Error too large" << std::endl; 181dc00e230Sjeremylt } 182dc00e230Sjeremylt } 183182fbe45STzanio 184182fbe45STzanio // 10. Open a socket connection to GLVis and send the mesh and solution for 185182fbe45STzanio // visualization. 186182fbe45STzanio if (visualization) { 187182fbe45STzanio char vishost[] = "localhost"; 188182fbe45STzanio int visport = 19916; 189182fbe45STzanio mfem::socketstream sol_sock(vishost, visport); 190182fbe45STzanio sol_sock.precision(8); 191182fbe45STzanio sol_sock << "solution\n" << *mesh << sol << std::flush; 192182fbe45STzanio } 193182fbe45STzanio 194182fbe45STzanio // 11. Free memory and exit. 195182fbe45STzanio delete fespace; 196182fbe45STzanio delete fec; 197182fbe45STzanio delete mesh; 198182fbe45STzanio CeedDestroy(&ceed); 199182fbe45STzanio return 0; 200182fbe45STzanio } 201