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