xref: /libCEED/examples/mfem/bp1.cpp (revision 47a7da020d3cc99cd8a1634c544023154b55d335)
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 
44 int main(int argc, char *argv[]) {
45   // 1. Parse command-line options.
46   const char *ceed_spec = "/cpu/self";
47   const char *mesh_file = "../../../mfem/data/star.mesh";
48   int order = 1;
49   bool visualization = true;
50 
51   mfem::OptionsParser args(argc, argv);
52   args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification.");
53   args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
54   args.AddOption(&order, "-o", "--order",
55                  "Finite element order (polynomial degree).");
56   args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
57                  "--no-visualization",
58                  "Enable or disable GLVis visualization.");
59   args.Parse();
60   if (!args.Good()) {
61     args.PrintUsage(std::cout);
62     return 1;
63   }
64   args.PrintOptions(std::cout);
65 
66   // 2. Initialize a Ceed device object using the given Ceed specification.
67   Ceed ceed;
68   CeedInit(ceed_spec, &ceed);
69 
70   // 3. Read the mesh from the given mesh file.
71   mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1);
72   int dim = mesh->Dimension();
73 
74   // 4. Refine the mesh to increase the resolution. In this example we do
75   //    'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
76   //    largest number that gives a final system with no more than 50,000
77   //    unknowns, approximately.
78   {
79     double max_dofs = 50000;
80     int ref_levels =
81       (int)floor((log(max_dofs/mesh->GetNE())-dim*log(order))/log(2.)/dim);
82     for (int l = 0; l < ref_levels; l++) {
83       mesh->UniformRefinement();
84     }
85   }
86   if (mesh->GetNodalFESpace() == NULL) {
87     mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES);
88   }
89   if (mesh->NURBSext) {
90     mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES);
91   }
92 
93   // 5. Define a finite element space on the mesh. Here we use continuous
94   //    Lagrange finite elements of the specified order.
95   MFEM_VERIFY(order > 0, "invalid order");
96   mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim);
97   mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec);
98   std::cout << "Number of finite element unknowns: "
99             << fespace->GetTrueVSize() << std::endl;
100 
101   // 6. Construct a rhs vector using the linear form f(v) = (solution, v), where
102   //    v is a test function.
103   mfem::LinearForm b(fespace);
104   mfem::FunctionCoefficient sol_coeff(solution);
105   b.AddDomainIntegrator(new mfem::DomainLFIntegrator(sol_coeff));
106   b.Assemble();
107 
108   // 7. Construct a CeedMassOperator utilizing the 'ceed' device and using the
109   //    'fespace' object to extract data needed by the Ceed objects.
110   CeedMassOperator mass(ceed, fespace);
111 
112   // 8. Solve the discrete system using the conjugate gradients (CG) method.
113   mfem::CGSolver cg;
114   cg.SetRelTol(1e-6);
115   cg.SetMaxIter(100);
116   cg.SetPrintLevel(3);
117   cg.SetOperator(mass);
118 
119   mfem::GridFunction sol(fespace);
120   sol = 0.0;
121   cg.Mult(b, sol);
122 
123   // 9. Compute and print the L2 projection error.
124   std::cout << "L2 projection error: " << sol.ComputeL2Error(sol_coeff)
125             << std::endl;
126 
127   // 10. Open a socket connection to GLVis and send the mesh and solution for
128   //     visualization.
129   if (visualization) {
130     char vishost[] = "localhost";
131     int  visport   = 19916;
132     mfem::socketstream sol_sock(vishost, visport);
133     sol_sock.precision(8);
134     sol_sock << "solution\n" << *mesh << sol << std::flush;
135   }
136 
137   // 11. Free memory and exit.
138   delete fespace;
139   delete fec;
140   delete mesh;
141   CeedDestroy(&ceed);
142   return 0;
143 }
144