xref: /libCEED/examples/mfem/bp3.cpp (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
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
2628688798Sjeremylt //     ./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 
34c0c38e35SVeselin Dobrev #include "bp3.hpp"
35182fbe45STzanio 
36*2b730f8bSJeremy L Thompson #include <ceed.h>
37*2b730f8bSJeremy L Thompson 
38*2b730f8bSJeremy L Thompson #include <mfem.hpp>
39*2b730f8bSJeremy L Thompson 
40182fbe45STzanio /// Exact solution
41182fbe45STzanio double solution(const mfem::Vector &pt) {
42182fbe45STzanio   static const double x[3] = {-0.32, 0.15, 0.24};
43182fbe45STzanio   static const double k[3] = {1.21, 1.45, 1.37};
44182fbe45STzanio   double              val  = sin(M_PI * (x[0] + k[0] * pt(0)));
45*2b730f8bSJeremy L Thompson   for (int d = 1; d < pt.Size(); d++) val *= sin(M_PI * (x[d] + k[d] * pt(d)));
46182fbe45STzanio   return val;
47182fbe45STzanio }
48182fbe45STzanio 
49182fbe45STzanio /// Right-hand side
50182fbe45STzanio double rhs(const mfem::Vector &pt) {
51182fbe45STzanio   static const double x[3] = {-0.32, 0.15, 0.24};
52182fbe45STzanio   static const double k[3] = {1.21, 1.45, 1.37};
53182fbe45STzanio   double              f[3], l[3], val, lap;
54182fbe45STzanio   f[0] = sin(M_PI * (x[0] + k[0] * pt(0)));
55182fbe45STzanio   l[0] = M_PI * M_PI * k[0] * k[0] * f[0];
56182fbe45STzanio   val  = f[0];
57182fbe45STzanio   lap  = l[0];
58182fbe45STzanio   for (int d = 1; d < pt.Size(); d++) {
59182fbe45STzanio     f[d] = sin(M_PI * (x[d] + k[d] * pt(d)));
60182fbe45STzanio     l[d] = M_PI * M_PI * k[d] * k[d] * f[d];
61182fbe45STzanio     lap  = lap * f[d] + val * l[d];
62182fbe45STzanio     val  = val * f[d];
63182fbe45STzanio   }
64182fbe45STzanio   return lap;
65182fbe45STzanio }
66182fbe45STzanio 
67f063656dSJed Brown //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000
68182fbe45STzanio int main(int argc, char *argv[]) {
69182fbe45STzanio   // 1. Parse command-line options.
70182fbe45STzanio   const char *ceed_spec = "/cpu/self";
71c0c38e35SVeselin Dobrev #ifndef MFEM_DIR
72182fbe45STzanio   const char *mesh_file = "../../../mfem/data/star.mesh";
73c0c38e35SVeselin Dobrev #else
74c0c38e35SVeselin Dobrev   const char *mesh_file = MFEM_DIR "/data/star.mesh";
75c0c38e35SVeselin Dobrev #endif
76182fbe45STzanio   int    order         = 2;
77182fbe45STzanio   bool   visualization = true;
78dc00e230Sjeremylt   bool   test          = false;
79e2b2c771Svaleria   double max_nnodes    = 50000;
80182fbe45STzanio 
81182fbe45STzanio   mfem::OptionsParser args(argc, argv);
82182fbe45STzanio   args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification.");
83182fbe45STzanio   args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
84*2b730f8bSJeremy L Thompson   args.AddOption(&order, "-o", "--order", "Finite element order (polynomial degree).");
85e2b2c771Svaleria   args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)");
86*2b730f8bSJeremy L Thompson   args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", "--no-visualization", "Enable or disable GLVis visualization.");
87*2b730f8bSJeremy L Thompson   args.AddOption(&test, "-t", "--test", "-no-test", "--no-test", "Enable or disable test mode.");
88182fbe45STzanio   args.Parse();
89182fbe45STzanio   if (!args.Good()) {
90182fbe45STzanio     args.PrintUsage(std::cout);
91182fbe45STzanio     return 1;
92182fbe45STzanio   }
93dc00e230Sjeremylt   if (!test) {
94182fbe45STzanio     args.PrintOptions(std::cout);
95dc00e230Sjeremylt   }
96182fbe45STzanio 
97182fbe45STzanio   // 2. Initialize a Ceed device object using the given Ceed specification.
98182fbe45STzanio   Ceed ceed;
99182fbe45STzanio   CeedInit(ceed_spec, &ceed);
100182fbe45STzanio 
101182fbe45STzanio   // 3. Read the mesh from the given mesh file.
102182fbe45STzanio   mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1);
103182fbe45STzanio   int         dim  = mesh->Dimension();
104182fbe45STzanio 
105182fbe45STzanio   // 4. Refine the mesh to increase the resolution. In this example we do
106182fbe45STzanio   //    'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
10731d4d2baSJed Brown   //    largest number that gives a final system with no more than 50,000
10831d4d2baSJed Brown   //    unknowns, approximately.
109182fbe45STzanio   {
110*2b730f8bSJeremy L Thompson     int ref_levels = (int)floor((log(max_nnodes / mesh->GetNE()) - dim * log(order)) / log(2.) / dim);
111182fbe45STzanio     for (int l = 0; l < ref_levels; l++) {
112182fbe45STzanio       mesh->UniformRefinement();
113182fbe45STzanio     }
114182fbe45STzanio   }
115182fbe45STzanio   if (mesh->GetNodalFESpace() == NULL) {
116182fbe45STzanio     mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES);
117182fbe45STzanio   }
118182fbe45STzanio   if (mesh->NURBSext) {
119182fbe45STzanio     mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES);
120182fbe45STzanio   }
121182fbe45STzanio 
122182fbe45STzanio   // 5. Define a finite element space on the mesh. Here we use continuous
123182fbe45STzanio   //    Lagrange finite elements of the specified order.
124182fbe45STzanio   MFEM_VERIFY(order > 0, "invalid order");
125182fbe45STzanio   mfem::FiniteElementCollection *fec     = new mfem::H1_FECollection(order, dim);
126182fbe45STzanio   mfem::FiniteElementSpace      *fespace = new mfem::FiniteElementSpace(mesh, fec);
127dc00e230Sjeremylt   if (!test) {
128*2b730f8bSJeremy L Thompson     std::cout << "Number of finite element unknowns: " << fespace->GetTrueVSize() << std::endl;
129dc00e230Sjeremylt   }
130182fbe45STzanio 
131182fbe45STzanio   mfem::FunctionCoefficient sol_coeff(solution);
132182fbe45STzanio   mfem::Array<int>          ess_tdof_list;
133182fbe45STzanio   mfem::GridFunction        sol(fespace);
134182fbe45STzanio   if (mesh->bdr_attributes.Size()) {
135182fbe45STzanio     mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max());
136182fbe45STzanio     ess_bdr = 1;
137182fbe45STzanio     fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
138182fbe45STzanio     sol.ProjectBdrCoefficient(sol_coeff, ess_bdr);
139182fbe45STzanio   }
140182fbe45STzanio 
141182fbe45STzanio   // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where
142182fbe45STzanio   //    v is a test function.
143182fbe45STzanio   mfem::LinearForm          b(fespace);
144182fbe45STzanio   mfem::FunctionCoefficient rhs_coeff(rhs);
145182fbe45STzanio   b.AddDomainIntegrator(new mfem::DomainLFIntegrator(rhs_coeff));
146182fbe45STzanio   b.Assemble();
147182fbe45STzanio 
148182fbe45STzanio   // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using
149182fbe45STzanio   //    the 'fespace' object to extract data needed by the Ceed objects.
150182fbe45STzanio   CeedDiffusionOperator diff(ceed, fespace);
151182fbe45STzanio 
152182fbe45STzanio   mfem::Operator *D;
153182fbe45STzanio   mfem::Vector    X, B;
154182fbe45STzanio   diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B);
155182fbe45STzanio 
156182fbe45STzanio   // 8. Solve the discrete system using the conjugate gradients (CG) method.
157182fbe45STzanio   mfem::CGSolver cg;
158182fbe45STzanio   cg.SetRelTol(1e-6);
159182fbe45STzanio   cg.SetMaxIter(1000);
160dc00e230Sjeremylt   if (test) {
161dc00e230Sjeremylt     cg.SetPrintLevel(0);
162dc00e230Sjeremylt   } else {
163182fbe45STzanio     cg.SetPrintLevel(3);
164dc00e230Sjeremylt   }
165182fbe45STzanio   cg.SetOperator(*D);
166182fbe45STzanio 
167182fbe45STzanio   cg.Mult(B, X);
168182fbe45STzanio 
169182fbe45STzanio   // 9. Compute and print the L2 norm of the error.
1709647a07eSDavid Medina   double err_l2 = sol.ComputeL2Error(sol_coeff);
171dc00e230Sjeremylt   if (!test) {
172*2b730f8bSJeremy L Thompson     std::cout << "L2 projection error: " << err_l2 << std::endl;
173dc00e230Sjeremylt   } else {
174f063656dSJed Brown     if (fabs(sol.ComputeL2Error(sol_coeff)) > 2e-3) {
1759647a07eSDavid Medina       std::cout << "Error too large: " << err_l2 << std::endl;
176dc00e230Sjeremylt     }
177dc00e230Sjeremylt   }
178182fbe45STzanio 
179182fbe45STzanio   // 10. Open a socket connection to GLVis and send the mesh and solution for
180182fbe45STzanio   //     visualization.
181182fbe45STzanio   if (visualization) {
182182fbe45STzanio     char               vishost[] = "localhost";
183182fbe45STzanio     int                visport   = 19916;
184182fbe45STzanio     mfem::socketstream sol_sock(vishost, visport);
185182fbe45STzanio     sol_sock.precision(8);
186182fbe45STzanio     sol_sock << "solution\n" << *mesh << sol << std::flush;
187182fbe45STzanio   }
188182fbe45STzanio 
189182fbe45STzanio   // 11. Free memory and exit.
190182fbe45STzanio   delete fespace;
191182fbe45STzanio   delete fec;
192182fbe45STzanio   delete mesh;
193182fbe45STzanio   CeedDestroy(&ceed);
194182fbe45STzanio   return 0;
195182fbe45STzanio }
196