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