xref: /libCEED/examples/ceed/ex1-volume.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
366087c08SValeria Barra //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
566087c08SValeria Barra //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
766087c08SValeria Barra 
866087c08SValeria Barra //                             libCEED Example 1
966087c08SValeria Barra //
1066087c08SValeria Barra // This example illustrates a simple usage of libCEED to compute the volume of a
1166087c08SValeria Barra // 3D body using matrix-free application of a mass operator.  Arbitrary mesh and
12ded9b81dSJeremy L Thompson // solution degrees in 1D, 2D and 3D are supported from the same code.
1366087c08SValeria Barra //
1466087c08SValeria Barra // The example has no dependencies, and is designed to be self-contained. For
1566087c08SValeria Barra // additional examples that use external discretization libraries (MFEM, PETSc,
1666087c08SValeria Barra // etc.) see the subdirectories in libceed/examples.
1766087c08SValeria Barra //
1866087c08SValeria Barra // All libCEED objects use a Ceed device object constructed based on a command
1966087c08SValeria Barra // line argument (-ceed).
2066087c08SValeria Barra //
2166087c08SValeria Barra // Build with:
2266087c08SValeria Barra //
2366087c08SValeria Barra //     make ex1-volume [CEED_DIR=</path/to/libceed>]
2466087c08SValeria Barra //
2566087c08SValeria Barra // Sample runs:
2666087c08SValeria Barra //
2766087c08SValeria Barra //     ./ex1-volume
2866087c08SValeria Barra //     ./ex1-volume -ceed /cpu/self
2928688798Sjeremylt //     ./ex1-volume -ceed /gpu/cuda
3066087c08SValeria Barra //
3166087c08SValeria Barra // Test in 1D-3D
32dc8efd83SLeila Ghaffari //TESTARGS(name="1D_User_QFunction") -ceed {ceed_resource} -d 1 -t
33dc8efd83SLeila Ghaffari //TESTARGS(name="2D_User_QFunction") -ceed {ceed_resource} -d 2 -t
34dc8efd83SLeila Ghaffari //TESTARGS(name="3D_User_QFunction") -ceed {ceed_resource} -d 3 -t
35dc8efd83SLeila Ghaffari //TESTARGS(name="1D_Gallery_QFunction") -ceed {ceed_resource} -d 1 -t -g
36dc8efd83SLeila Ghaffari //TESTARGS(name="2D_Gallery_QFunction") -ceed {ceed_resource} -d 2 -t -g
37dc8efd83SLeila Ghaffari //TESTARGS(name="3D_Gallery_QFunction") -ceed {ceed_resource} -d 3 -t -g
3866087c08SValeria Barra 
3966087c08SValeria Barra /// @file
4066087c08SValeria Barra /// libCEED example using mass operator to compute volume
4166087c08SValeria Barra 
42*2b730f8bSJeremy L Thompson #include "ex1-volume.h"
43*2b730f8bSJeremy L Thompson 
4466087c08SValeria Barra #include <ceed.h>
4566087c08SValeria Barra #include <math.h>
463d576824SJeremy L Thompson #include <stdlib.h>
4766087c08SValeria Barra #include <string.h>
4866087c08SValeria Barra 
49*2b730f8bSJeremy L Thompson // Auxiliary functions
50*2b730f8bSJeremy L Thompson int        GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]);
51*2b730f8bSJeremy L Thompson int        BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
52*2b730f8bSJeremy L Thompson                                      CeedElemRestriction *restr, CeedElemRestriction *restr_i);
53*2b730f8bSJeremy L Thompson int        SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords);
54*2b730f8bSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords);
5566087c08SValeria Barra 
56*2b730f8bSJeremy L Thompson // Main example
5766087c08SValeria Barra int main(int argc, const char *argv[]) {
5866087c08SValeria Barra   const char *ceed_spec   = "/cpu/self";
59990fdeb6SJeremy L Thompson   CeedInt     dim         = 3;               // dimension of the mesh
60990fdeb6SJeremy L Thompson   CeedInt     num_comp_x  = 3;               // number of x components
61990fdeb6SJeremy L Thompson   CeedInt     mesh_degree = 4;               // polynomial degree for the mesh
62990fdeb6SJeremy L Thompson   CeedInt     sol_degree  = 4;               // polynomial degree for the solution
63990fdeb6SJeremy L Thompson   CeedInt     num_qpts    = sol_degree + 2;  // number of 1D quadrature points
64990fdeb6SJeremy L Thompson   CeedInt     prob_size   = -1;              // approximate problem size
65990fdeb6SJeremy L Thompson   CeedInt     help = 0, test = 0, gallery = 0;
6666087c08SValeria Barra 
6766087c08SValeria Barra   // Process command line arguments.
6866087c08SValeria Barra   for (int ia = 1; ia < argc; ia++) {
69ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
7066087c08SValeria Barra     int next_arg = ((ia + 1) < argc), parse_error = 0;
7166087c08SValeria Barra     if (!strcmp(argv[ia], "-h")) {
7266087c08SValeria Barra       help = 1;
7366087c08SValeria Barra     } else if (!strcmp(argv[ia], "-c") || !strcmp(argv[ia], "-ceed")) {
7466087c08SValeria Barra       parse_error = next_arg ? ceed_spec = argv[++ia], 0 : 1;
7566087c08SValeria Barra     } else if (!strcmp(argv[ia], "-d")) {
7666087c08SValeria Barra       parse_error = next_arg ? dim = atoi(argv[++ia]), 0 : 1;
77d1d35e2fSjeremylt       num_comp_x                   = dim;
7866087c08SValeria Barra     } else if (!strcmp(argv[ia], "-m")) {
79ded9b81dSJeremy L Thompson       parse_error = next_arg ? mesh_degree = atoi(argv[++ia]), 0 : 1;
80ded9b81dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-p")) {
81ded9b81dSJeremy L Thompson       parse_error = next_arg ? sol_degree = atoi(argv[++ia]), 0 : 1;
8266087c08SValeria Barra     } else if (!strcmp(argv[ia], "-q")) {
8366087c08SValeria Barra       parse_error = next_arg ? num_qpts = atoi(argv[++ia]), 0 : 1;
8466087c08SValeria Barra     } else if (!strcmp(argv[ia], "-s")) {
8566087c08SValeria Barra       parse_error = next_arg ? prob_size = atoi(argv[++ia]), 0 : 1;
8666087c08SValeria Barra     } else if (!strcmp(argv[ia], "-t")) {
8766087c08SValeria Barra       test = 1;
8866087c08SValeria Barra     } else if (!strcmp(argv[ia], "-g")) {
8966087c08SValeria Barra       gallery = 1;
9066087c08SValeria Barra     }
9166087c08SValeria Barra     if (parse_error) {
9266087c08SValeria Barra       printf("Error parsing command line options.\n");
9366087c08SValeria Barra       return 1;
9466087c08SValeria Barra     }
95ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
9666087c08SValeria Barra   }
9766087c08SValeria Barra   if (prob_size < 0) prob_size = test ? 8 * 16 : 256 * 1024;
9866087c08SValeria Barra 
9966087c08SValeria Barra   // Print the values of all options:
10066087c08SValeria Barra   if (!test || help) {
101ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
10266087c08SValeria Barra     printf("Selected options: [command line option] : <current value>\n");
10366087c08SValeria Barra     printf("  Ceed specification [-c] : %s\n", ceed_spec);
104990fdeb6SJeremy L Thompson     printf("  Mesh dimension     [-d] : %" CeedInt_FMT "\n", dim);
105990fdeb6SJeremy L Thompson     printf("  Mesh degree        [-m] : %" CeedInt_FMT "\n", mesh_degree);
106990fdeb6SJeremy L Thompson     printf("  Solution degree    [-p] : %" CeedInt_FMT "\n", sol_degree);
107990fdeb6SJeremy L Thompson     printf("  Num. 1D quadr. pts [-q] : %" CeedInt_FMT "\n", num_qpts);
108990fdeb6SJeremy L Thompson     printf("  Approx. # unknowns [-s] : %" CeedInt_FMT "\n", prob_size);
10966087c08SValeria Barra     printf("  QFunction source   [-g] : %s\n", gallery ? "gallery" : "header");
11066087c08SValeria Barra     if (help) {
11166087c08SValeria Barra       printf("Test/quiet mode is %s\n", (test ? "ON" : "OFF (use -t to enable)"));
11266087c08SValeria Barra       return 0;
11366087c08SValeria Barra     }
11466087c08SValeria Barra     printf("\n");
115ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
11666087c08SValeria Barra   }
11766087c08SValeria Barra 
11866087c08SValeria Barra   // Select appropriate backend and logical device based on the <ceed-spec>
11966087c08SValeria Barra   // command line argument.
12066087c08SValeria Barra   Ceed ceed;
12166087c08SValeria Barra   CeedInit(ceed_spec, &ceed);
12266087c08SValeria Barra 
12366087c08SValeria Barra   // Construct the mesh and solution bases.
12466087c08SValeria Barra   CeedBasis mesh_basis, sol_basis;
125*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis);
126*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis);
12766087c08SValeria Barra 
12866087c08SValeria Barra   // Determine the mesh size based on the given approximate problem size.
129990fdeb6SJeremy L Thompson   CeedInt num_xyz[dim];
130d1d35e2fSjeremylt   GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz);
13166087c08SValeria Barra   if (!test) {
132ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
133990fdeb6SJeremy L Thompson     printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]);
134*2b730f8bSJeremy L Thompson     if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]);
135*2b730f8bSJeremy L Thompson     if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]);
13666087c08SValeria Barra     printf("\n");
137ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
13866087c08SValeria Barra   }
13966087c08SValeria Barra 
14066087c08SValeria Barra   // Build CeedElemRestriction objects describing the mesh and solution discrete
14166087c08SValeria Barra   // representations.
14266087c08SValeria Barra   CeedInt             mesh_size, sol_size;
14315910d16Sjeremylt   CeedElemRestriction mesh_restr, sol_restr, sol_restr_i;
144*2b730f8bSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restr, NULL);
145*2b730f8bSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restr, &sol_restr_i);
14666087c08SValeria Barra   if (!test) {
147e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
148990fdeb6SJeremy L Thompson     printf("Number of mesh nodes     : %" CeedInt_FMT "\n", mesh_size / dim);
149990fdeb6SJeremy L Thompson     printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size);
150e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
15166087c08SValeria Barra   }
15266087c08SValeria Barra 
15366087c08SValeria Barra   // Create a CeedVector with the mesh coordinates.
15466087c08SValeria Barra   CeedVector mesh_coords;
15566087c08SValeria Barra   CeedVectorCreate(ceed, mesh_size, &mesh_coords);
156d1d35e2fSjeremylt   SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords);
15766087c08SValeria Barra 
15866087c08SValeria Barra   // Apply a transformation to the mesh.
15966087c08SValeria Barra   CeedScalar exact_vol = TransformMeshCoords(dim, mesh_size, mesh_coords);
16066087c08SValeria Barra 
161ded9b81dSJeremy L Thompson   // Context data to be passed to the 'f_build_mass' QFunction.
162777ff853SJeremy L Thompson   CeedQFunctionContext build_ctx;
163777ff853SJeremy L Thompson   struct BuildContext  build_ctx_data;
164777ff853SJeremy L Thompson   build_ctx_data.dim = build_ctx_data.space_dim = dim;
165777ff853SJeremy L Thompson   CeedQFunctionContextCreate(ceed, &build_ctx);
166*2b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
16766087c08SValeria Barra 
168ded9b81dSJeremy L Thompson   // Create the QFunction that builds the mass operator (i.e. computes its
16966087c08SValeria Barra   // quadrature data) and set its context data.
170d1d35e2fSjeremylt   CeedQFunction qf_build;
17166087c08SValeria Barra   switch (gallery) {
17266087c08SValeria Barra     case 0:
17366087c08SValeria Barra       // This creates the QFunction directly.
174*2b730f8bSJeremy L Thompson       CeedQFunctionCreateInterior(ceed, 1, f_build_mass, f_build_mass_loc, &qf_build);
175d1d35e2fSjeremylt       CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
176d1d35e2fSjeremylt       CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
177d1d35e2fSjeremylt       CeedQFunctionAddOutput(qf_build, "qdata", 1, CEED_EVAL_NONE);
178d1d35e2fSjeremylt       CeedQFunctionSetContext(qf_build, build_ctx);
17966087c08SValeria Barra       break;
18066087c08SValeria Barra     case 1: {
18166087c08SValeria Barra       // This creates the QFunction via the gallery.
18266087c08SValeria Barra       char name[13] = "";
183990fdeb6SJeremy L Thompson       snprintf(name, sizeof name, "Mass%" CeedInt_FMT "DBuild", dim);
184d1d35e2fSjeremylt       CeedQFunctionCreateInteriorByName(ceed, name, &qf_build);
18566087c08SValeria Barra       break;
18666087c08SValeria Barra     }
18766087c08SValeria Barra   }
18866087c08SValeria Barra 
18966087c08SValeria Barra   // Create the operator that builds the quadrature data for the mass operator.
190d1d35e2fSjeremylt   CeedOperator op_build;
191*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
192*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_build, "dx", mesh_restr, mesh_basis, CEED_VECTOR_ACTIVE);
193*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE);
194*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_build, "qdata", sol_restr_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
19566087c08SValeria Barra 
19666087c08SValeria Barra   // Compute the quadrature data for the mass operator.
197d1d35e2fSjeremylt   CeedVector q_data;
19866087c08SValeria Barra   CeedInt    elem_qpts = CeedIntPow(num_qpts, dim);
19966087c08SValeria Barra   CeedInt    num_elem  = 1;
200*2b730f8bSJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d];
201d1d35e2fSjeremylt   CeedVectorCreate(ceed, num_elem * elem_qpts, &q_data);
202*2b730f8bSJeremy L Thompson   CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE);
20366087c08SValeria Barra 
204ded9b81dSJeremy L Thompson   // Create the QFunction that defines the action of the mass operator.
205d1d35e2fSjeremylt   CeedQFunction qf_apply;
20666087c08SValeria Barra   switch (gallery) {
20766087c08SValeria Barra     case 0:
20866087c08SValeria Barra       // This creates the QFunction directly.
209*2b730f8bSJeremy L Thompson       CeedQFunctionCreateInterior(ceed, 1, f_apply_mass, f_apply_mass_loc, &qf_apply);
210d1d35e2fSjeremylt       CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP);
211d1d35e2fSjeremylt       CeedQFunctionAddInput(qf_apply, "qdata", 1, CEED_EVAL_NONE);
212d1d35e2fSjeremylt       CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP);
21366087c08SValeria Barra       break;
21466087c08SValeria Barra     case 1:
21566087c08SValeria Barra       // This creates the QFunction via the gallery.
216d1d35e2fSjeremylt       CeedQFunctionCreateInteriorByName(ceed, "MassApply", &qf_apply);
21766087c08SValeria Barra       break;
21866087c08SValeria Barra   }
21966087c08SValeria Barra 
22066087c08SValeria Barra   // Create the mass operator.
221d1d35e2fSjeremylt   CeedOperator op_apply;
222*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
223d1d35e2fSjeremylt   CeedOperatorSetField(op_apply, "u", sol_restr, sol_basis, CEED_VECTOR_ACTIVE);
224*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply, "qdata", sol_restr_i, CEED_BASIS_COLLOCATED, q_data);
225d1d35e2fSjeremylt   CeedOperatorSetField(op_apply, "v", sol_restr, sol_basis, CEED_VECTOR_ACTIVE);
22666087c08SValeria Barra 
22766087c08SValeria Barra   // Create auxiliary solution-size vectors.
22866087c08SValeria Barra   CeedVector u, v;
22966087c08SValeria Barra   CeedVectorCreate(ceed, sol_size, &u);
23066087c08SValeria Barra   CeedVectorCreate(ceed, sol_size, &v);
23166087c08SValeria Barra 
23266087c08SValeria Barra   // Initialize 'u' and 'v' with ones.
23366087c08SValeria Barra   CeedVectorSetValue(u, 1.0);
23466087c08SValeria Barra 
235ded9b81dSJeremy L Thompson   // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1
236d1d35e2fSjeremylt   CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
23766087c08SValeria Barra 
23866087c08SValeria Barra   // Compute and print the sum of the entries of 'v' giving the mesh volume.
239d1d35e2fSjeremylt   const CeedScalar *v_array;
240d1d35e2fSjeremylt   CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
24166087c08SValeria Barra   CeedScalar vol = 0.;
242*2b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < sol_size; i++) vol += v_array[i];
243d1d35e2fSjeremylt   CeedVectorRestoreArrayRead(v, &v_array);
24466087c08SValeria Barra   if (!test) {
245ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
24666087c08SValeria Barra     printf(" done.\n");
24766087c08SValeria Barra     printf("Exact mesh volume    : % .14g\n", exact_vol);
24866087c08SValeria Barra     printf("Computed mesh volume : % .14g\n", vol);
24966087c08SValeria Barra     printf("Volume error         : % .14g\n", vol - exact_vol);
250ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
25166087c08SValeria Barra   } else {
25280a9ef05SNatalie Beams     CeedScalar tol = (dim == 1 ? 100. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5);
253*2b730f8bSJeremy L Thompson     if (fabs(vol - exact_vol) > tol) printf("Volume error : % .1e\n", vol - exact_vol);
25466087c08SValeria Barra   }
25566087c08SValeria Barra 
25666087c08SValeria Barra   // Free dynamically allocated memory.
25766087c08SValeria Barra   CeedVectorDestroy(&u);
25866087c08SValeria Barra   CeedVectorDestroy(&v);
259d1d35e2fSjeremylt   CeedVectorDestroy(&q_data);
26066087c08SValeria Barra   CeedVectorDestroy(&mesh_coords);
261d1d35e2fSjeremylt   CeedOperatorDestroy(&op_apply);
262d1d35e2fSjeremylt   CeedQFunctionDestroy(&qf_apply);
263777ff853SJeremy L Thompson   CeedQFunctionContextDestroy(&build_ctx);
264d1d35e2fSjeremylt   CeedOperatorDestroy(&op_build);
265d1d35e2fSjeremylt   CeedQFunctionDestroy(&qf_build);
26666087c08SValeria Barra   CeedElemRestrictionDestroy(&sol_restr);
26766087c08SValeria Barra   CeedElemRestrictionDestroy(&mesh_restr);
26866087c08SValeria Barra   CeedElemRestrictionDestroy(&sol_restr_i);
26966087c08SValeria Barra   CeedBasisDestroy(&sol_basis);
27066087c08SValeria Barra   CeedBasisDestroy(&mesh_basis);
27166087c08SValeria Barra   CeedDestroy(&ceed);
27266087c08SValeria Barra   return 0;
27366087c08SValeria Barra }
27466087c08SValeria Barra 
275*2b730f8bSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) {
27666087c08SValeria Barra   // Use the approximate formula:
277ded9b81dSJeremy L Thompson   //    prob_size ~ num_elem * degree^dim
278ded9b81dSJeremy L Thompson   CeedInt num_elem = prob_size / CeedIntPow(degree, dim);
27966087c08SValeria Barra   CeedInt s        = 0;  // find s: num_elem/2 < 2^s <= num_elem
28066087c08SValeria Barra   while (num_elem > 1) {
28166087c08SValeria Barra     num_elem /= 2;
28266087c08SValeria Barra     s++;
28366087c08SValeria Barra   }
28466087c08SValeria Barra   CeedInt r = s % dim;
285990fdeb6SJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
286990fdeb6SJeremy L Thompson     CeedInt sd = s / dim;
287*2b730f8bSJeremy L Thompson     if (r > 0) {
288*2b730f8bSJeremy L Thompson       sd++;
289*2b730f8bSJeremy L Thompson       r--;
290*2b730f8bSJeremy L Thompson     }
291d1d35e2fSjeremylt     num_xyz[d] = 1 << sd;
29266087c08SValeria Barra   }
29366087c08SValeria Barra   return 0;
29466087c08SValeria Barra }
29566087c08SValeria Barra 
296*2b730f8bSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
297*2b730f8bSJeremy L Thompson                               CeedElemRestriction *restr, CeedElemRestriction *restr_i) {
298ded9b81dSJeremy L Thompson   CeedInt p         = degree + 1;
299d1d35e2fSjeremylt   CeedInt num_nodes = CeedIntPow(p, dim);         // number of scalar nodes per element
30066087c08SValeria Barra   CeedInt elem_qpts = CeedIntPow(num_qpts, dim);  // number of qpts per element
30166087c08SValeria Barra   CeedInt nd[3], num_elem = 1, scalar_size = 1;
302990fdeb6SJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
303d1d35e2fSjeremylt     num_elem *= num_xyz[d];
304d1d35e2fSjeremylt     nd[d] = num_xyz[d] * (p - 1) + 1;
30566087c08SValeria Barra     scalar_size *= nd[d];
30666087c08SValeria Barra   }
307d1d35e2fSjeremylt   *size = scalar_size * num_comp;
30866087c08SValeria Barra   // elem:         0             1                 n-1
30966087c08SValeria Barra   //           |---*-...-*---|---*-...-*---|- ... -|--...--|
310d1d35e2fSjeremylt   // num_nodes:   0   1    p-1  p  p+1       2*p             n*p
311d1d35e2fSjeremylt   CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes);
31266087c08SValeria Barra   for (CeedInt e = 0; e < num_elem; e++) {
313d1d35e2fSjeremylt     CeedInt e_xyz[3] = {1, 1, 1}, re = e;
314*2b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
315*2b730f8bSJeremy L Thompson       e_xyz[d] = re % num_xyz[d];
316*2b730f8bSJeremy L Thompson       re /= num_xyz[d];
317*2b730f8bSJeremy L Thompson     }
318d1d35e2fSjeremylt     CeedInt *loc_el_nodes = elem_nodes + e * num_nodes;
319990fdeb6SJeremy L Thompson     for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) {
320d1d35e2fSjeremylt       CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes;
321990fdeb6SJeremy L Thompson       for (CeedInt d = 0; d < dim; d++) {
322d1d35e2fSjeremylt         g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride;
323d1d35e2fSjeremylt         g_nodes_stride *= nd[d];
324d1d35e2fSjeremylt         r_nodes /= p;
32566087c08SValeria Barra       }
326d1d35e2fSjeremylt       loc_el_nodes[l_nodes] = g_nodes;
32766087c08SValeria Barra     }
32866087c08SValeria Barra   }
329*2b730f8bSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes,
330*2b730f8bSJeremy L Thompson                             restr);
331*2b730f8bSJeremy L Thompson   if (restr_i) CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, restr_i);
332d1d35e2fSjeremylt   free(elem_nodes);
33366087c08SValeria Barra   return 0;
33466087c08SValeria Barra }
33566087c08SValeria Barra 
336*2b730f8bSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) {
337ded9b81dSJeremy L Thompson   CeedInt p = mesh_degree + 1;
33866087c08SValeria Barra   CeedInt nd[3], num_elem = 1, scalar_size = 1;
339990fdeb6SJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
340d1d35e2fSjeremylt     num_elem *= num_xyz[d];
341d1d35e2fSjeremylt     nd[d] = num_xyz[d] * (p - 1) + 1;
34266087c08SValeria Barra     scalar_size *= nd[d];
34366087c08SValeria Barra   }
34466087c08SValeria Barra   CeedScalar *coords;
3459c774eddSJeremy L Thompson   CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords);
346ded9b81dSJeremy L Thompson   CeedScalar *nodes = malloc(sizeof(CeedScalar) * p);
34766087c08SValeria Barra   // The H1 basis uses Lobatto quadrature points as nodes.
348ded9b81dSJeremy L Thompson   CeedLobattoQuadrature(p, nodes, NULL);  // nodes are in [-1,1]
349*2b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < p; i++) {
350*2b730f8bSJeremy L Thompson     nodes[i] = 0.5 + 0.5 * nodes[i];
351*2b730f8bSJeremy L Thompson   }
352d1d35e2fSjeremylt   for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) {
353d1d35e2fSjeremylt     CeedInt r_nodes = gs_nodes;
354990fdeb6SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
355d1d35e2fSjeremylt       CeedInt d_1d                       = r_nodes % nd[d];
356*2b730f8bSJeremy L Thompson       coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d];
357d1d35e2fSjeremylt       r_nodes /= nd[d];
35866087c08SValeria Barra     }
35966087c08SValeria Barra   }
36066087c08SValeria Barra   free(nodes);
36166087c08SValeria Barra   CeedVectorRestoreArray(mesh_coords, &coords);
36266087c08SValeria Barra   return 0;
36366087c08SValeria Barra }
36466087c08SValeria Barra 
36566087c08SValeria Barra #ifndef M_PI
36666087c08SValeria Barra #define M_PI 3.14159265358979323846
36766087c08SValeria Barra #define M_PI_2 1.57079632679489661923
36866087c08SValeria Barra #endif
36966087c08SValeria Barra 
370*2b730f8bSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) {
37166087c08SValeria Barra   CeedScalar  exact_volume;
37266087c08SValeria Barra   CeedScalar *coords;
37366087c08SValeria Barra   CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
37466087c08SValeria Barra   if (dim == 1) {
37566087c08SValeria Barra     for (CeedInt i = 0; i < mesh_size; i++) {
37666087c08SValeria Barra       // map [0,1] to [0,1] varying the mesh density
37766087c08SValeria Barra       coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5));
37866087c08SValeria Barra     }
37966087c08SValeria Barra     exact_volume = 1.;
38066087c08SValeria Barra   } else {
38166087c08SValeria Barra     CeedInt num_nodes = mesh_size / dim;
38266087c08SValeria Barra     for (CeedInt i = 0; i < num_nodes; i++) {
38366087c08SValeria Barra       // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar
38466087c08SValeria Barra       // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi
38566087c08SValeria Barra       CeedScalar u = coords[i], v = coords[i + num_nodes];
38666087c08SValeria Barra       u                     = 1. + u;
38766087c08SValeria Barra       v                     = M_PI_2 * v;
38866087c08SValeria Barra       coords[i]             = u * cos(v);
38966087c08SValeria Barra       coords[i + num_nodes] = u * sin(v);
39066087c08SValeria Barra     }
39166087c08SValeria Barra     exact_volume = 3. / 4. * M_PI;
39266087c08SValeria Barra   }
39366087c08SValeria Barra   CeedVectorRestoreArray(mesh_coords, &coords);
39466087c08SValeria Barra   return exact_volume;
39566087c08SValeria Barra }
396