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 // 24182fbe45STzanio // bp3 25182fbe45STzanio // bp3 -ceed /cpu/self 26182fbe45STzanio // bp3 -m ../../../mfem/data/fichera.mesh -o 4 27182fbe45STzanio // bp3 -m ../../../mfem/data/square-disc-nurbs.mesh -o 6 28182fbe45STzanio // bp3 -m ../../../mfem/data/inline-segment.mesh -o 8 29182fbe45STzanio 30182fbe45STzanio #include <ceed.h> 31182fbe45STzanio #include <mfem.hpp> 32*c0c38e35SVeselin Dobrev #include "bp3.hpp" 33182fbe45STzanio 34182fbe45STzanio /// Exact solution 35182fbe45STzanio double solution(const mfem::Vector &pt) { 36182fbe45STzanio static const double x[3] = { -0.32, 0.15, 0.24 }; 37182fbe45STzanio static const double k[3] = { 1.21, 1.45, 1.37 }; 38182fbe45STzanio double val = sin(M_PI*(x[0]+k[0]*pt(0))); 39182fbe45STzanio for (int d = 1; d < pt.Size(); d++) 40182fbe45STzanio val *= sin(M_PI*(x[d]+k[d]*pt(d))); 41182fbe45STzanio return val; 42182fbe45STzanio } 43182fbe45STzanio 44182fbe45STzanio /// Right-hand side 45182fbe45STzanio double rhs(const mfem::Vector &pt) { 46182fbe45STzanio static const double x[3] = { -0.32, 0.15, 0.24 }; 47182fbe45STzanio static const double k[3] = { 1.21, 1.45, 1.37 }; 48182fbe45STzanio double f[3], l[3], val, lap; 49182fbe45STzanio f[0] = sin(M_PI*(x[0]+k[0]*pt(0))); 50182fbe45STzanio l[0] = M_PI*M_PI*k[0]*k[0]*f[0]; 51182fbe45STzanio val = f[0]; 52182fbe45STzanio lap = l[0]; 53182fbe45STzanio for (int d = 1; d < pt.Size(); d++) { 54182fbe45STzanio f[d] = sin(M_PI*(x[d]+k[d]*pt(d))); 55182fbe45STzanio l[d] = M_PI*M_PI*k[d]*k[d]*f[d]; 56182fbe45STzanio lap = lap*f[d] + val*l[d]; 57182fbe45STzanio val = val*f[d]; 58182fbe45STzanio } 59182fbe45STzanio return lap; 60182fbe45STzanio } 61182fbe45STzanio 62182fbe45STzanio 63182fbe45STzanio int main(int argc, char *argv[]) { 64182fbe45STzanio // 1. Parse command-line options. 65182fbe45STzanio const char *ceed_spec = "/cpu/self"; 66*c0c38e35SVeselin Dobrev #ifndef MFEM_DIR 67182fbe45STzanio const char *mesh_file = "../../../mfem/data/star.mesh"; 68*c0c38e35SVeselin Dobrev #else 69*c0c38e35SVeselin Dobrev const char *mesh_file = MFEM_DIR "/data/star.mesh"; 70*c0c38e35SVeselin Dobrev #endif 71182fbe45STzanio int order = 2; 72182fbe45STzanio bool visualization = true; 73182fbe45STzanio 74182fbe45STzanio mfem::OptionsParser args(argc, argv); 75182fbe45STzanio args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 76182fbe45STzanio args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 77182fbe45STzanio args.AddOption(&order, "-o", "--order", 78182fbe45STzanio "Finite element order (polynomial degree)."); 79182fbe45STzanio args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", 80182fbe45STzanio "--no-visualization", 81182fbe45STzanio "Enable or disable GLVis visualization."); 82182fbe45STzanio args.Parse(); 83182fbe45STzanio if (!args.Good()) { 84182fbe45STzanio args.PrintUsage(std::cout); 85182fbe45STzanio return 1; 86182fbe45STzanio } 87182fbe45STzanio args.PrintOptions(std::cout); 88182fbe45STzanio 89182fbe45STzanio // 2. Initialize a Ceed device object using the given Ceed specification. 90182fbe45STzanio Ceed ceed; 91182fbe45STzanio CeedInit(ceed_spec, &ceed); 92182fbe45STzanio 93182fbe45STzanio // 3. Read the mesh from the given mesh file. 94182fbe45STzanio mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1); 95182fbe45STzanio int dim = mesh->Dimension(); 96182fbe45STzanio 97182fbe45STzanio // 4. Refine the mesh to increase the resolution. In this example we do 98182fbe45STzanio // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the 99182fbe45STzanio // largest number that gives a final system with no more than 50,000 (1,000 100182fbe45STzanio // in 1D) unknowns, approximately. 101182fbe45STzanio { 102182fbe45STzanio double max_dofs = (dim > 1) ? 50000 : 1000; 103182fbe45STzanio int ref_levels = 104182fbe45STzanio (int)floor((log(max_dofs/mesh->GetNE())-dim*log(order))/log(2.)/dim); 105182fbe45STzanio for (int l = 0; l < ref_levels; l++) { 106182fbe45STzanio mesh->UniformRefinement(); 107182fbe45STzanio } 108182fbe45STzanio } 109182fbe45STzanio if (mesh->GetNodalFESpace() == NULL) { 110182fbe45STzanio mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 111182fbe45STzanio } 112182fbe45STzanio if (mesh->NURBSext) { 113182fbe45STzanio mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 114182fbe45STzanio } 115182fbe45STzanio 116182fbe45STzanio // 5. Define a finite element space on the mesh. Here we use continuous 117182fbe45STzanio // Lagrange finite elements of the specified order. 118182fbe45STzanio MFEM_VERIFY(order > 0, "invalid order"); 119182fbe45STzanio mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 120182fbe45STzanio mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 121182fbe45STzanio std::cout << "Number of finite element unknowns: " 122182fbe45STzanio << fespace->GetTrueVSize() << std::endl; 123182fbe45STzanio 124182fbe45STzanio mfem::FunctionCoefficient sol_coeff(solution); 125182fbe45STzanio mfem::Array<int> ess_tdof_list; 126182fbe45STzanio mfem::GridFunction sol(fespace); 127182fbe45STzanio if (mesh->bdr_attributes.Size()) { 128182fbe45STzanio mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max()); 129182fbe45STzanio ess_bdr = 1; 130182fbe45STzanio fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list); 131182fbe45STzanio sol.ProjectBdrCoefficient(sol_coeff, ess_bdr); 132182fbe45STzanio } 133182fbe45STzanio 134182fbe45STzanio // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where 135182fbe45STzanio // v is a test function. 136182fbe45STzanio mfem::LinearForm b(fespace); 137182fbe45STzanio mfem::FunctionCoefficient rhs_coeff(rhs); 138182fbe45STzanio b.AddDomainIntegrator(new mfem::DomainLFIntegrator(rhs_coeff)); 139182fbe45STzanio b.Assemble(); 140182fbe45STzanio 141182fbe45STzanio // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using 142182fbe45STzanio // the 'fespace' object to extract data needed by the Ceed objects. 143182fbe45STzanio CeedDiffusionOperator diff(ceed, fespace); 144182fbe45STzanio 145182fbe45STzanio mfem::Operator *D; 146182fbe45STzanio mfem::Vector X, B; 147182fbe45STzanio diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B); 148182fbe45STzanio 149182fbe45STzanio // 8. Solve the discrete system using the conjugate gradients (CG) method. 150182fbe45STzanio mfem::CGSolver cg; 151182fbe45STzanio cg.SetRelTol(1e-6); 152182fbe45STzanio cg.SetMaxIter(1000); 153182fbe45STzanio cg.SetPrintLevel(3); 154182fbe45STzanio cg.SetOperator(*D); 155182fbe45STzanio 156182fbe45STzanio cg.Mult(B, X); 157182fbe45STzanio 158182fbe45STzanio // 9. Compute and print the L2 norm of the error. 159182fbe45STzanio std::cout << "L2 norm of the error: " << sol.ComputeL2Error(sol_coeff) 160182fbe45STzanio << std::endl; 161182fbe45STzanio 162182fbe45STzanio // 10. Open a socket connection to GLVis and send the mesh and solution for 163182fbe45STzanio // visualization. 164182fbe45STzanio if (visualization) { 165182fbe45STzanio char vishost[] = "localhost"; 166182fbe45STzanio int visport = 19916; 167182fbe45STzanio mfem::socketstream sol_sock(vishost, visport); 168182fbe45STzanio sol_sock.precision(8); 169182fbe45STzanio sol_sock << "solution\n" << *mesh << sol << std::flush; 170182fbe45STzanio } 171182fbe45STzanio 172182fbe45STzanio // 11. Free memory and exit. 173182fbe45STzanio delete fespace; 174182fbe45STzanio delete fec; 175182fbe45STzanio delete mesh; 176182fbe45STzanio CeedDestroy(&ceed); 177182fbe45STzanio return 0; 178182fbe45STzanio } 179