xref: /libCEED/examples/mfem/bp1.cpp (revision dc00e230834c021cd7c5b556bbaba87dc4209676)
1182fbe45STzanio //                         libCEED + MFEM Example: BP1
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 simple linear system with a
7182fbe45STzanio // mass matrix (L2-projection of a given analytic function provided by
8182fbe45STzanio // 'solution'). The mass matrix required for performing the projection is
9182fbe45STzanio // expressed as a new class, CeedMassOperator, derived from mfem::Operator.
10182fbe45STzanio // Internally, CeedMassOperator uses a CeedOperator object constructed based on
11182fbe45STzanio // an mfem::FiniteElementSpace. All libCEED objects use a Ceed device object
12182fbe45STzanio // constructed based on a command line argument (-ceed).
13182fbe45STzanio //
14182fbe45STzanio // The mass matrix is inverted using a simple conjugate gradient algorithm
15182fbe45STzanio // corresponding to CEED BP1, 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 bp1 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>]
21182fbe45STzanio //
22182fbe45STzanio // Sample runs:
23182fbe45STzanio //
24182fbe45STzanio //     bp1
25182fbe45STzanio //     bp1 -ceed /cpu/self
26182fbe45STzanio //     bp1 -ceed /gpu/occa
27182fbe45STzanio //     bp1 -ceed /cpu/occa
28182fbe45STzanio //     bp1 -ceed /omp/occa
29182fbe45STzanio //     bp1 -ceed /ocl/occa
30182fbe45STzanio //     bp1 -m ../../../mfem/data/fichera.mesh
31182fbe45STzanio //     bp1 -m ../../../mfem/data/star.vtk -o 3
32182fbe45STzanio //     bp1 -m ../../../mfem/data/inline-segment.mesh -o 8
33182fbe45STzanio 
34182fbe45STzanio #include <ceed.h>
35182fbe45STzanio #include <mfem.hpp>
36c0c38e35SVeselin Dobrev #include "bp1.hpp"
37182fbe45STzanio 
38182fbe45STzanio /// Continuous function to project on the discrete FE space
39182fbe45STzanio double solution(const mfem::Vector &pt) {
40182fbe45STzanio   return pt.Norml2(); // distance to the origin
41182fbe45STzanio }
42182fbe45STzanio 
43182fbe45STzanio 
44182fbe45STzanio int main(int argc, char *argv[]) {
45182fbe45STzanio   // 1. Parse command-line options.
46182fbe45STzanio   const char *ceed_spec = "/cpu/self";
47c0c38e35SVeselin Dobrev #ifndef MFEM_DIR
48182fbe45STzanio   const char *mesh_file = "../../../mfem/data/star.mesh";
49c0c38e35SVeselin Dobrev #else
50c0c38e35SVeselin Dobrev   const char *mesh_file = MFEM_DIR "/data/star.mesh";
51c0c38e35SVeselin Dobrev #endif
52182fbe45STzanio   int order = 1;
53182fbe45STzanio   bool visualization = true;
54*dc00e230Sjeremylt   bool test = false;
55182fbe45STzanio 
56182fbe45STzanio   mfem::OptionsParser args(argc, argv);
57182fbe45STzanio   args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification.");
58182fbe45STzanio   args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
59182fbe45STzanio   args.AddOption(&order, "-o", "--order",
60182fbe45STzanio                  "Finite element order (polynomial degree).");
61182fbe45STzanio   args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
62182fbe45STzanio                  "--no-visualization",
63182fbe45STzanio                  "Enable or disable GLVis visualization.");
64*dc00e230Sjeremylt   args.AddOption(&test, "-t", "--test", "-no-test",
65*dc00e230Sjeremylt                  "--no-test",
66*dc00e230Sjeremylt                  "Enable or disable test mode.");
67182fbe45STzanio   args.Parse();
68182fbe45STzanio   if (!args.Good()) {
69182fbe45STzanio     args.PrintUsage(std::cout);
70182fbe45STzanio     return 1;
71182fbe45STzanio   }
72*dc00e230Sjeremylt   if (!test) {
73182fbe45STzanio     args.PrintOptions(std::cout);
74*dc00e230Sjeremylt   }
75182fbe45STzanio 
76182fbe45STzanio   // 2. Initialize a Ceed device object using the given Ceed specification.
77182fbe45STzanio   Ceed ceed;
78182fbe45STzanio   CeedInit(ceed_spec, &ceed);
79182fbe45STzanio 
80182fbe45STzanio   // 3. Read the mesh from the given mesh file.
81182fbe45STzanio   mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1);
82182fbe45STzanio   int dim = mesh->Dimension();
83182fbe45STzanio 
84182fbe45STzanio   // 4. Refine the mesh to increase the resolution. In this example we do
85182fbe45STzanio   //    'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
86182fbe45STzanio   //    largest number that gives a final system with no more than 50,000
87182fbe45STzanio   //    unknowns, approximately.
88182fbe45STzanio   {
89182fbe45STzanio     double max_dofs = 50000;
90182fbe45STzanio     int ref_levels =
91182fbe45STzanio       (int)floor((log(max_dofs/mesh->GetNE())-dim*log(order))/log(2.)/dim);
92182fbe45STzanio     for (int l = 0; l < ref_levels; l++) {
93182fbe45STzanio       mesh->UniformRefinement();
94182fbe45STzanio     }
95182fbe45STzanio   }
96182fbe45STzanio   if (mesh->GetNodalFESpace() == NULL) {
97182fbe45STzanio     mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES);
98182fbe45STzanio   }
99182fbe45STzanio   if (mesh->NURBSext) {
100182fbe45STzanio     mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES);
101182fbe45STzanio   }
102182fbe45STzanio 
103182fbe45STzanio   // 5. Define a finite element space on the mesh. Here we use continuous
104182fbe45STzanio   //    Lagrange finite elements of the specified order.
105182fbe45STzanio   MFEM_VERIFY(order > 0, "invalid order");
106182fbe45STzanio   mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim);
107182fbe45STzanio   mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec);
108*dc00e230Sjeremylt   if (!test) {
109182fbe45STzanio     std::cout << "Number of finite element unknowns: "
110182fbe45STzanio               << fespace->GetTrueVSize() << std::endl;
111*dc00e230Sjeremylt   }
112182fbe45STzanio 
113182fbe45STzanio   // 6. Construct a rhs vector using the linear form f(v) = (solution, v), where
114182fbe45STzanio   //    v is a test function.
115182fbe45STzanio   mfem::LinearForm b(fespace);
116182fbe45STzanio   mfem::FunctionCoefficient sol_coeff(solution);
117182fbe45STzanio   b.AddDomainIntegrator(new mfem::DomainLFIntegrator(sol_coeff));
118182fbe45STzanio   b.Assemble();
119182fbe45STzanio 
120182fbe45STzanio   // 7. Construct a CeedMassOperator utilizing the 'ceed' device and using the
121182fbe45STzanio   //    'fespace' object to extract data needed by the Ceed objects.
122182fbe45STzanio   CeedMassOperator mass(ceed, fespace);
123182fbe45STzanio 
124182fbe45STzanio   // 8. Solve the discrete system using the conjugate gradients (CG) method.
125182fbe45STzanio   mfem::CGSolver cg;
126182fbe45STzanio   cg.SetRelTol(1e-6);
127182fbe45STzanio   cg.SetMaxIter(100);
128*dc00e230Sjeremylt   if (test) {
129*dc00e230Sjeremylt     cg.SetPrintLevel(0);
130*dc00e230Sjeremylt   } else {
131182fbe45STzanio     cg.SetPrintLevel(3);
132*dc00e230Sjeremylt   }
133182fbe45STzanio   cg.SetOperator(mass);
134182fbe45STzanio 
135182fbe45STzanio   mfem::GridFunction sol(fespace);
136182fbe45STzanio   sol = 0.0;
137182fbe45STzanio   cg.Mult(b, sol);
138182fbe45STzanio 
139182fbe45STzanio   // 9. Compute and print the L2 projection error.
140*dc00e230Sjeremylt   if (!test) {
141182fbe45STzanio     std::cout << "L2 projection error: " << sol.ComputeL2Error(sol_coeff)
142182fbe45STzanio               << std::endl;
143*dc00e230Sjeremylt   } else {
1449b872752Sjeremylt     if (fabs(sol.ComputeL2Error(sol_coeff))>1e-4) {
1459b872752Sjeremylt       std::cout << "Error too large" << std::endl;
146*dc00e230Sjeremylt     }
1479b872752Sjeremylt   }
148182fbe45STzanio 
149182fbe45STzanio   // 10. Open a socket connection to GLVis and send the mesh and solution for
150182fbe45STzanio   //     visualization.
151182fbe45STzanio   if (visualization) {
152182fbe45STzanio     char vishost[] = "localhost";
153182fbe45STzanio     int  visport   = 19916;
154182fbe45STzanio     mfem::socketstream sol_sock(vishost, visport);
155182fbe45STzanio     sol_sock.precision(8);
156182fbe45STzanio     sol_sock << "solution\n" << *mesh << sol << std::flush;
157182fbe45STzanio   }
158182fbe45STzanio 
159182fbe45STzanio   // 11. Free memory and exit.
160182fbe45STzanio   delete fespace;
161182fbe45STzanio   delete fec;
162182fbe45STzanio   delete mesh;
163182fbe45STzanio   CeedDestroy(&ceed);
164182fbe45STzanio   return 0;
165182fbe45STzanio }
166