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