xref: /libCEED/examples/petsc/area.c (revision 15910d16b955338d1102d4e730fc58bca8f202b9)
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
2032d2ee49SValeria Barra // the surface area of a simple closed surface, such as the one of a cube or a
2132d2ee49SValeria Barra // tensor-product discrete sphere 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 //
3232d2ee49SValeria Barra //     ./area -problem cube -petscspace_degree 3 -dm_refine 2
3332d2ee49SValeria Barra //
3432d2ee49SValeria Barra //     ./area -problem sphere -petscspace_degree 3 -dm_refine 2
35cb32e2e7SValeria Barra //
36cb32e2e7SValeria Barra //   In parallel:
37cb32e2e7SValeria Barra //
3832d2ee49SValeria Barra //     mpiexec -n 4 ./area -probelm cube -petscspace_degree 3 -dm_refine 2
3932d2ee49SValeria Barra //
4032d2ee49SValeria Barra //     mpiexec -n 4 ./area -problem sphere -petscspace_degree 3 -dm_refine 2
4132d2ee49SValeria Barra //
4232d2ee49SValeria Barra //   The above example runs use 2 levels of refinement for the mesh.
4332d2ee49SValeria Barra //   Use -dm_refine k, for k levels of uniform refinement.
44cb32e2e7SValeria Barra //
45cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -petscspace_degree 3
46cb32e2e7SValeria Barra 
47cb32e2e7SValeria Barra /// @file
4832d2ee49SValeria Barra /// libCEED example using the mass operator to compute a cube or a cubed-sphere surface area using PETSc with DMPlex
49ccf0fe6fSjeremylt static const char help[] =
5032d2ee49SValeria Barra   "Compute surface area of a cube or a cubed-sphere using DMPlex in PETSc\n";
51cb32e2e7SValeria Barra 
52cb32e2e7SValeria Barra #include <string.h>
53cb32e2e7SValeria Barra #include <petscdmplex.h>
54cb32e2e7SValeria Barra #include <ceed.h>
5532d2ee49SValeria Barra #include "qfunctions/area/areacube.h"
5632d2ee49SValeria Barra #include "qfunctions/area/areasphere.h"
5732d2ee49SValeria Barra 
5832d2ee49SValeria Barra #ifndef M_PI
5932d2ee49SValeria Barra #  define M_PI    3.14159265358979323846
6032d2ee49SValeria Barra #endif
61cb32e2e7SValeria Barra 
62cb32e2e7SValeria Barra // Auxiliary function to define CEED restrictions from DMPlex data
6361dbc9d2Sjeremylt static int CreateRestrictionPlex(Ceed ceed, CeedInterlaceMode imode, CeedInt P,
64a8d32208Sjeremylt                                  CeedInt ncomp, CeedElemRestriction *Erestrict,
65a8d32208Sjeremylt                                  DM dm) {
66cb32e2e7SValeria Barra   PetscInt ierr;
67cb32e2e7SValeria Barra   PetscInt c, cStart, cEnd, nelem, nnodes, *erestrict, eoffset;
68cb32e2e7SValeria Barra   PetscSection section;
69cb32e2e7SValeria Barra   Vec Uloc;
70cb32e2e7SValeria Barra 
71cb32e2e7SValeria Barra   PetscFunctionBegin;
72cb32e2e7SValeria Barra 
73cb32e2e7SValeria Barra   // Get Nelem
74cb32e2e7SValeria Barra   ierr = DMGetSection(dm, &section); CHKERRQ(ierr);
75cb32e2e7SValeria Barra   ierr = DMPlexGetHeightStratum(dm, 0, &cStart,& cEnd); CHKERRQ(ierr);
76cb32e2e7SValeria Barra   nelem = cEnd - cStart;
77cb32e2e7SValeria Barra 
78cb32e2e7SValeria Barra   // Get indices
79cb32e2e7SValeria Barra   ierr = PetscMalloc1(nelem*P*P, &erestrict); CHKERRQ(ierr);
80cb32e2e7SValeria Barra   for (c=cStart, eoffset = 0; c<cEnd; c++) {
81cb32e2e7SValeria Barra     PetscInt numindices, *indices, i;
82cb32e2e7SValeria Barra     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
83cb32e2e7SValeria Barra                                    &indices, NULL); CHKERRQ(ierr);
84cb32e2e7SValeria Barra     for (i=0; i<numindices; i+=ncomp) {
85cb32e2e7SValeria Barra       for (PetscInt j=0; j<ncomp; j++) {
86cb32e2e7SValeria Barra         if (indices[i+j] != indices[i] + (PetscInt)(copysign(j, indices[i])))
87cb32e2e7SValeria Barra           SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
88cb32e2e7SValeria Barra                    "Cell %D closure indices not interlaced", c);
89cb32e2e7SValeria Barra       }
90cb32e2e7SValeria Barra       // NO BC on closed surfaces
91cb32e2e7SValeria Barra       PetscInt loc = indices[i];
92cb32e2e7SValeria Barra       erestrict[eoffset++] = loc/ncomp;
93cb32e2e7SValeria Barra     }
94cb32e2e7SValeria Barra     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
95cb32e2e7SValeria Barra                                        &indices, NULL); CHKERRQ(ierr);
96cb32e2e7SValeria Barra   }
97cb32e2e7SValeria Barra 
98cb32e2e7SValeria Barra   // Setup CEED restriction
99cb32e2e7SValeria Barra   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
100cb32e2e7SValeria Barra   ierr = VecGetLocalSize(Uloc, &nnodes); CHKERRQ(ierr);
101cb32e2e7SValeria Barra 
102cb32e2e7SValeria Barra   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
10361dbc9d2Sjeremylt   CeedElemRestrictionCreate(ceed, imode, nelem, P*P, nnodes/ncomp, ncomp,
104cb32e2e7SValeria Barra                             CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
105cb32e2e7SValeria Barra                             Erestrict);
106cb32e2e7SValeria Barra   ierr = PetscFree(erestrict); CHKERRQ(ierr);
107cb32e2e7SValeria Barra 
108cb32e2e7SValeria Barra   PetscFunctionReturn(0);
109cb32e2e7SValeria Barra }
110cb32e2e7SValeria Barra 
11132d2ee49SValeria Barra // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
11232d2ee49SValeria Barra static PetscErrorCode ProjectToUnitSphere(DM dm) {
11332d2ee49SValeria Barra   Vec            coordinates;
11432d2ee49SValeria Barra   PetscScalar   *coords;
11532d2ee49SValeria Barra   PetscInt       Nv, v, dim, d;
11632d2ee49SValeria Barra   PetscErrorCode ierr;
11732d2ee49SValeria Barra 
11832d2ee49SValeria Barra   PetscFunctionBeginUser;
11932d2ee49SValeria Barra   ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
12032d2ee49SValeria Barra   ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
12132d2ee49SValeria Barra   ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
12232d2ee49SValeria Barra   Nv  /= dim;
12332d2ee49SValeria Barra   ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
12432d2ee49SValeria Barra   for (v = 0; v < Nv; ++v) {
12532d2ee49SValeria Barra     PetscReal r = 0.0;
12632d2ee49SValeria Barra 
12732d2ee49SValeria Barra     for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
12832d2ee49SValeria Barra     r = PetscSqrtReal(r);
12932d2ee49SValeria Barra     for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
13032d2ee49SValeria Barra   }
13132d2ee49SValeria Barra   ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
13232d2ee49SValeria Barra   PetscFunctionReturn(0);
13332d2ee49SValeria Barra }
13432d2ee49SValeria Barra 
135cb32e2e7SValeria Barra int main(int argc, char **argv) {
136cb32e2e7SValeria Barra   PetscInt ierr;
137cb32e2e7SValeria Barra   MPI_Comm comm;
138cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
139cb0b5415Sjeremylt        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
140cb32e2e7SValeria Barra   PetscInt lsize, gsize, xlsize,
141cb32e2e7SValeria Barra            qextra  = 1, // default number of extra quadrature points
142cb32e2e7SValeria Barra            ncompx  = 3, // number of components of 3D physical coordinates
143cb32e2e7SValeria Barra            ncompu  = 1, // dimension of field to which apply mass operator
144cb32e2e7SValeria Barra            topodim = 2, // topological dimension of manifold
145cb32e2e7SValeria Barra            degree  = 3; // default degree for finite element bases
146cb32e2e7SValeria Barra   PetscBool read_mesh = PETSC_FALSE,
147cb32e2e7SValeria Barra             test_mode = PETSC_FALSE;
148cb32e2e7SValeria Barra   PetscSpace sp;
149cb32e2e7SValeria Barra   PetscFE fe;
150cb32e2e7SValeria Barra   Vec X, Xloc, V, Vloc;
151cb32e2e7SValeria Barra   DM  dm, dmcoord;
152cb32e2e7SValeria Barra   Ceed ceed;
153cb32e2e7SValeria Barra   CeedInt P, Q;
15432d2ee49SValeria Barra   CeedOperator op_setupgeo, op_apply;
15532d2ee49SValeria Barra   CeedQFunction qf_setupgeo, qf_apply;
15632d2ee49SValeria Barra   CeedBasis basisx, basisu;
157*15910d16Sjeremylt   CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi;
15832d2ee49SValeria Barra   PetscFunctionList geomfactorlist = NULL;
15932d2ee49SValeria Barra   char problemtype[PETSC_MAX_PATH_LEN] = "sphere";
160cb32e2e7SValeria Barra 
161cb32e2e7SValeria Barra   ierr = PetscInitialize(&argc, &argv, NULL, help);
162cb32e2e7SValeria Barra   if (ierr) return ierr;
163cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
164cb32e2e7SValeria Barra 
16532d2ee49SValeria Barra   // Set up problem type command line option
16632d2ee49SValeria Barra   ierr = PetscFunctionListAdd(&geomfactorlist, "cube", &SetupMassGeoCube);
16732d2ee49SValeria Barra   CHKERRQ(ierr);
16832d2ee49SValeria Barra   ierr = PetscFunctionListAdd(&geomfactorlist, "sphere", &SetupMassGeoSphere);
16932d2ee49SValeria Barra   CHKERRQ(ierr);
17032d2ee49SValeria Barra 
17132d2ee49SValeria Barra   // Read command line options
172ccf0fe6fSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc",
173ccf0fe6fSjeremylt                            NULL);
174cb32e2e7SValeria Barra   CHKERRQ(ierr);
17532d2ee49SValeria Barra   ierr = PetscOptionsFList("-problem", "Problem to solve", NULL, geomfactorlist,
17632d2ee49SValeria Barra                            problemtype, problemtype, sizeof problemtype, NULL);
17732d2ee49SValeria Barra   CHKERRQ(ierr);
178cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
179cb32e2e7SValeria Barra                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
180cb32e2e7SValeria Barra   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
181cb32e2e7SValeria Barra                             NULL, ceedresource, ceedresource,
182cb32e2e7SValeria Barra                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
183cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-test",
184cb32e2e7SValeria Barra                           "Testing mode (do not print unless error is large)",
185cb32e2e7SValeria Barra                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
186cb32e2e7SValeria Barra   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
187cb32e2e7SValeria Barra                             filename, filename, sizeof(filename), &read_mesh);
188cb32e2e7SValeria Barra   CHKERRQ(ierr);
189cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
190cb32e2e7SValeria Barra 
19132d2ee49SValeria Barra   // Setup function pointer for geometric factors
19232d2ee49SValeria Barra   int (*geomfp)(void *ctx, const CeedInt Q, const CeedScalar *const *in,
19332d2ee49SValeria Barra                 CeedScalar *const *out);
19432d2ee49SValeria Barra   ierr = PetscFunctionListFind(geomfactorlist, problemtype,
19532d2ee49SValeria Barra                                (void(* *)(void))&geomfp); CHKERRQ(ierr);
19632d2ee49SValeria Barra   const char *str;
19732d2ee49SValeria Barra   if (geomfp == SetupMassGeoCube)
19832d2ee49SValeria Barra     str = SetupMassGeoCube_loc;
19932d2ee49SValeria Barra   else if (geomfp == SetupMassGeoSphere)
20032d2ee49SValeria Barra     str = SetupMassGeoSphere_loc;
20132d2ee49SValeria Barra   else
20232d2ee49SValeria Barra     return CeedError(ceed, 1, "Function not found in the list");
20332d2ee49SValeria Barra 
204cb32e2e7SValeria Barra   // Setup DM
205cb32e2e7SValeria Barra   if (read_mesh) {
206cb32e2e7SValeria Barra     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
207cb32e2e7SValeria Barra     CHKERRQ(ierr);
208cb32e2e7SValeria Barra   } else {
209cb32e2e7SValeria Barra     // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box.
210cb32e2e7SValeria Barra     PetscBool simplex = PETSC_FALSE;
211cb32e2e7SValeria Barra     ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm);
212cb32e2e7SValeria Barra     CHKERRQ(ierr);
213cb32e2e7SValeria Barra     // Set the object name
21432d2ee49SValeria Barra     ierr = PetscObjectSetName((PetscObject)dm, problemtype); CHKERRQ(ierr);
215cb32e2e7SValeria Barra     // Distribute mesh over processes
216cb32e2e7SValeria Barra     {
217cb32e2e7SValeria Barra       DM dmDist = NULL;
218cb32e2e7SValeria Barra       PetscPartitioner part;
219cb32e2e7SValeria Barra 
220cb32e2e7SValeria Barra       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
221cb32e2e7SValeria Barra       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
222cb32e2e7SValeria Barra       ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
223cb32e2e7SValeria Barra       if (dmDist) {
224cb32e2e7SValeria Barra         ierr = DMDestroy(&dm); CHKERRQ(ierr);
225cb32e2e7SValeria Barra         dm  = dmDist;
226cb32e2e7SValeria Barra       }
227cb32e2e7SValeria Barra     }
22832d2ee49SValeria Barra     // Refine DMPlex with uniform refinement using runtime option -dm_refine
22932d2ee49SValeria Barra     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
23032d2ee49SValeria Barra     ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
23132d2ee49SValeria Barra     if (!strcmp(problemtype, "sphere"))
23232d2ee49SValeria Barra       ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);
233cb32e2e7SValeria Barra     // View DMPlex via runtime option
234cb32e2e7SValeria Barra     ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
235cb32e2e7SValeria Barra   }
236cb32e2e7SValeria Barra 
237cb32e2e7SValeria Barra   // Create FE
238cb32e2e7SValeria Barra   ierr = PetscFECreateDefault(PETSC_COMM_SELF, topodim, ncompu, PETSC_FALSE, NULL,
239cb32e2e7SValeria Barra                               PETSC_DETERMINE, &fe);
240cb32e2e7SValeria Barra   CHKERRQ(ierr);
241cb32e2e7SValeria Barra   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
242cb32e2e7SValeria Barra   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
243cb32e2e7SValeria Barra   ierr = DMCreateDS(dm); CHKERRQ(ierr);
244cb32e2e7SValeria Barra   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
245cb32e2e7SValeria Barra   CHKERRQ(ierr);
246cb32e2e7SValeria Barra 
247cb32e2e7SValeria Barra   // Get basis space degree
248cb32e2e7SValeria Barra   ierr = PetscFEGetBasisSpace(fe, &sp); CHKERRQ(ierr);
249cb32e2e7SValeria Barra   ierr = PetscSpaceGetDegree(sp, &degree, NULL); CHKERRQ(ierr);
250cb32e2e7SValeria Barra   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
251cb32e2e7SValeria Barra   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
252cb32e2e7SValeria Barra                              "-petscspace_degree %D must be at least 1", degree);
253cb32e2e7SValeria Barra 
254cb32e2e7SValeria Barra   // Create vectors
255cb32e2e7SValeria Barra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
256cb32e2e7SValeria Barra   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
257cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
258cb32e2e7SValeria Barra   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
259cb32e2e7SValeria Barra   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
260cb32e2e7SValeria Barra   ierr = VecDuplicate(X, &V); CHKERRQ(ierr);
261cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &Vloc); CHKERRQ(ierr);
262cb32e2e7SValeria Barra 
263cb32e2e7SValeria Barra   // Set up libCEED
264cb32e2e7SValeria Barra   CeedInit(ceedresource, &ceed);
265cb32e2e7SValeria Barra 
266cb32e2e7SValeria Barra   // Print summary
267cb32e2e7SValeria Barra   P = degree + 1;
268cb32e2e7SValeria Barra   Q = P + qextra;
269cb32e2e7SValeria Barra   const char *usedresource;
270cb32e2e7SValeria Barra   CeedGetResource(ceed, &usedresource);
271cb32e2e7SValeria Barra   if (!test_mode) {
272cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
27332d2ee49SValeria Barra                        "\n-- libCEED + PETSc Surface Area of a Manifold --\n"
274cb32e2e7SValeria Barra                        "  libCEED:\n"
275cb32e2e7SValeria Barra                        "    libCEED Backend                    : %s\n"
276cb32e2e7SValeria Barra                        "  Mesh:\n"
277cb32e2e7SValeria Barra                        "    Number of 1D Basis Nodes (p)       : %d\n"
278cb32e2e7SValeria Barra                        "    Number of 1D Quadrature Points (q) : %d\n"
279db419314Sjeremylt                        "    Global nodes                       : %D\n"
280db419314Sjeremylt                        "    DoF per node                       : %D\n",
281db419314Sjeremylt                        usedresource, P, Q,  gsize/ncompu, ncompu);
282cb32e2e7SValeria Barra     CHKERRQ(ierr);
283cb32e2e7SValeria Barra   }
284cb32e2e7SValeria Barra 
28532d2ee49SValeria Barra   // Setup libCEED's objects:
286cb32e2e7SValeria Barra   // Create bases
287cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q,
288cb32e2e7SValeria Barra                                   CEED_GAUSS, &basisu);
289cb32e2e7SValeria Barra   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q,
290cb32e2e7SValeria Barra                                   CEED_GAUSS, &basisx);
291cb32e2e7SValeria Barra 
292cb32e2e7SValeria Barra   // CEED restrictions
293cb32e2e7SValeria Barra   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
294cb32e2e7SValeria Barra   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
295cb32e2e7SValeria Barra   CHKERRQ(ierr);
296cb32e2e7SValeria Barra 
29774312f09Sjeremylt   CreateRestrictionPlex(ceed, CEED_INTERLACED, 2, ncompx, &Erestrictx, dmcoord);
298a8d32208Sjeremylt   CHKERRQ(ierr);
29974312f09Sjeremylt   CreateRestrictionPlex(ceed, CEED_INTERLACED, P, ncompu, &Erestrictu, dm);
300a8d32208Sjeremylt   CHKERRQ(ierr);
301cb32e2e7SValeria Barra 
302cb32e2e7SValeria Barra   CeedInt cStart, cEnd;
303cb32e2e7SValeria Barra   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
304cb32e2e7SValeria Barra   const CeedInt nelem = cEnd - cStart;
305cb32e2e7SValeria Barra 
3067509a596Sjeremylt   // CEED strided restrictions
307cb32e2e7SValeria Barra   const CeedInt qdatasize = 1;
3087509a596Sjeremylt   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q, nelem*Q*Q, qdatasize,
309523b8ea0Sjeremylt                                    CEED_STRIDES_BACKEND, &Erestrictqdi);
310cb32e2e7SValeria Barra 
311cb32e2e7SValeria Barra   // Element coordinates
312cb32e2e7SValeria Barra   Vec coords;
313cb32e2e7SValeria Barra   const PetscScalar *coordArray;
314cb32e2e7SValeria Barra   PetscSection section;
315cb32e2e7SValeria Barra   ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr);
316cb32e2e7SValeria Barra   ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr);
317cb32e2e7SValeria Barra   ierr = DMGetSection(dmcoord, &section); CHKERRQ(ierr);
318cb32e2e7SValeria Barra 
319cb32e2e7SValeria Barra   CeedVector xcoord;
320cb32e2e7SValeria Barra   CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL);
321cb32e2e7SValeria Barra   CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES,
322cb32e2e7SValeria Barra                      (PetscScalar *)coordArray);
323cb32e2e7SValeria Barra   ierr = VecRestoreArrayRead(coords, &coordArray);
324cb32e2e7SValeria Barra 
325cb32e2e7SValeria Barra   // Create the vectors that will be needed in setup and apply
326cb32e2e7SValeria Barra   CeedVector uceed, vceed, qdata;
327cb32e2e7SValeria Barra   CeedInt nqpts;
328cb32e2e7SValeria Barra   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
329cb32e2e7SValeria Barra   CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata);
330cb32e2e7SValeria Barra   CeedVectorCreate(ceed, xlsize, &uceed);
331cb32e2e7SValeria Barra   CeedVectorCreate(ceed, xlsize, &vceed);
332cb32e2e7SValeria Barra 
33332d2ee49SValeria Barra   // Create the Q-function that builds the operator for the geomteric factors
33432d2ee49SValeria Barra   //   (i.e., the quadrature data)
33532d2ee49SValeria Barra   CeedQFunctionCreateInterior(ceed, 1, geomfp, str, &qf_setupgeo);
33632d2ee49SValeria Barra   CeedQFunctionAddInput(qf_setupgeo, "x", ncompx, CEED_EVAL_INTERP);
337cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*topodim, CEED_EVAL_GRAD);
338cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT);
339cb32e2e7SValeria Barra   CeedQFunctionAddOutput(qf_setupgeo, "qdata", qdatasize, CEED_EVAL_NONE);
340cb32e2e7SValeria Barra 
341cb32e2e7SValeria Barra   // Set up the mass operator
342cb32e2e7SValeria Barra   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_apply);
343cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_apply, "u", ncompu, CEED_EVAL_INTERP);
344cb32e2e7SValeria Barra   CeedQFunctionAddInput(qf_apply, "qdata", qdatasize, CEED_EVAL_NONE);
345cb32e2e7SValeria Barra   CeedQFunctionAddOutput(qf_apply, "v", ncompu, CEED_EVAL_INTERP);
346cb32e2e7SValeria Barra 
347cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the operator
34832d2ee49SValeria Barra   CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo);
349a8d32208Sjeremylt   CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx,
350a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
351a8d32208Sjeremylt   CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx,
352a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
353*15910d16Sjeremylt   CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
354a8d32208Sjeremylt                        CEED_VECTOR_NONE);
355a8d32208Sjeremylt   CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi,
356cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
357cb32e2e7SValeria Barra 
358cb32e2e7SValeria Barra   // Create the mass operator
35932d2ee49SValeria Barra   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
360a8d32208Sjeremylt   CeedOperatorSetField(op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
361a8d32208Sjeremylt   CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED,
362a8d32208Sjeremylt                        qdata);
363a8d32208Sjeremylt   CeedOperatorSetField(op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
364cb32e2e7SValeria Barra 
365cb32e2e7SValeria Barra   // Compute the quadrature data for the mass operator
366cb32e2e7SValeria Barra   CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
367cb32e2e7SValeria Barra 
368cb32e2e7SValeria Barra   PetscScalar *v;
369cb32e2e7SValeria Barra   ierr = VecZeroEntries(Vloc); CHKERRQ(ierr);
370cb32e2e7SValeria Barra   ierr = VecGetArray(Vloc, &v);
371cb32e2e7SValeria Barra   CeedVectorSetArray(vceed, CEED_MEM_HOST, CEED_USE_POINTER, v);
372cb32e2e7SValeria Barra 
373cb32e2e7SValeria Barra   // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1
374cb32e2e7SValeria Barra   if (!test_mode) {
375cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
37632d2ee49SValeria Barra                        "Computing the mesh area using the formula: area = 1^T M 1\n");
377cb32e2e7SValeria Barra     CHKERRQ(ierr);
378cb32e2e7SValeria Barra   }
379cb32e2e7SValeria Barra 
380cb32e2e7SValeria Barra   // Initialize u and v with ones
381cb32e2e7SValeria Barra   CeedVectorSetValue(uceed, 1.0);
382cb32e2e7SValeria Barra 
383cb32e2e7SValeria Barra   // Apply the mass operator: 'u' -> 'v'
384cb32e2e7SValeria Barra   CeedOperatorApply(op_apply, uceed, vceed, CEED_REQUEST_IMMEDIATE);
385cb32e2e7SValeria Barra   CeedVectorSyncArray(vceed, CEED_MEM_HOST);
386cb32e2e7SValeria Barra 
387cb32e2e7SValeria Barra   // Gather output vector
388cb32e2e7SValeria Barra   ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr);
389cb32e2e7SValeria Barra   ierr = VecZeroEntries(V); CHKERRQ(ierr);
390cb32e2e7SValeria Barra   ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
391cb32e2e7SValeria Barra   ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
392cb32e2e7SValeria Barra 
393cb32e2e7SValeria Barra   // Compute and print the sum of the entries of 'v' giving the mesh surface area
394cb32e2e7SValeria Barra   PetscScalar area;
395cb32e2e7SValeria Barra   ierr = VecSum(V, &area); CHKERRQ(ierr);
396cb32e2e7SValeria Barra 
397cb32e2e7SValeria Barra   // Compute the exact surface area and print the result
39832d2ee49SValeria Barra   CeedScalar exact_surfarea = 4 * M_PI;
39932d2ee49SValeria Barra   if (!strcmp(problemtype, "cube")) {
40032d2ee49SValeria Barra     PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube
40132d2ee49SValeria Barra     exact_surfarea = 6 * (2*l) * (2*l);
40232d2ee49SValeria Barra   }
40332d2ee49SValeria Barra 
404cb32e2e7SValeria Barra   if (!test_mode) {
405cb32e2e7SValeria Barra     ierr = PetscPrintf(comm, "Exact mesh surface area    : % .14g\n",
406cb32e2e7SValeria Barra                        exact_surfarea); CHKERRQ(ierr);
407cb32e2e7SValeria Barra     ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area);
408cb32e2e7SValeria Barra     CHKERRQ(ierr);
409cb32e2e7SValeria Barra     ierr = PetscPrintf(comm, "Area error                 : % .14g\n",
410cb32e2e7SValeria Barra                        fabs(area - exact_surfarea)); CHKERRQ(ierr);
411cb32e2e7SValeria Barra   }
412cb32e2e7SValeria Barra 
413cb32e2e7SValeria Barra   // PETSc cleanup
414cb32e2e7SValeria Barra   ierr = DMDestroy(&dm); CHKERRQ(ierr);
415cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
416cb32e2e7SValeria Barra   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
417cb32e2e7SValeria Barra   ierr = VecDestroy(&V); CHKERRQ(ierr);
418cb32e2e7SValeria Barra   ierr = VecDestroy(&Vloc); CHKERRQ(ierr);
419cb32e2e7SValeria Barra 
420cb32e2e7SValeria Barra   // libCEED cleanup
421cb32e2e7SValeria Barra   CeedQFunctionDestroy(&qf_setupgeo);
422cb32e2e7SValeria Barra   CeedOperatorDestroy(&op_setupgeo);
423cb32e2e7SValeria Barra   CeedVectorDestroy(&xcoord);
424cb32e2e7SValeria Barra   CeedVectorDestroy(&uceed);
425cb32e2e7SValeria Barra   CeedVectorDestroy(&vceed);
426cb32e2e7SValeria Barra   CeedVectorDestroy(&qdata);
427cb32e2e7SValeria Barra   CeedBasisDestroy(&basisx);
428cb32e2e7SValeria Barra   CeedBasisDestroy(&basisu);
429cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictu);
430cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictx);
431cb32e2e7SValeria Barra   CeedElemRestrictionDestroy(&Erestrictqdi);
432cb32e2e7SValeria Barra   CeedQFunctionDestroy(&qf_apply);
433cb32e2e7SValeria Barra   CeedOperatorDestroy(&op_apply);
434cb32e2e7SValeria Barra   CeedDestroy(&ceed);
435cb32e2e7SValeria Barra   return PetscFinalize();
436cb32e2e7SValeria Barra }
437