1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2ea61e9acSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ea61e9acSJeremy L Thompson //
4ea61e9acSJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5ea61e9acSJeremy L Thompson //
6ea61e9acSJeremy L Thompson // This file is part of CEED: http://github.com/ceed
7ea61e9acSJeremy L Thompson
8182fbe45STzanio // libCEED + MFEM Example: BP3
9182fbe45STzanio //
10ea61e9acSJeremy L Thompson // This example illustrates a simple usage of libCEED with the MFEM (mfem.org) finite element library.
11182fbe45STzanio //
12ea61e9acSJeremy L Thompson // The example reads a mesh from a file and solves a linear system with a diffusion stiffness matrix (with a prescribed analytic solution, provided by
13ea61e9acSJeremy L Thompson // the function 'solution'). The diffusion matrix is expressed as a new class, CeedDiffusionOperator, derived from mfem::Operator. Internally,
14ea61e9acSJeremy L Thompson // CeedDiffusionOperator uses a CeedOperator object constructed based on an mfem::FiniteElementSpace. All libCEED objects use a Ceed logical device
15182fbe45STzanio // object constructed based on a command line argument. (-ceed).
16182fbe45STzanio //
17ea61e9acSJeremy L Thompson // The linear system is inverted using the conjugate gradients algorithm corresponding to CEED BP3, see http://ceed.exascaleproject.org/bps.
18ea61e9acSJeremy L Thompson // Arbitrary mesh and solution orders in 1D, 2D and 3D are supported from the same code.
19182fbe45STzanio //
20182fbe45STzanio // Build with:
21182fbe45STzanio //
22182fbe45STzanio // make bp3 [MFEM_DIR=</path/to/mfem>] [CEED_DIR=</path/to/libceed>]
23182fbe45STzanio //
24182fbe45STzanio // Sample runs:
25182fbe45STzanio //
2666087c08SValeria Barra // ./bp3
2766087c08SValeria Barra // ./bp3 -ceed /cpu/self
2828688798Sjeremylt // ./bp3 -ceed /gpu/cuda
2966087c08SValeria Barra // ./bp3 -m ../../../mfem/data/fichera.mesh -o 4
3066087c08SValeria Barra // ./bp3 -m ../../../mfem/data/square-disc-nurbs.mesh -o 6
3166087c08SValeria Barra // ./bp3 -m ../../../mfem/data/inline-segment.mesh -o 8
32182fbe45STzanio
335d6bafb2Sjeremylt /// @file
345d6bafb2Sjeremylt /// MFEM diffusion operator based on libCEED
355d6bafb2Sjeremylt
36c0c38e35SVeselin Dobrev #include "bp3.hpp"
37182fbe45STzanio
382b730f8bSJeremy L Thompson #include <ceed.h>
392b730f8bSJeremy L Thompson
402b730f8bSJeremy L Thompson #include <mfem.hpp>
412b730f8bSJeremy L Thompson
42182fbe45STzanio /// Exact solution
solution(const mfem::Vector & pt)43182fbe45STzanio double solution(const mfem::Vector &pt) {
44182fbe45STzanio static const double x[3] = {-0.32, 0.15, 0.24};
45182fbe45STzanio static const double k[3] = {1.21, 1.45, 1.37};
46182fbe45STzanio double val = sin(M_PI * (x[0] + k[0] * pt(0)));
472b730f8bSJeremy L Thompson for (int d = 1; d < pt.Size(); d++) val *= sin(M_PI * (x[d] + k[d] * pt(d)));
48182fbe45STzanio return val;
49182fbe45STzanio }
50182fbe45STzanio
51182fbe45STzanio /// Right-hand side
rhs(const mfem::Vector & pt)52182fbe45STzanio double rhs(const mfem::Vector &pt) {
53182fbe45STzanio static const double x[3] = {-0.32, 0.15, 0.24};
54182fbe45STzanio static const double k[3] = {1.21, 1.45, 1.37};
55182fbe45STzanio double f[3], l[3], val, lap;
56182fbe45STzanio f[0] = sin(M_PI * (x[0] + k[0] * pt(0)));
57182fbe45STzanio l[0] = M_PI * M_PI * k[0] * k[0] * f[0];
58182fbe45STzanio val = f[0];
59182fbe45STzanio lap = l[0];
60182fbe45STzanio for (int d = 1; d < pt.Size(); d++) {
61182fbe45STzanio f[d] = sin(M_PI * (x[d] + k[d] * pt(d)));
62182fbe45STzanio l[d] = M_PI * M_PI * k[d] * k[d] * f[d];
63182fbe45STzanio lap = lap * f[d] + val * l[d];
64182fbe45STzanio val = val * f[d];
65182fbe45STzanio }
66182fbe45STzanio return lap;
67182fbe45STzanio }
68182fbe45STzanio
69f063656dSJed Brown //TESTARGS -ceed {ceed_resource} -t -no-vis --size 2000
main(int argc,char * argv[])70182fbe45STzanio int main(int argc, char *argv[]) {
71182fbe45STzanio // 1. Parse command-line options.
72182fbe45STzanio const char *ceed_spec = "/cpu/self";
73c0c38e35SVeselin Dobrev #ifndef MFEM_DIR
74182fbe45STzanio const char *mesh_file = "../../../mfem/data/star.mesh";
75c0c38e35SVeselin Dobrev #else
76c0c38e35SVeselin Dobrev const char *mesh_file = MFEM_DIR "/data/star.mesh";
77c0c38e35SVeselin Dobrev #endif
78182fbe45STzanio int order = 2;
79182fbe45STzanio bool visualization = true;
80dc00e230Sjeremylt bool test = false;
81e2b2c771Svaleria double max_nnodes = 50000;
82182fbe45STzanio
83182fbe45STzanio mfem::OptionsParser args(argc, argv);
84182fbe45STzanio args.AddOption(&ceed_spec, "-c", "-ceed", "Ceed specification.");
85182fbe45STzanio args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
862b730f8bSJeremy L Thompson args.AddOption(&order, "-o", "--order", "Finite element order (polynomial degree).");
87e2b2c771Svaleria args.AddOption(&max_nnodes, "-s", "--size", "Maximum size (number of DoFs)");
882b730f8bSJeremy L Thompson args.AddOption(&visualization, "-vis", "--visualization", "-no-vis", "--no-visualization", "Enable or disable GLVis visualization.");
892b730f8bSJeremy L Thompson args.AddOption(&test, "-t", "--test", "-no-test", "--no-test", "Enable or disable test mode.");
90182fbe45STzanio args.Parse();
91182fbe45STzanio if (!args.Good()) {
92182fbe45STzanio args.PrintUsage(std::cout);
93182fbe45STzanio return 1;
94182fbe45STzanio }
95dc00e230Sjeremylt if (!test) {
96182fbe45STzanio args.PrintOptions(std::cout);
97dc00e230Sjeremylt }
98182fbe45STzanio
99182fbe45STzanio // 2. Initialize a Ceed device object using the given Ceed specification.
100182fbe45STzanio Ceed ceed;
101182fbe45STzanio CeedInit(ceed_spec, &ceed);
102182fbe45STzanio
103182fbe45STzanio // 3. Read the mesh from the given mesh file.
104182fbe45STzanio mfem::Mesh *mesh = new mfem::Mesh(mesh_file, 1, 1);
105182fbe45STzanio int dim = mesh->Dimension();
106182fbe45STzanio
107ea61e9acSJeremy L Thompson // 4. Refine the mesh to increase the resolution.
108ea61e9acSJeremy L Thompson // In this example we do 'ref_levels' of uniform refinement.
109ea61e9acSJeremy L Thompson // We choose 'ref_levels' to be the largest number that gives a final system with no more than 50,000 unknowns, approximately.
110182fbe45STzanio {
1112b730f8bSJeremy L Thompson int ref_levels = (int)floor((log(max_nnodes / mesh->GetNE()) - dim * log(order)) / log(2.) / dim);
112182fbe45STzanio for (int l = 0; l < ref_levels; l++) {
113182fbe45STzanio mesh->UniformRefinement();
114182fbe45STzanio }
115182fbe45STzanio }
116182fbe45STzanio if (mesh->GetNodalFESpace() == NULL) {
117182fbe45STzanio mesh->SetCurvature(1, false, -1, mfem::Ordering::byNODES);
118182fbe45STzanio }
119182fbe45STzanio if (mesh->NURBSext) {
120182fbe45STzanio mesh->SetCurvature(order, false, -1, mfem::Ordering::byNODES);
121182fbe45STzanio }
122182fbe45STzanio
123ea61e9acSJeremy L Thompson // 5. Define a finite element space on the mesh.
124ea61e9acSJeremy L Thompson // Here we use continuous Lagrange finite elements of the specified order.
125182fbe45STzanio MFEM_VERIFY(order > 0, "invalid order");
126182fbe45STzanio mfem::FiniteElementCollection *fec = new mfem::H1_FECollection(order, dim);
127182fbe45STzanio mfem::FiniteElementSpace *fespace = new mfem::FiniteElementSpace(mesh, fec);
128dc00e230Sjeremylt if (!test) {
1292b730f8bSJeremy L Thompson std::cout << "Number of finite element unknowns: " << fespace->GetTrueVSize() << std::endl;
130dc00e230Sjeremylt }
131182fbe45STzanio
132182fbe45STzanio mfem::FunctionCoefficient sol_coeff(solution);
133182fbe45STzanio mfem::Array<int> ess_tdof_list;
134182fbe45STzanio mfem::GridFunction sol(fespace);
135182fbe45STzanio if (mesh->bdr_attributes.Size()) {
136182fbe45STzanio mfem::Array<int> ess_bdr(mesh->bdr_attributes.Max());
137182fbe45STzanio ess_bdr = 1;
138182fbe45STzanio fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
139182fbe45STzanio sol.ProjectBdrCoefficient(sol_coeff, ess_bdr);
140182fbe45STzanio }
141182fbe45STzanio
142ea61e9acSJeremy L Thompson // 6. Construct a rhs vector using the linear form f(v) = (rhs, v), where 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
148ea61e9acSJeremy L Thompson // 7. Construct a CeedDiffusionOperator utilizing the 'ceed' device and using the 'fespace' object to extract data needed by the Ceed objects.
149182fbe45STzanio CeedDiffusionOperator diff(ceed, fespace);
150182fbe45STzanio
151182fbe45STzanio mfem::Operator *D;
152182fbe45STzanio mfem::Vector X, B;
153182fbe45STzanio diff.FormLinearSystem(ess_tdof_list, sol, b, D, X, B);
154182fbe45STzanio
155182fbe45STzanio // 8. Solve the discrete system using the conjugate gradients (CG) method.
156182fbe45STzanio mfem::CGSolver cg;
157182fbe45STzanio cg.SetRelTol(1e-6);
158182fbe45STzanio cg.SetMaxIter(1000);
159dc00e230Sjeremylt if (test) {
160dc00e230Sjeremylt cg.SetPrintLevel(0);
161dc00e230Sjeremylt } else {
162182fbe45STzanio cg.SetPrintLevel(3);
163dc00e230Sjeremylt }
164182fbe45STzanio cg.SetOperator(*D);
165182fbe45STzanio
166182fbe45STzanio cg.Mult(B, X);
167182fbe45STzanio
168182fbe45STzanio // 9. Compute and print the L2 norm of the error.
1699647a07eSDavid Medina double err_l2 = sol.ComputeL2Error(sol_coeff);
170dc00e230Sjeremylt if (!test) {
1712b730f8bSJeremy L Thompson std::cout << "L2 projection error: " << err_l2 << std::endl;
172dc00e230Sjeremylt } else {
173f063656dSJed Brown if (fabs(sol.ComputeL2Error(sol_coeff)) > 2e-3) {
1749647a07eSDavid Medina std::cout << "Error too large: " << err_l2 << std::endl;
175dc00e230Sjeremylt }
176dc00e230Sjeremylt }
177182fbe45STzanio
178ea61e9acSJeremy L Thompson // 10. Open a socket connection to GLVis and send the mesh and solution for visualization.
179182fbe45STzanio if (visualization) {
180182fbe45STzanio char vishost[] = "localhost";
181182fbe45STzanio int visport = 19916;
182182fbe45STzanio mfem::socketstream sol_sock(vishost, visport);
183182fbe45STzanio sol_sock.precision(8);
184182fbe45STzanio sol_sock << "solution\n" << *mesh << sol << std::flush;
185182fbe45STzanio }
186182fbe45STzanio
187182fbe45STzanio // 11. Free memory and exit.
188182fbe45STzanio delete fespace;
189182fbe45STzanio delete fec;
190182fbe45STzanio delete mesh;
191637f263aSJeremy L Thompson delete D;
192182fbe45STzanio CeedDestroy(&ceed);
193182fbe45STzanio return 0;
194182fbe45STzanio }
195