xref: /libCEED/examples/petsc/area.c (revision ccf0fe6fd2ec6c692fac5dfe411b8e3ec625937d)
1cb32e2e7SValeria Barra // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2cb32e2e7SValeria Barra // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3cb32e2e7SValeria Barra // reserved. See files LICENSE and NOTICE for details.
4cb32e2e7SValeria Barra //
5cb32e2e7SValeria Barra // This file is part of CEED, a collection of benchmarks, miniapps, software
6cb32e2e7SValeria Barra // libraries and APIs for efficient high-order finite element and spectral
7cb32e2e7SValeria Barra // element discretizations for exascale applications. For more information and
8cb32e2e7SValeria Barra // source code availability see http://github.com/ceed.
9cb32e2e7SValeria Barra //
10cb32e2e7SValeria Barra // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11cb32e2e7SValeria Barra // a collaborative effort of two U.S. Department of Energy organizations (Office
12cb32e2e7SValeria Barra // of Science and the National Nuclear Security Administration) responsible for
13cb32e2e7SValeria Barra // the planning and preparation of a capable exascale ecosystem, including
14cb32e2e7SValeria Barra // software, applications, hardware, advanced system engineering and early
15cb32e2e7SValeria Barra // testbed platforms, in support of the nation's exascale computing imperative.
16cb32e2e7SValeria Barra 
17cb32e2e7SValeria Barra //                        libCEED + PETSc Example: Surface Area
18cb32e2e7SValeria Barra //
19cb32e2e7SValeria Barra // This example demonstrates a simple usage of libCEED with PETSc to calculate
20cb32e2e7SValeria Barra // the surface area of a simple closed surface, such as the one of a cube
21cb32e2e7SValeria Barra // via the mass operator.
22cb32e2e7SValeria Barra //
23cb32e2e7SValeria Barra // The code uses higher level communication protocols in DMPlex.
24cb32e2e7SValeria Barra //
25cb32e2e7SValeria Barra // Build with:
26cb32e2e7SValeria Barra //
27cb32e2e7SValeria Barra //     make area [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
28cb32e2e7SValeria Barra //
29cb32e2e7SValeria Barra // Sample runs:
30cb32e2e7SValeria Barra //   Sequential:
31cb32e2e7SValeria Barra //
32cb32e2e7SValeria Barra //     area -petscspace_degree 3
33cb32e2e7SValeria Barra //
34cb32e2e7SValeria Barra //   In parallel:
35cb32e2e7SValeria Barra //
36cb32e2e7SValeria Barra //     mpiexec -n 4 area -petscspace_degree 3
37cb32e2e7SValeria Barra //
38cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -petscspace_degree 3
39cb32e2e7SValeria Barra 
40cb32e2e7SValeria Barra /// @file
41cb32e2e7SValeria Barra /// libCEED example using the mass operator to compute surface area using PETSc with DMPlex
42*ccf0fe6fSjeremylt static const char help[] =
43*ccf0fe6fSjeremylt   "Compute surface area of a cube using DMPlex in PETSc\n";
44cb32e2e7SValeria Barra 
45cb32e2e7SValeria Barra #include <string.h>
46cb32e2e7SValeria Barra #include <petscdmplex.h>
47cb32e2e7SValeria Barra #include <ceed.h>
48cb32e2e7SValeria Barra #include "qfunctions/area/area.h"
49cb32e2e7SValeria Barra 
50cb32e2e7SValeria Barra // Auxiliary function to define CEED restrictions from DMPlex data
51cb32e2e7SValeria Barra static int CreateRestrictionPlex(Ceed ceed, CeedInt P, CeedInt ncomp,
52cb32e2e7SValeria Barra                                  CeedElemRestriction *Erestrict, DM dm) {
53cb32e2e7SValeria Barra   PetscInt ierr;
54cb32e2e7SValeria Barra   PetscInt c, cStart, cEnd, nelem, nnodes, *erestrict, eoffset;
55cb32e2e7SValeria Barra   PetscSection section;
56cb32e2e7SValeria Barra   Vec Uloc;
57cb32e2e7SValeria Barra 
58cb32e2e7SValeria Barra   PetscFunctionBegin;
59cb32e2e7SValeria Barra 
60cb32e2e7SValeria Barra   // Get Nelem
61cb32e2e7SValeria Barra   ierr = DMGetSection(dm, &section); CHKERRQ(ierr);
62cb32e2e7SValeria Barra   ierr = DMPlexGetHeightStratum(dm, 0, &cStart,& cEnd); CHKERRQ(ierr);
63cb32e2e7SValeria Barra   nelem = cEnd - cStart;
64cb32e2e7SValeria Barra 
65cb32e2e7SValeria Barra   // Get indices
66cb32e2e7SValeria Barra   ierr = PetscMalloc1(nelem*P*P, &erestrict); CHKERRQ(ierr);
67cb32e2e7SValeria Barra   for (c=cStart, eoffset = 0; c<cEnd; c++) {
68cb32e2e7SValeria Barra     PetscInt numindices, *indices, i;
69cb32e2e7SValeria Barra     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
70cb32e2e7SValeria Barra                                    &indices, NULL); CHKERRQ(ierr);
71cb32e2e7SValeria Barra     for (i=0; i<numindices; i+=ncomp) {
72cb32e2e7SValeria Barra       for (PetscInt j=0; j<ncomp; j++) {
73cb32e2e7SValeria Barra         if (indices[i+j] != indices[i] + (PetscInt)(copysign(j, indices[i])))
74cb32e2e7SValeria Barra           SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
75cb32e2e7SValeria Barra                    "Cell %D closure indices not interlaced", c);
76cb32e2e7SValeria Barra       }
77cb32e2e7SValeria Barra       // NO BC on closed surfaces
78cb32e2e7SValeria Barra       PetscInt loc = indices[i];
79cb32e2e7SValeria Barra       erestrict[eoffset++] = loc/ncomp;
80cb32e2e7SValeria Barra     }
81cb32e2e7SValeria Barra     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
82cb32e2e7SValeria Barra                                        &indices, NULL); CHKERRQ(ierr);
83cb32e2e7SValeria Barra   }
84cb32e2e7SValeria Barra 
85cb32e2e7SValeria Barra   // Setup CEED restriction
86cb32e2e7SValeria Barra   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
87cb32e2e7SValeria Barra   ierr = VecGetLocalSize(Uloc, &nnodes); CHKERRQ(ierr);
88cb32e2e7SValeria Barra 
89cb32e2e7SValeria Barra   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
90cb32e2e7SValeria Barra   CeedElemRestrictionCreate(ceed, nelem, P*P, nnodes/ncomp, ncomp,
91cb32e2e7SValeria Barra                             CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
92cb32e2e7SValeria Barra                             Erestrict);
93cb32e2e7SValeria Barra   ierr = PetscFree(erestrict); CHKERRQ(ierr);
94cb32e2e7SValeria Barra 
95cb32e2e7SValeria Barra   PetscFunctionReturn(0);
96cb32e2e7SValeria Barra }
97cb32e2e7SValeria Barra 
98cb32e2e7SValeria Barra int main(int argc, char **argv) {
99cb32e2e7SValeria Barra   PetscInt ierr;
100cb32e2e7SValeria Barra   MPI_Comm comm;
101cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
102cb0b5415Sjeremylt        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
103cb32e2e7SValeria Barra   PetscInt lsize, gsize, xlsize,
104cb32e2e7SValeria Barra            qextra  = 1, // default number of extra quadrature points
105cb32e2e7SValeria Barra            ncompx  = 3, // number of components of 3D physical coordinates
106cb32e2e7SValeria Barra            ncompu  = 1, // dimension of field to which apply mass operator
107cb32e2e7SValeria Barra            topodim = 2, // topological dimension of manifold
108cb32e2e7SValeria Barra            degree  = 3; // default degree for finite element bases
109cb32e2e7SValeria Barra   PetscBool read_mesh = PETSC_FALSE,
110cb32e2e7SValeria Barra             test_mode = PETSC_FALSE;
111cb32e2e7SValeria Barra   PetscSpace sp;
112cb32e2e7SValeria Barra   PetscFE fe;
113cb32e2e7SValeria Barra   Vec X, Xloc, V, Vloc;
114cb32e2e7SValeria Barra   DM  dm, dmcoord;
115cb32e2e7SValeria Barra   Ceed ceed;
116cb32e2e7SValeria Barra   CeedInt P, Q;
117cb32e2e7SValeria Barra 
118cb32e2e7SValeria Barra   ierr = PetscInitialize(&argc, &argv, NULL, help);
119cb32e2e7SValeria Barra   if (ierr) return ierr;
120cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
121cb32e2e7SValeria Barra 
122cb32e2e7SValeria Barra   // Read CL options
123*ccf0fe6fSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc",
124*ccf0fe6fSjeremylt                            NULL);
125cb32e2e7SValeria Barra   CHKERRQ(ierr);
126cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
127cb32e2e7SValeria Barra                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
128cb32e2e7SValeria Barra   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
129cb32e2e7SValeria Barra                             NULL, ceedresource, ceedresource,
130cb32e2e7SValeria Barra                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
131cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-test",
132cb32e2e7SValeria Barra                           "Testing mode (do not print unless error is large)",
133cb32e2e7SValeria Barra                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
134cb32e2e7SValeria Barra   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
135cb32e2e7SValeria Barra                             filename, filename, sizeof(filename), &read_mesh);
136cb32e2e7SValeria Barra   CHKERRQ(ierr);
137cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
138cb32e2e7SValeria Barra 
139cb32e2e7SValeria Barra   // Setup DM
140cb32e2e7SValeria Barra   PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube
141cb32e2e7SValeria Barra   if (read_mesh) {
142cb32e2e7SValeria Barra     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
143cb32e2e7SValeria Barra     CHKERRQ(ierr);
144cb32e2e7SValeria Barra   } else {
145cb32e2e7SValeria Barra     // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box.
146cb32e2e7SValeria Barra     PetscBool simplex = PETSC_FALSE;
147cb32e2e7SValeria Barra     ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm);
148cb32e2e7SValeria Barra     CHKERRQ(ierr);
149cb32e2e7SValeria Barra     // Set the object name
150cb32e2e7SValeria Barra     ierr = PetscObjectSetName((PetscObject) dm, "Cube"); CHKERRQ(ierr);
151cb32e2e7SValeria Barra     // Distribute mesh over processes
152cb32e2e7SValeria Barra     {
153cb32e2e7SValeria Barra       DM dmDist = NULL;
154cb32e2e7SValeria Barra       PetscPartitioner part;
155cb32e2e7SValeria Barra 
156cb32e2e7SValeria Barra       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
157cb32e2e7SValeria Barra       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
158cb32e2e7SValeria Barra       ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
159cb32e2e7SValeria Barra       if (dmDist) {
160cb32e2e7SValeria Barra         ierr = DMDestroy(&dm); CHKERRQ(ierr);
161cb32e2e7SValeria Barra         dm  = dmDist;
162cb32e2e7SValeria Barra       }
163cb32e2e7SValeria Barra     }
164cb32e2e7SValeria Barra     // View DMPlex via runtime option
165cb32e2e7SValeria Barra     ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
166cb32e2e7SValeria Barra   }
167cb32e2e7SValeria Barra 
168cb32e2e7SValeria Barra   // Create FE
169cb32e2e7SValeria Barra   ierr = PetscFECreateDefault(PETSC_COMM_SELF, topodim, ncompu, PETSC_FALSE, NULL,
170cb32e2e7SValeria Barra                               PETSC_DETERMINE, &fe);
171cb32e2e7SValeria Barra   CHKERRQ(ierr);
172cb32e2e7SValeria Barra   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
173cb32e2e7SValeria Barra   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
174cb32e2e7SValeria Barra   ierr = DMCreateDS(dm); CHKERRQ(ierr);
175cb32e2e7SValeria Barra   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
176cb32e2e7SValeria Barra   CHKERRQ(ierr);
177cb32e2e7SValeria Barra 
178cb32e2e7SValeria Barra   // Get basis space degree
179cb32e2e7SValeria Barra   ierr = PetscFEGetBasisSpace(fe, &sp); CHKERRQ(ierr);
180cb32e2e7SValeria Barra   ierr = PetscSpaceGetDegree(sp, &degree, NULL); CHKERRQ(ierr);
181cb32e2e7SValeria Barra   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
182cb32e2e7SValeria Barra   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
183cb32e2e7SValeria Barra                              "-petscspace_degree %D must be at least 1", degree);
184cb32e2e7SValeria Barra 
185cb32e2e7SValeria Barra   // Create vectors
186cb32e2e7SValeria Barra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
187cb32e2e7SValeria Barra   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
188cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
189cb32e2e7SValeria Barra   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
190cb32e2e7SValeria Barra   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
191cb32e2e7SValeria Barra   ierr = VecDuplicate(X, &V); CHKERRQ(ierr);
192cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &Vloc); CHKERRQ(ierr);
193cb32e2e7SValeria Barra 
194cb32e2e7SValeria Barra   // Set up libCEED
195cb32e2e7SValeria Barra   CeedInit(ceedresource, &ceed);
196cb32e2e7SValeria Barra 
197cb32e2e7SValeria Barra   // Print summary
198cb32e2e7SValeria Barra   P = degree + 1;
199cb32e2e7SValeria Barra   Q = P + qextra;
200cb32e2e7SValeria Barra   const char *usedresource;
201cb32e2e7SValeria Barra   CeedGetResource(ceed, &usedresource);
202cb32e2e7SValeria Barra   if (!test_mode) {
203cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
204cb32e2e7SValeria Barra                        "\n-- libCEED + PETSc Surface Area problem --\n"
205cb32e2e7SValeria Barra                        "  libCEED:\n"
206cb32e2e7SValeria Barra                        "    libCEED Backend                    : %s\n"
207cb32e2e7SValeria Barra                        "  Mesh:\n"
208cb32e2e7SValeria Barra                        "    Number of 1D Basis Nodes (p)       : %d\n"
209cb32e2e7SValeria Barra                        "    Number of 1D Quadrature Points (q) : %d\n"
210cb32e2e7SValeria Barra                        "    Global nodes                       : %D\n",
211cb32e2e7SValeria Barra                        usedresource, P, Q,  gsize/ncompu);
212cb32e2e7SValeria Barra     CHKERRQ(ierr);
213cb32e2e7SValeria Barra   }
214cb32e2e7SValeria Barra 
215cb32e2e7SValeria Barra   // Setup libCEED's objects
216cb32e2e7SValeria Barra   // Create CEED operators and API objects they need
217cb32e2e7SValeria Barra   CeedOperator op_setupgeo, op_apply;
218cb32e2e7SValeria Barra   CeedQFunction qf_setupgeo, qf_apply;
219cb32e2e7SValeria Barra   CeedBasis basisx, basisu;
220cb32e2e7SValeria Barra   CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi,
221cb32e2e7SValeria Barra                       Erestrictqdi;
222cb32e2e7SValeria Barra 
223cb32e2e7SValeria Barra   // Create bases
224cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q,
225cb32e2e7SValeria Barra                                   CEED_GAUSS, &basisu);
226cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q,
227cb32e2e7SValeria Barra                                   CEED_GAUSS, &basisx);
228cb32e2e7SValeria Barra 
229cb32e2e7SValeria Barra   // CEED restrictions
230cb32e2e7SValeria Barra   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
231cb32e2e7SValeria Barra   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
232cb32e2e7SValeria Barra   CHKERRQ(ierr);
233cb32e2e7SValeria Barra 
234cb32e2e7SValeria Barra   CreateRestrictionPlex(ceed, 2, ncompx, &Erestrictx, dmcoord); CHKERRQ(ierr);
235cb32e2e7SValeria Barra   CreateRestrictionPlex(ceed, P, ncompu, &Erestrictu, dm); CHKERRQ(ierr);
236cb32e2e7SValeria Barra 
237cb32e2e7SValeria Barra   CeedInt cStart, cEnd;
238cb32e2e7SValeria Barra   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
239cb32e2e7SValeria Barra   const CeedInt nelem = cEnd - cStart;
240cb32e2e7SValeria Barra 
241cb32e2e7SValeria Barra   // CEED identity restrictions
242cb32e2e7SValeria Barra   const CeedInt qdatasize = 1;
243cb32e2e7SValeria Barra   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q, nelem*Q*Q,
244cb32e2e7SValeria Barra                                     qdatasize, &Erestrictqdi);
245cb32e2e7SValeria Barra   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q, nelem*Q*Q, 1,
246cb32e2e7SValeria Barra                                     &Erestrictxi);
247cb32e2e7SValeria Barra 
248cb32e2e7SValeria Barra   // Element coordinates
249cb32e2e7SValeria Barra   Vec coords;
250cb32e2e7SValeria Barra   const PetscScalar *coordArray;
251cb32e2e7SValeria Barra   PetscSection section;
252cb32e2e7SValeria Barra   ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr);
253cb32e2e7SValeria Barra   ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr);
254cb32e2e7SValeria Barra   ierr = DMGetSection(dmcoord, &section); CHKERRQ(ierr);
255cb32e2e7SValeria Barra 
256cb32e2e7SValeria Barra   CeedVector xcoord;
257cb32e2e7SValeria Barra   CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL);
258cb32e2e7SValeria Barra   CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES,
259cb32e2e7SValeria Barra                      (PetscScalar *)coordArray);
260cb32e2e7SValeria Barra   ierr = VecRestoreArrayRead(coords, &coordArray);
261cb32e2e7SValeria Barra 
262cb32e2e7SValeria Barra   // Create the vectors that will be needed in setup and apply
263cb32e2e7SValeria Barra   CeedVector uceed, vceed, qdata;
264cb32e2e7SValeria Barra   CeedInt nqpts;
265cb32e2e7SValeria Barra   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
266cb32e2e7SValeria Barra   CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata);
267cb32e2e7SValeria Barra   CeedVectorCreate(ceed, xlsize, &uceed);
268cb32e2e7SValeria Barra   CeedVectorCreate(ceed, xlsize, &vceed);
269cb32e2e7SValeria Barra 
270cb32e2e7SValeria Barra   /* Create the Q-function that builds the operator for the geomteric factors
271cb32e2e7SValeria Barra      (i.e., the quadrature data) */
272cb32e2e7SValeria Barra   CeedQFunctionCreateInterior(ceed, 1, SetupMassGeo,
273cb32e2e7SValeria Barra                               SetupMassGeo_loc, &qf_setupgeo);
274cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*topodim, CEED_EVAL_GRAD);
275cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT);
276cb32e2e7SValeria Barra   CeedQFunctionAddOutput(qf_setupgeo, "qdata", qdatasize, CEED_EVAL_NONE);
277cb32e2e7SValeria Barra 
278cb32e2e7SValeria Barra   // Set up the mass operator
279cb32e2e7SValeria Barra   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_apply);
280cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_apply, "u", ncompu, CEED_EVAL_INTERP);
281cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_apply, "qdata", qdatasize, CEED_EVAL_NONE);
282cb32e2e7SValeria Barra   CeedQFunctionAddOutput(qf_apply, "v", ncompu, CEED_EVAL_INTERP);
283cb32e2e7SValeria Barra 
284cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the operator
285cb32e2e7SValeria Barra   CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo);
286cb32e2e7SValeria Barra   CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, CEED_TRANSPOSE,
287cb32e2e7SValeria Barra                        basisx, CEED_VECTOR_ACTIVE);
288cb32e2e7SValeria Barra   CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, CEED_NOTRANSPOSE,
289cb32e2e7SValeria Barra                        basisx, CEED_VECTOR_NONE);
290cb32e2e7SValeria Barra   CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, CEED_NOTRANSPOSE,
291cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
292cb32e2e7SValeria Barra 
293cb32e2e7SValeria Barra   // Create the mass operator
294cb32e2e7SValeria Barra   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
295cb32e2e7SValeria Barra   CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE,
296cb32e2e7SValeria Barra                        basisu, CEED_VECTOR_ACTIVE);
297cb32e2e7SValeria Barra   CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_NOTRANSPOSE,
298cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, qdata);
299cb32e2e7SValeria Barra   CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE,
300cb32e2e7SValeria Barra                        basisu, CEED_VECTOR_ACTIVE);
301cb32e2e7SValeria Barra 
302cb32e2e7SValeria Barra   // Compute the quadrature data for the mass operator
303cb32e2e7SValeria Barra   CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
304cb32e2e7SValeria Barra 
305cb32e2e7SValeria Barra   PetscScalar *v;
306cb32e2e7SValeria Barra   ierr = VecZeroEntries(Vloc); CHKERRQ(ierr);
307cb32e2e7SValeria Barra   ierr = VecGetArray(Vloc, &v);
308cb32e2e7SValeria Barra   CeedVectorSetArray(vceed, CEED_MEM_HOST, CEED_USE_POINTER, v);
309cb32e2e7SValeria Barra 
310cb32e2e7SValeria Barra   // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1
311cb32e2e7SValeria Barra   if (!test_mode) {
312cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
313cb32e2e7SValeria Barra                        "Computing the mesh volume using the formula: vol = 1^T M 1\n");
314cb32e2e7SValeria Barra     CHKERRQ(ierr);
315cb32e2e7SValeria Barra   }
316cb32e2e7SValeria Barra 
317cb32e2e7SValeria Barra   // Initialize u and v with ones
318cb32e2e7SValeria Barra   CeedVectorSetValue(uceed, 1.0);
319cb32e2e7SValeria Barra   CeedVectorSetValue(vceed, 1.0);
320cb32e2e7SValeria Barra 
321cb32e2e7SValeria Barra   // Apply the mass operator: 'u' -> 'v'
322cb32e2e7SValeria Barra   CeedOperatorApply(op_apply, uceed, vceed, CEED_REQUEST_IMMEDIATE);
323cb32e2e7SValeria Barra   CeedVectorSyncArray(vceed, CEED_MEM_HOST);
324cb32e2e7SValeria Barra 
325cb32e2e7SValeria Barra   // Gather output vector
326cb32e2e7SValeria Barra   ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr);
327cb32e2e7SValeria Barra   ierr = VecZeroEntries(V); CHKERRQ(ierr);
328cb32e2e7SValeria Barra   ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
329cb32e2e7SValeria Barra   ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
330cb32e2e7SValeria Barra 
331cb32e2e7SValeria Barra   // Compute and print the sum of the entries of 'v' giving the mesh surface area
332cb32e2e7SValeria Barra   PetscScalar area;
333cb32e2e7SValeria Barra   ierr = VecSum(V, &area); CHKERRQ(ierr);
334cb32e2e7SValeria Barra 
335cb32e2e7SValeria Barra   // Compute the exact surface area and print the result
336cb32e2e7SValeria Barra   CeedScalar exact_surfarea = 6 * (2*l) * (2*l);
337cb32e2e7SValeria Barra   if (!test_mode) {
338cb32e2e7SValeria Barra     ierr = PetscPrintf(comm, "Exact mesh surface area    : % .14g\n",
339cb32e2e7SValeria Barra                        exact_surfarea); CHKERRQ(ierr);
340cb32e2e7SValeria Barra     ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area);
341cb32e2e7SValeria Barra     CHKERRQ(ierr);
342cb32e2e7SValeria Barra     ierr = PetscPrintf(comm, "Area error                 : % .14g\n",
343cb32e2e7SValeria Barra                        fabs(area - exact_surfarea)); CHKERRQ(ierr);
344cb32e2e7SValeria Barra   }
345cb32e2e7SValeria Barra 
346cb32e2e7SValeria Barra   // PETSc cleanup
347cb32e2e7SValeria Barra   ierr = DMDestroy(&dm); CHKERRQ(ierr);
348cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
349cb32e2e7SValeria Barra   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
350cb32e2e7SValeria Barra   ierr = VecDestroy(&V); CHKERRQ(ierr);
351cb32e2e7SValeria Barra   ierr = VecDestroy(&Vloc); CHKERRQ(ierr);
352cb32e2e7SValeria Barra 
353cb32e2e7SValeria Barra   // libCEED cleanup
354cb32e2e7SValeria Barra   CeedQFunctionDestroy(&qf_setupgeo);
355cb32e2e7SValeria Barra   CeedOperatorDestroy(&op_setupgeo);
356cb32e2e7SValeria Barra   CeedVectorDestroy(&xcoord);
357cb32e2e7SValeria Barra   CeedVectorDestroy(&uceed);
358cb32e2e7SValeria Barra   CeedVectorDestroy(&vceed);
359cb32e2e7SValeria Barra   CeedVectorDestroy(&qdata);
360cb32e2e7SValeria Barra   CeedBasisDestroy(&basisx);
361cb32e2e7SValeria Barra   CeedBasisDestroy(&basisu);
362cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictu);
363cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictx);
364cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictxi);
365cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictqdi);
366cb32e2e7SValeria Barra   CeedQFunctionDestroy(&qf_apply);
367cb32e2e7SValeria Barra   CeedOperatorDestroy(&op_apply);
368cb32e2e7SValeria Barra   CeedDestroy(&ceed);
369cb32e2e7SValeria Barra   return PetscFinalize();
370cb32e2e7SValeria Barra }
371