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