xref: /libCEED/examples/mfem/bp3.cpp (revision 286887987a5e31fe81687a1110e0e2b7aea371d1)
1182fbe45STzanio //                         libCEED + MFEM Example: BP3
2182fbe45STzanio //
3182fbe45STzanio // This example illustrates a simple usage of libCEED with the MFEM (mfem.org)
4182fbe45STzanio // finite element library.
5182fbe45STzanio //
6182fbe45STzanio // The example reads a mesh from a file and solves a linear system with a
7182fbe45STzanio // diffusion stiffness matrix (with a prescribed analytic solution, provided by
8182fbe45STzanio // the function 'solution'). The diffusion matrix is expressed as a new class,
9182fbe45STzanio // CeedDiffusionOperator, derived from mfem::Operator. Internally,
10182fbe45STzanio // CeedDiffusionOperator uses a CeedOperator object constructed based on an
11182fbe45STzanio // mfem::FiniteElementSpace. All libCEED objects use a Ceed logical device
12182fbe45STzanio // object constructed based on a command line argument. (-ceed).
13182fbe45STzanio //
14182fbe45STzanio // The linear system is inverted using the conjugate gradients algorithm
15182fbe45STzanio // corresponding to CEED BP3, see http://ceed.exascaleproject.org/bps. Arbitrary
16182fbe45STzanio // mesh and solution orders in 1D, 2D and 3D are supported from the same code.
17182fbe45STzanio //
18182fbe45STzanio // Build with:
19182fbe45STzanio //
20182fbe45STzanio //     make bp3 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>]
21182fbe45STzanio //
22182fbe45STzanio // Sample runs:
23182fbe45STzanio //
2466087c08SValeria Barra //     ./bp3
2566087c08SValeria Barra //     ./bp3 -ceed /cpu/self
26*28688798Sjeremylt //     ./bp3 -ceed /gpu/cuda
2766087c08SValeria Barra //     ./bp3 -m ../../../mfem/data/fichera.mesh -o 4
2866087c08SValeria Barra //     ./bp3 -m ../../../mfem/data/square-disc-nurbs.mesh -o 6
2966087c08SValeria Barra //     ./bp3 -m ../../../mfem/data/inline-segment.mesh -o 8
30182fbe45STzanio 
315d6bafb2Sjeremylt /// @file
325d6bafb2Sjeremylt /// MFEM diffusion operator based on libCEED
335d6bafb2Sjeremylt 
34182fbe45STzanio #include <ceed.h>
35182fbe45STzanio #include <mfem.hpp>
36c0c38e35SVeselin Dobrev #include "bp3.hpp"
37182fbe45STzanio 
38182fbe45STzanio /// Exact solution
39182fbe45STzanio double solution(const mfem::Vector &pt) {
40182fbe45STzanio   static const double x[3] = { -0.32, 0.15, 0.24 };
41182fbe45STzanio   static const double k[3] = { 1.21, 1.45, 1.37 };
42182fbe45STzanio   double val = sin(M_PI*(x[0]+k[0]*pt(0)));
43182fbe45STzanio   for (int d = 1; d < pt.Size(); d++)
44182fbe45STzanio     val *= sin(M_PI*(x[d]+k[d]*pt(d)));
45182fbe45STzanio   return val;
46182fbe45STzanio }
47182fbe45STzanio 
48182fbe45STzanio /// Right-hand side
49182fbe45STzanio double rhs(const mfem::Vector &pt) {
50182fbe45STzanio   static const double x[3] = { -0.32, 0.15, 0.24 };
51182fbe45STzanio   static const double k[3] = { 1.21, 1.45, 1.37 };
52182fbe45STzanio   double f[3], l[3], val, lap;
53182fbe45STzanio   f[0] = sin(M_PI*(x[0]+k[0]*pt(0)));
54182fbe45STzanio   l[0] = M_PI*M_PI*k[0]*k[0]*f[0];
55182fbe45STzanio   val = f[0];
56182fbe45STzanio   lap = l[0];
57182fbe45STzanio   for (int d = 1; d < pt.Size(); d++) {
58182fbe45STzanio     f[d] = sin(M_PI*(x[d]+k[d]*pt(d)));
59182fbe45STzanio     l[d] = M_PI*M_PI*k[d]*k[d]*f[d];
60182fbe45STzanio     lap = lap*f[d] + val*l[d];
61182fbe45STzanio     val = val*f[d];
62182fbe45STzanio   }
63182fbe45STzanio   return lap;
64182fbe45STzanio }
65182fbe45STzanio 
66f063656dSJed Brown //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000
67182fbe45STzanio int main(int argc, char *argv[]) {
68182fbe45STzanio   // 1. Parse command-line options.
69182fbe45STzanio   const char *ceed_spec = "/cpu/self";
70c0c38e35SVeselin Dobrev   #ifndef MFEM_DIR
71182fbe45STzanio   const char *mesh_file = "../../../mfem/data/star.mesh";
72c0c38e35SVeselin Dobrev   #else
73c0c38e35SVeselin Dobrev   const char *mesh_file = MFEM_DIR "/data/star.mesh";
74c0c38e35SVeselin Dobrev   #endif
75182fbe45STzanio   int order = 2;
76182fbe45STzanio   bool visualization = true;
77dc00e230Sjeremylt   bool test = false;
78e2b2c771Svaleria   double max_nnodes = 50000;
79182fbe45STzanio 
80182fbe45STzanio   mfem::OptionsParser args(argc, argv);
81182fbe45STzanio   args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification.");
82182fbe45STzanio   args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
83182fbe45STzanio   args.AddOption(&order, "-o", "--order",
84182fbe45STzanio                  "Finite element order (polynomial degree).");
85e2b2c771Svaleria   args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)");
86182fbe45STzanio   args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
87182fbe45STzanio                  "--no-visualization",
88182fbe45STzanio                  "Enable or disable GLVis visualization.");
89dc00e230Sjeremylt   args.AddOption(&test, "-t", "--test", "-no-test",
90dc00e230Sjeremylt                  "--no-test",
91dc00e230Sjeremylt                  "Enable or disable test mode.");
92182fbe45STzanio   args.Parse();
93182fbe45STzanio   if (!args.Good()) {
94182fbe45STzanio     args.PrintUsage(std::cout);
95182fbe45STzanio     return 1;
96182fbe45STzanio   }
97dc00e230Sjeremylt   if (!test) {
98182fbe45STzanio     args.PrintOptions(std::cout);
99dc00e230Sjeremylt   }
100182fbe45STzanio 
101182fbe45STzanio   // 2. Initialize a Ceed device object using the given Ceed specification.
102182fbe45STzanio   Ceed ceed;
103182fbe45STzanio   CeedInit(ceed_spec, &ceed);
104182fbe45STzanio 
105182fbe45STzanio   // 3. Read the mesh from the given mesh file.
106182fbe45STzanio   mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1);
107182fbe45STzanio   int dim = mesh->Dimension();
108182fbe45STzanio 
109182fbe45STzanio   // 4. Refine the mesh to increase the resolution. In this example we do
110182fbe45STzanio   //    'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
11131d4d2baSJed Brown   //    largest number that gives a final system with no more than 50,000
11231d4d2baSJed Brown   //    unknowns, approximately.
113182fbe45STzanio   {
114182fbe45STzanio     int ref_levels =
115e2b2c771Svaleria       (int)floor((log(max_nnodes/mesh->GetNE())-dim*log(order))/log(2.)/dim);
116182fbe45STzanio     for (int l = 0; l < ref_levels; l++) {
117182fbe45STzanio       mesh->UniformRefinement();
118182fbe45STzanio     }
119182fbe45STzanio   }
120182fbe45STzanio   if (mesh->GetNodalFESpace() == NULL) {
121182fbe45STzanio     mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES);
122182fbe45STzanio   }
123182fbe45STzanio   if (mesh->NURBSext) {
124182fbe45STzanio     mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES);
125182fbe45STzanio   }
126182fbe45STzanio 
127182fbe45STzanio   // 5. Define a finite element space on the mesh. Here we use continuous
128182fbe45STzanio   //    Lagrange finite elements of the specified order.
129182fbe45STzanio   MFEM_VERIFY(order > 0, "invalid order");
130182fbe45STzanio   mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim);
131182fbe45STzanio   mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec);
132dc00e230Sjeremylt   if (!test) {
133182fbe45STzanio     std::cout << "Number of finite element unknowns: "
134182fbe45STzanio               << fespace->GetTrueVSize() << std::endl;
135dc00e230Sjeremylt   }
136182fbe45STzanio 
137182fbe45STzanio   mfem::FunctionCoefficient sol_coeff(solution);
138182fbe45STzanio   mfem::Array<int> ess_tdof_list;
139182fbe45STzanio   mfem::GridFunction sol(fespace);
140182fbe45STzanio   if (mesh->bdr_attributes.Size()) {
141182fbe45STzanio     mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max());
142182fbe45STzanio     ess_bdr = 1;
143182fbe45STzanio     fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
144182fbe45STzanio     sol.ProjectBdrCoefficient(sol_coeff, ess_bdr);
145182fbe45STzanio   }
146182fbe45STzanio 
147182fbe45STzanio   // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where
148182fbe45STzanio   //    v is a test function.
149182fbe45STzanio   mfem::LinearForm b(fespace);
150182fbe45STzanio   mfem::FunctionCoefficient rhs_coeff(rhs);
151182fbe45STzanio   b.AddDomainIntegrator(new mfem::DomainLFIntegrator(rhs_coeff));
152182fbe45STzanio   b.Assemble();
153182fbe45STzanio 
154182fbe45STzanio   // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using
155182fbe45STzanio   //    the 'fespace' object to extract data needed by the Ceed objects.
156182fbe45STzanio   CeedDiffusionOperator diff(ceed, fespace);
157182fbe45STzanio 
158182fbe45STzanio   mfem::Operator *D;
159182fbe45STzanio   mfem::Vector X, B;
160182fbe45STzanio   diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B);
161182fbe45STzanio 
162182fbe45STzanio   // 8. Solve the discrete system using the conjugate gradients (CG) method.
163182fbe45STzanio   mfem::CGSolver cg;
164182fbe45STzanio   cg.SetRelTol(1e-6);
165182fbe45STzanio   cg.SetMaxIter(1000);
166dc00e230Sjeremylt   if (test) {
167dc00e230Sjeremylt     cg.SetPrintLevel(0);
168dc00e230Sjeremylt   } else {
169182fbe45STzanio     cg.SetPrintLevel(3);
170dc00e230Sjeremylt   }
171182fbe45STzanio   cg.SetOperator(*D);
172182fbe45STzanio 
173182fbe45STzanio   cg.Mult(B, X);
174182fbe45STzanio 
175182fbe45STzanio   // 9. Compute and print the L2 norm of the error.
1769647a07eSDavid Medina   double err_l2 = sol.ComputeL2Error(sol_coeff);
177dc00e230Sjeremylt   if (!test) {
1789647a07eSDavid Medina     std::cout << "L2 projection error: " << err_l2
179182fbe45STzanio               << std::endl;
180dc00e230Sjeremylt   } else {
181f063656dSJed Brown     if (fabs(sol.ComputeL2Error(sol_coeff))>2e-3) {
1829647a07eSDavid Medina       std::cout << "Error too large: " << err_l2 << std::endl;
183dc00e230Sjeremylt     }
184dc00e230Sjeremylt   }
185182fbe45STzanio 
186182fbe45STzanio   // 10. Open a socket connection to GLVis and send the mesh and solution for
187182fbe45STzanio   //     visualization.
188182fbe45STzanio   if (visualization) {
189182fbe45STzanio     char vishost[] = "localhost";
190182fbe45STzanio     int  visport   = 19916;
191182fbe45STzanio     mfem::socketstream sol_sock(vishost, visport);
192182fbe45STzanio     sol_sock.precision(8);
193182fbe45STzanio     sol_sock << "solution\n" << *mesh << sol << std::flush;
194182fbe45STzanio   }
195182fbe45STzanio 
196182fbe45STzanio   // 11. Free memory and exit.
197182fbe45STzanio   delete fespace;
198182fbe45STzanio   delete fec;
199182fbe45STzanio   delete mesh;
200182fbe45STzanio   CeedDestroy(&ceed);
201182fbe45STzanio   return 0;
202182fbe45STzanio }
203