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