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 -m ../../../mfem/data/fichera.mesh -o 4 27 // bp3 -m ../../../mfem/data/square-disc-nurbs.mesh -o 6 28 // bp3 -m ../../../mfem/data/inline-segment.mesh -o 8 29 30 #include <ceed.h> 31 #include <mfem.hpp> 32 #include "bp3.hpp" 33 34 /// Exact solution 35 double solution(const mfem::Vector &pt) { 36 static const double x[3] = { -0.32, 0.15, 0.24 }; 37 static const double k[3] = { 1.21, 1.45, 1.37 }; 38 double val = sin(M_PI*(x[0]+k[0]*pt(0))); 39 for (int d = 1; d < pt.Size(); d++) 40 val *= sin(M_PI*(x[d]+k[d]*pt(d))); 41 return val; 42 } 43 44 /// Right-hand side 45 double rhs(const mfem::Vector &pt) { 46 static const double x[3] = { -0.32, 0.15, 0.24 }; 47 static const double k[3] = { 1.21, 1.45, 1.37 }; 48 double f[3], l[3], val, lap; 49 f[0] = sin(M_PI*(x[0]+k[0]*pt(0))); 50 l[0] = M_PI*M_PI*k[0]*k[0]*f[0]; 51 val = f[0]; 52 lap = l[0]; 53 for (int d = 1; d < pt.Size(); d++) { 54 f[d] = sin(M_PI*(x[d]+k[d]*pt(d))); 55 l[d] = M_PI*M_PI*k[d]*k[d]*f[d]; 56 lap = lap*f[d] + val*l[d]; 57 val = val*f[d]; 58 } 59 return lap; 60 } 61 62 //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000 63 int main(int argc, char *argv[]) { 64 // 1. Parse command-line options. 65 const char *ceed_spec = "/cpu/self"; 66 #ifndef MFEM_DIR 67 const char *mesh_file = "../../../mfem/data/star.mesh"; 68 #else 69 const char *mesh_file = MFEM_DIR "/data/star.mesh"; 70 #endif 71 int order = 2; 72 bool visualization = true; 73 bool test = false; 74 double max_dofs = 50000; 75 76 mfem::OptionsParser args(argc, argv); 77 args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification."); 78 args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use."); 79 args.AddOption(&order, "-o", "--order", 80 "Finite element order (polynomial degree)."); 81 args.AddOption(&max_dofs, "-s", "--size", "Maximum size (number of DoFs)"); 82 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", 83 "--no-visualization", 84 "Enable or disable GLVis visualization."); 85 args.AddOption(&test, "-t", "--test", "-no-test", 86 "--no-test", 87 "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 = 111 (int)floor((log(max_dofs/mesh->GetNE())-dim*log(order))/log(2.)/dim); 112 for (int l = 0; l < ref_levels; l++) { 113 mesh->UniformRefinement(); 114 } 115 } 116 if (mesh->GetNodalFESpace() == NULL) { 117 mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES); 118 } 119 if (mesh->NURBSext) { 120 mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES); 121 } 122 123 // 5. Define a finite element space on the mesh. Here we use continuous 124 // Lagrange finite elements of the specified order. 125 MFEM_VERIFY(order > 0, "invalid order"); 126 mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim); 127 mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec); 128 if (!test) { 129 std::cout << "Number of finite element unknowns: " 130 << fespace->GetTrueVSize() << std::endl; 131 } 132 133 mfem::FunctionCoefficient sol_coeff(solution); 134 mfem::Array<int> ess_tdof_list; 135 mfem::GridFunction sol(fespace); 136 if (mesh->bdr_attributes.Size()) { 137 mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max()); 138 ess_bdr = 1; 139 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list); 140 sol.ProjectBdrCoefficient(sol_coeff, ess_bdr); 141 } 142 143 // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where 144 // v is a test function. 145 mfem::LinearForm b(fespace); 146 mfem::FunctionCoefficient rhs_coeff(rhs); 147 b.AddDomainIntegrator(new mfem::DomainLFIntegrator(rhs_coeff)); 148 b.Assemble(); 149 150 // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using 151 // the 'fespace' object to extract data needed by the Ceed objects. 152 CeedDiffusionOperator diff(ceed, fespace); 153 154 mfem::Operator *D; 155 mfem::Vector X, B; 156 diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B); 157 158 // 8. Solve the discrete system using the conjugate gradients (CG) method. 159 mfem::CGSolver cg; 160 cg.SetRelTol(1e-6); 161 cg.SetMaxIter(1000); 162 if (test) { 163 cg.SetPrintLevel(0); 164 } else { 165 cg.SetPrintLevel(3); 166 } 167 cg.SetOperator(*D); 168 169 cg.Mult(B, X); 170 171 // 9. Compute and print the L2 norm of the error. 172 if (!test) { 173 std::cout << "L2 projection error: " << sol.ComputeL2Error(sol_coeff) 174 << std::endl; 175 } else { 176 if (fabs(sol.ComputeL2Error(sol_coeff))>2e-3) { 177 std::cout << "Error too large" << std::endl; 178 } 179 } 180 181 // 10. Open a socket connection to GLVis and send the mesh and solution for 182 // visualization. 183 if (visualization) { 184 char vishost[] = "localhost"; 185 int visport = 19916; 186 mfem::socketstream sol_sock(vishost, visport); 187 sol_sock.precision(8); 188 sol_sock << "solution\n" << *mesh << sol << std::flush; 189 } 190 191 // 11. Free memory and exit. 192 delete fespace; 193 delete fec; 194 delete mesh; 195 CeedDestroy(&ceed); 196 return 0; 197 } 198