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