1*182fbe45STzanio // libCEED + MFEM Example: BP3 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 linear system with a 7*182fbe45STzanio // diffusion stiffness matrix (with a prescribed analytic solution, provided by 8*182fbe45STzanio // the function 'solution'). The diffusion matrix is expressed as a new class, 9*182fbe45STzanio // CeedDiffusionOperator, derived from mfem::Operator. Internally, 10*182fbe45STzanio // CeedDiffusionOperator uses a CeedOperator object constructed based on an 11*182fbe45STzanio // mfem::FiniteElementSpace. All libCEED objects use a Ceed logical device 12*182fbe45STzanio // object constructed based on a command line argument. (-ceed). 13*182fbe45STzanio // 14*182fbe45STzanio // The linear system is inverted using the conjugate gradients algorithm 15*182fbe45STzanio // corresponding to CEED BP3, 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 bp3 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>] 21*182fbe45STzanio // 22*182fbe45STzanio // Sample runs: 23*182fbe45STzanio // 24*182fbe45STzanio // bp3 25*182fbe45STzanio // bp3 -ceed /cpu/self 26*182fbe45STzanio // bp3 -m ../../../mfem/data/fichera.mesh -o 4 27*182fbe45STzanio // bp3 -m ../../../mfem/data/square-disc-nurbs.mesh -o 6 28*182fbe45STzanio // bp3 -m ../../../mfem/data/inline-segment.mesh -o 8 29*182fbe45STzanio 30*182fbe45STzanio #include <ceed.h> 31*182fbe45STzanio #include <mfem.hpp> 32*182fbe45STzanio #include <bp3.hpp> 33*182fbe45STzanio 34*182fbe45STzanio /// Exact solution 35*182fbe45STzanio double solution(const mfem::Vector &pt) { 36*182fbe45STzanio static const double x[3] = { -0.32, 0.15, 0.24 }; 37*182fbe45STzanio static const double k[3] = { 1.21, 1.45, 1.37 }; 38*182fbe45STzanio double val = sin(M_PI*(x[0]+k[0]*pt(0))); 39*182fbe45STzanio for (int d = 1; d < pt.Size(); d++) 40*182fbe45STzanio val *= sin(M_PI*(x[d]+k[d]*pt(d))); 41*182fbe45STzanio return val; 42*182fbe45STzanio } 43*182fbe45STzanio 44*182fbe45STzanio /// Right-hand side 45*182fbe45STzanio double rhs(const mfem::Vector &pt) { 46*182fbe45STzanio static const double x[3] = { -0.32, 0.15, 0.24 }; 47*182fbe45STzanio static const double k[3] = { 1.21, 1.45, 1.37 }; 48*182fbe45STzanio double f[3], l[3], val, lap; 49*182fbe45STzanio f[0] = sin(M_PI*(x[0]+k[0]*pt(0))); 50*182fbe45STzanio l[0] = M_PI*M_PI*k[0]*k[0]*f[0]; 51*182fbe45STzanio val = f[0]; 52*182fbe45STzanio lap = l[0]; 53*182fbe45STzanio for (int d = 1; d < pt.Size(); d++) { 54*182fbe45STzanio f[d] = sin(M_PI*(x[d]+k[d]*pt(d))); 55*182fbe45STzanio l[d] = M_PI*M_PI*k[d]*k[d]*f[d]; 56*182fbe45STzanio lap = lap*f[d] + val*l[d]; 57*182fbe45STzanio val = val*f[d]; 58*182fbe45STzanio } 59*182fbe45STzanio return lap; 60*182fbe45STzanio } 61*182fbe45STzanio 62*182fbe45STzanio 63*182fbe45STzanio int main(int argc, char *argv[]) { 64*182fbe45STzanio // 1. Parse command-line options. 65*182fbe45STzanio const char *ceed_spec = "/cpu/self"; 66*182fbe45STzanio const char *mesh_file = "../../../mfem/data/star.mesh"; 67*182fbe45STzanio int order = 2; 68*182fbe45STzanio bool visualization = true; 69*182fbe45STzanio 70*182fbe45STzanio mfem::OptionsParser args(argc, argv); 71*182fbe45STzanio args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 72*182fbe45STzanio args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 73*182fbe45STzanio args.AddOption(&order, "-o", "--order", 74*182fbe45STzanio "Finite element order (polynomial degree)."); 75*182fbe45STzanio args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", 76*182fbe45STzanio "--no-visualization", 77*182fbe45STzanio "Enable or disable GLVis visualization."); 78*182fbe45STzanio args.Parse(); 79*182fbe45STzanio if (!args.Good()) { 80*182fbe45STzanio args.PrintUsage(std::cout); 81*182fbe45STzanio return 1; 82*182fbe45STzanio } 83*182fbe45STzanio args.PrintOptions(std::cout); 84*182fbe45STzanio 85*182fbe45STzanio // 2. Initialize a Ceed device object using the given Ceed specification. 86*182fbe45STzanio Ceed ceed; 87*182fbe45STzanio CeedInit(ceed_spec, &ceed); 88*182fbe45STzanio 89*182fbe45STzanio // 3. Read the mesh from the given mesh file. 90*182fbe45STzanio mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1); 91*182fbe45STzanio int dim = mesh->Dimension(); 92*182fbe45STzanio 93*182fbe45STzanio // 4. Refine the mesh to increase the resolution. In this example we do 94*182fbe45STzanio // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the 95*182fbe45STzanio // largest number that gives a final system with no more than 50,000 (1,000 96*182fbe45STzanio // in 1D) unknowns, approximately. 97*182fbe45STzanio { 98*182fbe45STzanio double max_dofs = (dim > 1) ? 50000 : 1000; 99*182fbe45STzanio int ref_levels = 100*182fbe45STzanio (int)floor((log(max_dofs/mesh->GetNE())-dim*log(order))/log(2.)/dim); 101*182fbe45STzanio for (int l = 0; l < ref_levels; l++) { 102*182fbe45STzanio mesh->UniformRefinement(); 103*182fbe45STzanio } 104*182fbe45STzanio } 105*182fbe45STzanio if (mesh->GetNodalFESpace() == NULL) { 106*182fbe45STzanio mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 107*182fbe45STzanio } 108*182fbe45STzanio if (mesh->NURBSext) { 109*182fbe45STzanio mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 110*182fbe45STzanio } 111*182fbe45STzanio 112*182fbe45STzanio // 5. Define a finite element space on the mesh. Here we use continuous 113*182fbe45STzanio // Lagrange finite elements of the specified order. 114*182fbe45STzanio MFEM_VERIFY(order > 0, "invalid order"); 115*182fbe45STzanio mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 116*182fbe45STzanio mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 117*182fbe45STzanio std::cout << "Number of finite element unknowns: " 118*182fbe45STzanio << fespace->GetTrueVSize() << std::endl; 119*182fbe45STzanio 120*182fbe45STzanio mfem::FunctionCoefficient sol_coeff(solution); 121*182fbe45STzanio mfem::Array<int> ess_tdof_list; 122*182fbe45STzanio mfem::GridFunction sol(fespace); 123*182fbe45STzanio if (mesh->bdr_attributes.Size()) { 124*182fbe45STzanio mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max()); 125*182fbe45STzanio ess_bdr = 1; 126*182fbe45STzanio fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list); 127*182fbe45STzanio sol.ProjectBdrCoefficient(sol_coeff, ess_bdr); 128*182fbe45STzanio } 129*182fbe45STzanio 130*182fbe45STzanio // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where 131*182fbe45STzanio // v is a test function. 132*182fbe45STzanio mfem::LinearForm b(fespace); 133*182fbe45STzanio mfem::FunctionCoefficient rhs_coeff(rhs); 134*182fbe45STzanio b.AddDomainIntegrator(new mfem::DomainLFIntegrator(rhs_coeff)); 135*182fbe45STzanio b.Assemble(); 136*182fbe45STzanio 137*182fbe45STzanio // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using 138*182fbe45STzanio // the 'fespace' object to extract data needed by the Ceed objects. 139*182fbe45STzanio CeedDiffusionOperator diff(ceed, fespace); 140*182fbe45STzanio 141*182fbe45STzanio mfem::Operator *D; 142*182fbe45STzanio mfem::Vector X, B; 143*182fbe45STzanio diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B); 144*182fbe45STzanio 145*182fbe45STzanio // 8. Solve the discrete system using the conjugate gradients (CG) method. 146*182fbe45STzanio mfem::CGSolver cg; 147*182fbe45STzanio cg.SetRelTol(1e-6); 148*182fbe45STzanio cg.SetMaxIter(1000); 149*182fbe45STzanio cg.SetPrintLevel(3); 150*182fbe45STzanio cg.SetOperator(*D); 151*182fbe45STzanio 152*182fbe45STzanio cg.Mult(B, X); 153*182fbe45STzanio 154*182fbe45STzanio // 9. Compute and print the L2 norm of the error. 155*182fbe45STzanio std::cout << "L2 norm of the error: " << sol.ComputeL2Error(sol_coeff) 156*182fbe45STzanio << std::endl; 157*182fbe45STzanio 158*182fbe45STzanio // 10. Open a socket connection to GLVis and send the mesh and solution for 159*182fbe45STzanio // visualization. 160*182fbe45STzanio if (visualization) { 161*182fbe45STzanio char vishost[] = "localhost"; 162*182fbe45STzanio int visport = 19916; 163*182fbe45STzanio mfem::socketstream sol_sock(vishost, visport); 164*182fbe45STzanio sol_sock.precision(8); 165*182fbe45STzanio sol_sock << "solution\n" << *mesh << sol << std::flush; 166*182fbe45STzanio } 167*182fbe45STzanio 168*182fbe45STzanio // 11. Free memory and exit. 169*182fbe45STzanio delete fespace; 170*182fbe45STzanio delete fec; 171*182fbe45STzanio delete mesh; 172*182fbe45STzanio CeedDestroy(&ceed); 173*182fbe45STzanio return 0; 174*182fbe45STzanio } 175