xref: /libCEED/examples/ceed/ex1-volume.c (revision d1d35e2f02dc969aee8debf3fd943dd784aa847a)
166087c08SValeria Barra // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC.
266087c08SValeria Barra // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707.
366087c08SValeria Barra // All Rights reserved. See files LICENSE and NOTICE for details.
466087c08SValeria Barra //
566087c08SValeria Barra // This file is part of CEED, a collection of benchmarks, miniapps, software
666087c08SValeria Barra // libraries and APIs for efficient high-order finite element and spectral
766087c08SValeria Barra // element discretizations for exascale applications. For more information and
866087c08SValeria Barra // source code availability see http://github.com/ceed.
966087c08SValeria Barra //
1066087c08SValeria Barra // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
1166087c08SValeria Barra // a collaborative effort of two U.S. Department of Energy organizations (Office
1266087c08SValeria Barra // of Science and the National Nuclear Security Administration) responsible for
1366087c08SValeria Barra // the planning and preparation of a capable exascale ecosystem, including
1466087c08SValeria Barra // software, applications, hardware, advanced system engineering and early
1566087c08SValeria Barra // testbed platforms, in support of the nation's exascale computing imperative.
1666087c08SValeria Barra 
1766087c08SValeria Barra //                             libCEED Example 1
1866087c08SValeria Barra //
1966087c08SValeria Barra // This example illustrates a simple usage of libCEED to compute the volume of a
2066087c08SValeria Barra // 3D body using matrix-free application of a mass operator.  Arbitrary mesh and
21ded9b81dSJeremy L Thompson // solution degrees in 1D, 2D and 3D are supported from the same code.
2266087c08SValeria Barra //
2366087c08SValeria Barra // The example has no dependencies, and is designed to be self-contained. For
2466087c08SValeria Barra // additional examples that use external discretization libraries (MFEM, PETSc,
2566087c08SValeria Barra // etc.) see the subdirectories in libceed/examples.
2666087c08SValeria Barra //
2766087c08SValeria Barra // All libCEED objects use a Ceed device object constructed based on a command
2866087c08SValeria Barra // line argument (-ceed).
2966087c08SValeria Barra //
3066087c08SValeria Barra // Build with:
3166087c08SValeria Barra //
3266087c08SValeria Barra //     make ex1-volume [CEED_DIR=</path/to/libceed>]
3366087c08SValeria Barra //
3466087c08SValeria Barra // Sample runs:
3566087c08SValeria Barra //
3666087c08SValeria Barra //     ./ex1-volume
3766087c08SValeria Barra //     ./ex1-volume -ceed /cpu/self
3828688798Sjeremylt //     ./ex1-volume -ceed /gpu/cuda
3966087c08SValeria Barra //
4066087c08SValeria Barra // Next line is grep'd from tap.sh to set its arguments
4166087c08SValeria Barra // Test in 1D-3D
42dc8efd83SLeila Ghaffari //TESTARGS(name="1D_User_QFunction") -ceed {ceed_resource} -d 1 -t
43dc8efd83SLeila Ghaffari //TESTARGS(name="2D_User_QFunction") -ceed {ceed_resource} -d 2 -t
44dc8efd83SLeila Ghaffari //TESTARGS(name="3D_User_QFunction") -ceed {ceed_resource} -d 3 -t
45dc8efd83SLeila Ghaffari //TESTARGS(name="1D_Gallery_QFunction") -ceed {ceed_resource} -d 1 -t -g
46dc8efd83SLeila Ghaffari //TESTARGS(name="2D_Gallery_QFunction") -ceed {ceed_resource} -d 2 -t -g
47dc8efd83SLeila Ghaffari //TESTARGS(name="3D_Gallery_QFunction") -ceed {ceed_resource} -d 3 -t -g
4866087c08SValeria Barra 
4966087c08SValeria Barra /// @file
5066087c08SValeria Barra /// libCEED example using mass operator to compute volume
5166087c08SValeria Barra 
5266087c08SValeria Barra #include <ceed.h>
5366087c08SValeria Barra #include <math.h>
543d576824SJeremy L Thompson #include <stdlib.h>
5566087c08SValeria Barra #include <string.h>
5666087c08SValeria Barra #include "ex1-volume.h"
5766087c08SValeria Barra 
5866087c08SValeria Barra // Auxiliary functions.
59*d1d35e2fSjeremylt int GetCartesianMeshSize(int dim, int degree, int prob_size, int num_xyz[3]);
60*d1d35e2fSjeremylt int BuildCartesianRestriction(Ceed ceed, int dim, int num_xyz[3], int degree,
61*d1d35e2fSjeremylt                               int num_comp, CeedInt *size, CeedInt num_qpts,
6266087c08SValeria Barra                               CeedElemRestriction *restr,
6366087c08SValeria Barra                               CeedElemRestriction *restr_i);
64*d1d35e2fSjeremylt int SetCartesianMeshCoords(int dim, int num_xyz[3], int mesh_degree,
6566087c08SValeria Barra                            CeedVector mesh_coords);
6666087c08SValeria Barra CeedScalar TransformMeshCoords(int dim, int mesh_size, CeedVector mesh_coords);
6766087c08SValeria Barra 
6866087c08SValeria Barra int main(int argc, const char *argv[]) {
6966087c08SValeria Barra   const char *ceed_spec = "/cpu/self";
7066087c08SValeria Barra   int dim         = 3;              // dimension of the mesh
71*d1d35e2fSjeremylt   int num_comp_x  = 3;              // number of x components
72ded9b81dSJeremy L Thompson   int mesh_degree = 4;              // polynomial degree for the mesh
73ded9b81dSJeremy L Thompson   int sol_degree  = 4;              // polynomial degree for the solution
74ded9b81dSJeremy L Thompson   int num_qpts    = sol_degree + 2; // number of 1D quadrature points
7566087c08SValeria Barra   int prob_size   = -1;             // approximate problem size
7666087c08SValeria Barra   int help = 0, test = 0, gallery = 0;
7766087c08SValeria Barra 
7866087c08SValeria Barra   // Process command line arguments.
7966087c08SValeria Barra   for (int ia = 1; ia < argc; ia++) {
80ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
8166087c08SValeria Barra     int next_arg = ((ia+1) < argc), parse_error = 0;
8266087c08SValeria Barra     if (!strcmp(argv[ia],"-h")) {
8366087c08SValeria Barra       help = 1;
8466087c08SValeria Barra     } else if (!strcmp(argv[ia],"-c") || !strcmp(argv[ia],"-ceed")) {
8566087c08SValeria Barra       parse_error = next_arg ? ceed_spec = argv[++ia], 0 : 1;
8666087c08SValeria Barra     } else if (!strcmp(argv[ia],"-d")) {
8766087c08SValeria Barra       parse_error = next_arg ? dim = atoi(argv[++ia]), 0 : 1;
88*d1d35e2fSjeremylt       num_comp_x = dim;
8966087c08SValeria Barra     } else if (!strcmp(argv[ia],"-m")) {
90ded9b81dSJeremy L Thompson       parse_error = next_arg ? mesh_degree = atoi(argv[++ia]), 0 : 1;
91ded9b81dSJeremy L Thompson     } else if (!strcmp(argv[ia],"-p")) {
92ded9b81dSJeremy L Thompson       parse_error = next_arg ? sol_degree= atoi(argv[++ia]), 0 : 1;
9366087c08SValeria Barra     } else if (!strcmp(argv[ia],"-q")) {
9466087c08SValeria Barra       parse_error = next_arg ? num_qpts = atoi(argv[++ia]), 0 : 1;
9566087c08SValeria Barra     } else if (!strcmp(argv[ia],"-s")) {
9666087c08SValeria Barra       parse_error = next_arg ? prob_size = atoi(argv[++ia]), 0 : 1;
9766087c08SValeria Barra     } else if (!strcmp(argv[ia],"-t")) {
9866087c08SValeria Barra       test = 1;
9966087c08SValeria Barra     } else if (!strcmp(argv[ia],"-g")) {
10066087c08SValeria Barra       gallery = 1;
10166087c08SValeria Barra     }
10266087c08SValeria Barra     if (parse_error) {
10366087c08SValeria Barra       printf("Error parsing command line options.\n");
10466087c08SValeria Barra       return 1;
10566087c08SValeria Barra     }
106ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
10766087c08SValeria Barra   }
10866087c08SValeria Barra   if (prob_size < 0) prob_size = test ? 8*16 : 256*1024;
10966087c08SValeria Barra 
11066087c08SValeria Barra   // Print the values of all options:
11166087c08SValeria Barra   if (!test || help) {
112ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
11366087c08SValeria Barra     printf("Selected options: [command line option] : <current value>\n");
11466087c08SValeria Barra     printf("  Ceed specification [-c] : %s\n", ceed_spec);
11566087c08SValeria Barra     printf("  Mesh dimension     [-d] : %d\n", dim);
116ded9b81dSJeremy L Thompson     printf("  Mesh degree        [-m] : %d\n", mesh_degree);
117ded9b81dSJeremy L Thompson     printf("  Solution degree    [-p] : %d\n", sol_degree);
11866087c08SValeria Barra     printf("  Num. 1D quadr. pts [-q] : %d\n", num_qpts);
11966087c08SValeria Barra     printf("  Approx. # unknowns [-s] : %d\n", prob_size);
12066087c08SValeria Barra     printf("  QFunction source   [-g] : %s\n", gallery?"gallery":"header");
12166087c08SValeria Barra     if (help) {
12266087c08SValeria Barra       printf("Test/quiet mode is %s\n", (test?"ON":"OFF (use -t to enable)"));
12366087c08SValeria Barra       return 0;
12466087c08SValeria Barra     }
12566087c08SValeria Barra     printf("\n");
126ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
12766087c08SValeria Barra   }
12866087c08SValeria Barra 
12966087c08SValeria Barra   // Select appropriate backend and logical device based on the <ceed-spec>
13066087c08SValeria Barra   // command line argument.
13166087c08SValeria Barra   Ceed ceed;
13266087c08SValeria Barra   CeedInit(ceed_spec, &ceed);
13366087c08SValeria Barra 
13466087c08SValeria Barra   // Construct the mesh and solution bases.
13566087c08SValeria Barra   CeedBasis mesh_basis, sol_basis;
136*d1d35e2fSjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1,
137*d1d35e2fSjeremylt                                   num_qpts, CEED_GAUSS, &mesh_basis);
138ded9b81dSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts,
13966087c08SValeria Barra                                   CEED_GAUSS, &sol_basis);
14066087c08SValeria Barra 
14166087c08SValeria Barra   // Determine the mesh size based on the given approximate problem size.
142*d1d35e2fSjeremylt   int num_xyz[dim];
143*d1d35e2fSjeremylt   GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz);
14466087c08SValeria Barra   if (!test) {
145ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
146*d1d35e2fSjeremylt     printf("Mesh size: nx = %d", num_xyz[0]);
147*d1d35e2fSjeremylt     if (dim > 1) { printf(", ny = %d", num_xyz[1]); }
148*d1d35e2fSjeremylt     if (dim > 2) { printf(", nz = %d", num_xyz[2]); }
14966087c08SValeria Barra     printf("\n");
150ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
15166087c08SValeria Barra   }
15266087c08SValeria Barra 
15366087c08SValeria Barra   // Build CeedElemRestriction objects describing the mesh and solution discrete
15466087c08SValeria Barra   // representations.
15566087c08SValeria Barra   CeedInt mesh_size, sol_size;
15615910d16Sjeremylt   CeedElemRestriction mesh_restr, sol_restr, sol_restr_i;
157*d1d35e2fSjeremylt   BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x,
158*d1d35e2fSjeremylt                             &mesh_size, num_qpts, &mesh_restr, NULL);
159*d1d35e2fSjeremylt   BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size,
16066087c08SValeria Barra                             num_qpts, &sol_restr, &sol_restr_i);
16166087c08SValeria Barra   if (!test) {
162e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
16366087c08SValeria Barra     printf("Number of mesh nodes     : %d\n", mesh_size/dim);
16466087c08SValeria Barra     printf("Number of solution nodes : %d\n", sol_size);
165e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
16666087c08SValeria Barra   }
16766087c08SValeria Barra 
16866087c08SValeria Barra   // Create a CeedVector with the mesh coordinates.
16966087c08SValeria Barra   CeedVector mesh_coords;
17066087c08SValeria Barra   CeedVectorCreate(ceed, mesh_size, &mesh_coords);
171*d1d35e2fSjeremylt   SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords);
17266087c08SValeria Barra 
17366087c08SValeria Barra   // Apply a transformation to the mesh.
17466087c08SValeria Barra   CeedScalar exact_vol = TransformMeshCoords(dim, mesh_size, mesh_coords);
17566087c08SValeria Barra 
176ded9b81dSJeremy L Thompson   // Context data to be passed to the 'f_build_mass' QFunction.
177777ff853SJeremy L Thompson   CeedQFunctionContext build_ctx;
178777ff853SJeremy L Thompson   struct BuildContext build_ctx_data;
179777ff853SJeremy L Thompson   build_ctx_data.dim = build_ctx_data.space_dim = dim;
180777ff853SJeremy L Thompson   CeedQFunctionContextCreate(ceed, &build_ctx);
181777ff853SJeremy L Thompson   CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER,
182777ff853SJeremy L Thompson                               sizeof(build_ctx_data), &build_ctx_data);
18366087c08SValeria Barra 
184ded9b81dSJeremy L Thompson   // Create the QFunction that builds the mass operator (i.e. computes its
18566087c08SValeria Barra   // quadrature data) and set its context data.
186*d1d35e2fSjeremylt   CeedQFunction qf_build;
18766087c08SValeria Barra   switch (gallery) {
18866087c08SValeria Barra   case 0:
18966087c08SValeria Barra     // This creates the QFunction directly.
19066087c08SValeria Barra     CeedQFunctionCreateInterior(ceed, 1, f_build_mass,
191*d1d35e2fSjeremylt                                 f_build_mass_loc, &qf_build);
192*d1d35e2fSjeremylt     CeedQFunctionAddInput(qf_build, "dx", num_comp_x*dim, CEED_EVAL_GRAD);
193*d1d35e2fSjeremylt     CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
194*d1d35e2fSjeremylt     CeedQFunctionAddOutput(qf_build, "qdata", 1, CEED_EVAL_NONE);
195*d1d35e2fSjeremylt     CeedQFunctionSetContext(qf_build, build_ctx);
19666087c08SValeria Barra     break;
19766087c08SValeria Barra   case 1: {
19866087c08SValeria Barra     // This creates the QFunction via the gallery.
19966087c08SValeria Barra     char name[13] = "";
20066087c08SValeria Barra     snprintf(name, sizeof name, "Mass%dDBuild", dim);
201*d1d35e2fSjeremylt     CeedQFunctionCreateInteriorByName(ceed, name, &qf_build);
20266087c08SValeria Barra     break;
20366087c08SValeria Barra   }
20466087c08SValeria Barra   }
20566087c08SValeria Barra 
20666087c08SValeria Barra   // Create the operator that builds the quadrature data for the mass operator.
207*d1d35e2fSjeremylt   CeedOperator op_build;
208*d1d35e2fSjeremylt   CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE,
209*d1d35e2fSjeremylt                      CEED_QFUNCTION_NONE, &op_build);
210*d1d35e2fSjeremylt   CeedOperatorSetField(op_build, "dx", mesh_restr, mesh_basis,
211a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
212*d1d35e2fSjeremylt   CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE,
21315910d16Sjeremylt                        mesh_basis, CEED_VECTOR_NONE);
214*d1d35e2fSjeremylt   CeedOperatorSetField(op_build, "qdata", sol_restr_i, CEED_BASIS_COLLOCATED,
215a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
21666087c08SValeria Barra 
21766087c08SValeria Barra   // Compute the quadrature data for the mass operator.
218*d1d35e2fSjeremylt   CeedVector q_data;
21966087c08SValeria Barra   CeedInt elem_qpts = CeedIntPow(num_qpts, dim);
22066087c08SValeria Barra   CeedInt num_elem = 1;
22166087c08SValeria Barra   for (int d = 0; d < dim; d++)
222*d1d35e2fSjeremylt     num_elem *= num_xyz[d];
223*d1d35e2fSjeremylt   CeedVectorCreate(ceed, num_elem*elem_qpts, &q_data);
224*d1d35e2fSjeremylt   CeedOperatorApply(op_build, mesh_coords, q_data,
22566087c08SValeria Barra                     CEED_REQUEST_IMMEDIATE);
22666087c08SValeria Barra 
227ded9b81dSJeremy L Thompson   // Create the QFunction that defines the action of the mass operator.
228*d1d35e2fSjeremylt   CeedQFunction qf_apply;
22966087c08SValeria Barra   switch (gallery) {
23066087c08SValeria Barra   case 0:
23166087c08SValeria Barra     // This creates the QFunction directly.
23266087c08SValeria Barra     CeedQFunctionCreateInterior(ceed, 1, f_apply_mass,
233*d1d35e2fSjeremylt                                 f_apply_mass_loc, &qf_apply);
234*d1d35e2fSjeremylt     CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP);
235*d1d35e2fSjeremylt     CeedQFunctionAddInput(qf_apply, "qdata", 1, CEED_EVAL_NONE);
236*d1d35e2fSjeremylt     CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP);
23766087c08SValeria Barra     break;
23866087c08SValeria Barra   case 1:
23966087c08SValeria Barra     // This creates the QFunction via the gallery.
240*d1d35e2fSjeremylt     CeedQFunctionCreateInteriorByName(ceed, "MassApply", &qf_apply);
24166087c08SValeria Barra     break;
24266087c08SValeria Barra   }
24366087c08SValeria Barra 
24466087c08SValeria Barra   // Create the mass operator.
245*d1d35e2fSjeremylt   CeedOperator op_apply;
246*d1d35e2fSjeremylt   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE,
247*d1d35e2fSjeremylt                      CEED_QFUNCTION_NONE, &op_apply);
248*d1d35e2fSjeremylt   CeedOperatorSetField(op_apply, "u", sol_restr, sol_basis, CEED_VECTOR_ACTIVE);
249*d1d35e2fSjeremylt   CeedOperatorSetField(op_apply, "qdata", sol_restr_i, CEED_BASIS_COLLOCATED,
250*d1d35e2fSjeremylt                        q_data);
251*d1d35e2fSjeremylt   CeedOperatorSetField(op_apply, "v", sol_restr, sol_basis, CEED_VECTOR_ACTIVE);
25266087c08SValeria Barra 
25366087c08SValeria Barra   // Create auxiliary solution-size vectors.
25466087c08SValeria Barra   CeedVector u, v;
25566087c08SValeria Barra   CeedVectorCreate(ceed, sol_size, &u);
25666087c08SValeria Barra   CeedVectorCreate(ceed, sol_size, &v);
25766087c08SValeria Barra 
25866087c08SValeria Barra   // Initialize 'u' and 'v' with ones.
25966087c08SValeria Barra   CeedVectorSetValue(u, 1.0);
26066087c08SValeria Barra 
261ded9b81dSJeremy L Thompson   // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1
262*d1d35e2fSjeremylt   CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
26366087c08SValeria Barra 
26466087c08SValeria Barra   // Compute and print the sum of the entries of 'v' giving the mesh volume.
265*d1d35e2fSjeremylt   const CeedScalar *v_array;
266*d1d35e2fSjeremylt   CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
26766087c08SValeria Barra   CeedScalar vol = 0.;
26866087c08SValeria Barra   for (CeedInt i = 0; i < sol_size; i++) {
269*d1d35e2fSjeremylt     vol += v_array[i];
27066087c08SValeria Barra   }
271*d1d35e2fSjeremylt   CeedVectorRestoreArrayRead(v, &v_array);
27266087c08SValeria Barra   if (!test) {
273ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
27466087c08SValeria Barra     printf(" done.\n");
27566087c08SValeria Barra     printf("Exact mesh volume    : % .14g\n", exact_vol);
27666087c08SValeria Barra     printf("Computed mesh volume : % .14g\n", vol);
27766087c08SValeria Barra     printf("Volume error         : % .14g\n", vol-exact_vol);
278ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
27966087c08SValeria Barra   } else {
28066087c08SValeria Barra     CeedScalar tol = (dim==1 ? 1E-14 : dim==2 ? 1E-7 : 1E-5);
28166087c08SValeria Barra     if (fabs(vol-exact_vol)>tol)
282e15f9bd0SJeremy L Thompson       // LCOV_EXCL_START
28366087c08SValeria Barra       printf("Volume error : % .1e\n", vol-exact_vol);
284e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
28566087c08SValeria Barra   }
28666087c08SValeria Barra 
28766087c08SValeria Barra   // Free dynamically allocated memory.
28866087c08SValeria Barra   CeedVectorDestroy(&u);
28966087c08SValeria Barra   CeedVectorDestroy(&v);
290*d1d35e2fSjeremylt   CeedVectorDestroy(&q_data);
29166087c08SValeria Barra   CeedVectorDestroy(&mesh_coords);
292*d1d35e2fSjeremylt   CeedOperatorDestroy(&op_apply);
293*d1d35e2fSjeremylt   CeedQFunctionDestroy(&qf_apply);
294777ff853SJeremy L Thompson   CeedQFunctionContextDestroy(&build_ctx);
295*d1d35e2fSjeremylt   CeedOperatorDestroy(&op_build);
296*d1d35e2fSjeremylt   CeedQFunctionDestroy(&qf_build);
29766087c08SValeria Barra   CeedElemRestrictionDestroy(&sol_restr);
29866087c08SValeria Barra   CeedElemRestrictionDestroy(&mesh_restr);
29966087c08SValeria Barra   CeedElemRestrictionDestroy(&sol_restr_i);
30066087c08SValeria Barra   CeedBasisDestroy(&sol_basis);
30166087c08SValeria Barra   CeedBasisDestroy(&mesh_basis);
30266087c08SValeria Barra   CeedDestroy(&ceed);
30366087c08SValeria Barra   return 0;
30466087c08SValeria Barra }
30566087c08SValeria Barra 
306*d1d35e2fSjeremylt int GetCartesianMeshSize(int dim, int degree, int prob_size, int num_xyz[dim]) {
30766087c08SValeria Barra   // Use the approximate formula:
308ded9b81dSJeremy L Thompson   //    prob_size ~ num_elem * degree^dim
309ded9b81dSJeremy L Thompson   CeedInt num_elem = prob_size / CeedIntPow(degree, dim);
31066087c08SValeria Barra   CeedInt s = 0;  // find s: num_elem/2 < 2^s <= num_elem
31166087c08SValeria Barra   while (num_elem > 1) {
31266087c08SValeria Barra     num_elem /= 2;
31366087c08SValeria Barra     s++;
31466087c08SValeria Barra   }
31566087c08SValeria Barra   CeedInt r = s%dim;
31666087c08SValeria Barra   for (int d = 0; d < dim; d++) {
31766087c08SValeria Barra     int sd = s/dim;
31866087c08SValeria Barra     if (r > 0) { sd++; r--; }
319*d1d35e2fSjeremylt     num_xyz[d] = 1 << sd;
32066087c08SValeria Barra   }
32166087c08SValeria Barra   return 0;
32266087c08SValeria Barra }
32366087c08SValeria Barra 
324*d1d35e2fSjeremylt int BuildCartesianRestriction(Ceed ceed, int dim, int num_xyz[dim], int degree,
325*d1d35e2fSjeremylt                               int num_comp, CeedInt *size, CeedInt num_qpts,
32666087c08SValeria Barra                               CeedElemRestriction *restr,
32766087c08SValeria Barra                               CeedElemRestriction *restr_i) {
328ded9b81dSJeremy L Thompson   CeedInt p = degree + 1;
329*d1d35e2fSjeremylt   CeedInt num_nodes = CeedIntPow(p, dim); // number of scalar nodes per element
33066087c08SValeria Barra   CeedInt elem_qpts = CeedIntPow(num_qpts, dim); // number of qpts per element
33166087c08SValeria Barra   CeedInt nd[3], num_elem = 1, scalar_size = 1;
33266087c08SValeria Barra   for (int d = 0; d < dim; d++) {
333*d1d35e2fSjeremylt     num_elem *= num_xyz[d];
334*d1d35e2fSjeremylt     nd[d] = num_xyz[d] * (p - 1) + 1;
33566087c08SValeria Barra     scalar_size *= nd[d];
33666087c08SValeria Barra   }
337*d1d35e2fSjeremylt   *size = scalar_size*num_comp;
33866087c08SValeria Barra   // elem:         0             1                 n-1
33966087c08SValeria Barra   //           |---*-...-*---|---*-...-*---|- ... -|--...--|
340*d1d35e2fSjeremylt   // num_nodes:   0   1    p-1  p  p+1       2*p             n*p
341*d1d35e2fSjeremylt   CeedInt *elem_nodes = malloc(sizeof(CeedInt)*num_elem*num_nodes);
34266087c08SValeria Barra   for (CeedInt e = 0; e < num_elem; e++) {
343*d1d35e2fSjeremylt     CeedInt e_xyz[3] = {1, 1, 1}, re = e;
344*d1d35e2fSjeremylt     for (int d = 0; d < dim; d++) { e_xyz[d] = re % num_xyz[d]; re /= num_xyz[d]; }
345*d1d35e2fSjeremylt     CeedInt *loc_el_nodes = elem_nodes + e*num_nodes;
346*d1d35e2fSjeremylt     for (int l_nodes = 0; l_nodes < num_nodes; l_nodes++) {
347*d1d35e2fSjeremylt       CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes;
34866087c08SValeria Barra       for (int d = 0; d < dim; d++) {
349*d1d35e2fSjeremylt         g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride;
350*d1d35e2fSjeremylt         g_nodes_stride *= nd[d];
351*d1d35e2fSjeremylt         r_nodes /= p;
35266087c08SValeria Barra       }
353*d1d35e2fSjeremylt       loc_el_nodes[l_nodes] = g_nodes;
35466087c08SValeria Barra     }
35566087c08SValeria Barra   }
356*d1d35e2fSjeremylt   CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size,
357*d1d35e2fSjeremylt                             num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES,
358*d1d35e2fSjeremylt                             elem_nodes, restr);
35915910d16Sjeremylt   if (restr_i)
3607509a596Sjeremylt     CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts,
361*d1d35e2fSjeremylt                                      num_comp, num_comp * elem_qpts * num_elem,
362523b8ea0Sjeremylt                                      CEED_STRIDES_BACKEND, restr_i);
363*d1d35e2fSjeremylt   free(elem_nodes);
36466087c08SValeria Barra   return 0;
36566087c08SValeria Barra }
36666087c08SValeria Barra 
367*d1d35e2fSjeremylt int SetCartesianMeshCoords(int dim, int num_xyz[dim], int mesh_degree,
36866087c08SValeria Barra                            CeedVector mesh_coords) {
369ded9b81dSJeremy L Thompson   CeedInt p = mesh_degree + 1;
37066087c08SValeria Barra   CeedInt nd[3], num_elem = 1, scalar_size = 1;
37166087c08SValeria Barra   for (int d = 0; d < dim; d++) {
372*d1d35e2fSjeremylt     num_elem *= num_xyz[d];
373*d1d35e2fSjeremylt     nd[d] = num_xyz[d] * (p - 1) + 1;
37466087c08SValeria Barra     scalar_size *= nd[d];
37566087c08SValeria Barra   }
37666087c08SValeria Barra   CeedScalar *coords;
37766087c08SValeria Barra   CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
378ded9b81dSJeremy L Thompson   CeedScalar *nodes = malloc(sizeof(CeedScalar) * p);
37966087c08SValeria Barra   // The H1 basis uses Lobatto quadrature points as nodes.
380ded9b81dSJeremy L Thompson   CeedLobattoQuadrature(p, nodes, NULL); // nodes are in [-1,1]
381ded9b81dSJeremy L Thompson   for (CeedInt i = 0; i < p; i++) { nodes[i] = 0.5 + 0.5 * nodes[i]; }
382*d1d35e2fSjeremylt   for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) {
383*d1d35e2fSjeremylt     CeedInt r_nodes = gs_nodes;
38466087c08SValeria Barra     for (int d = 0; d < dim; d++) {
385*d1d35e2fSjeremylt       CeedInt d_1d = r_nodes % nd[d];
386*d1d35e2fSjeremylt       coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d %
387*d1d35e2fSjeremylt                                             (p - 1)]) / num_xyz[d];
388*d1d35e2fSjeremylt       r_nodes /= nd[d];
38966087c08SValeria Barra     }
39066087c08SValeria Barra   }
39166087c08SValeria Barra   free(nodes);
39266087c08SValeria Barra   CeedVectorRestoreArray(mesh_coords, &coords);
39366087c08SValeria Barra   return 0;
39466087c08SValeria Barra }
39566087c08SValeria Barra 
39666087c08SValeria Barra #ifndef M_PI
39766087c08SValeria Barra #define M_PI    3.14159265358979323846
39866087c08SValeria Barra #define M_PI_2  1.57079632679489661923
39966087c08SValeria Barra #endif
40066087c08SValeria Barra 
40166087c08SValeria Barra CeedScalar TransformMeshCoords(int dim, int mesh_size, CeedVector mesh_coords) {
40266087c08SValeria Barra   CeedScalar exact_volume;
40366087c08SValeria Barra   CeedScalar *coords;
40466087c08SValeria Barra   CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
40566087c08SValeria Barra   if (dim == 1) {
40666087c08SValeria Barra     for (CeedInt i = 0; i < mesh_size; i++) {
40766087c08SValeria Barra       // map [0,1] to [0,1] varying the mesh density
40866087c08SValeria Barra       coords[i] = 0.5 + 1./sqrt(3.) * sin((2./3.) * M_PI*(coords[i] - 0.5));
40966087c08SValeria Barra     }
41066087c08SValeria Barra     exact_volume = 1.;
41166087c08SValeria Barra   } else {
41266087c08SValeria Barra     CeedInt num_nodes = mesh_size/dim;
41366087c08SValeria Barra     for (CeedInt i = 0; i < num_nodes; i++) {
41466087c08SValeria Barra       // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar
41566087c08SValeria Barra       // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi
41666087c08SValeria Barra       CeedScalar u = coords[i], v = coords[i+num_nodes];
41766087c08SValeria Barra       u = 1. + u;
41866087c08SValeria Barra       v = M_PI_2 * v;
41966087c08SValeria Barra       coords[i] = u * cos(v);
42066087c08SValeria Barra       coords[i+num_nodes] = u * sin(v);
42166087c08SValeria Barra     }
42266087c08SValeria Barra     exact_volume = 3./4. * M_PI;
42366087c08SValeria Barra   }
42466087c08SValeria Barra   CeedVectorRestoreArray(mesh_coords, &coords);
42566087c08SValeria Barra   return exact_volume;
42666087c08SValeria Barra }
427