xref: /libCEED/examples/ceed/ex1-volume.c (revision 82138112808ac45c6722ef2bfe52ea5cd96df80f)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 //
10ea61e9acSJeremy L Thompson // This example illustrates a simple usage of libCEED to compute the volume of a 3D body using matrix-free application of a mass operator.
11ea61e9acSJeremy L Thompson // Arbitrary mesh and solution degrees in 1D, 2D and 3D are supported from the same code.
1266087c08SValeria Barra //
13ea61e9acSJeremy L Thompson // The example has no dependencies, and is designed to be self-contained.
14ea61e9acSJeremy L Thompson // For additional examples that use external discretization libraries (MFEM, PETSc, etc.) see the subdirectories in libceed/examples.
1566087c08SValeria Barra //
16ea61e9acSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command line argument (-ceed).
1766087c08SValeria Barra //
1866087c08SValeria Barra // Build with:
1966087c08SValeria Barra //
2066087c08SValeria Barra //     make ex1-volume [CEED_DIR=</path/to/libceed>]
2166087c08SValeria Barra //
2266087c08SValeria Barra // Sample runs:
2366087c08SValeria Barra //
2466087c08SValeria Barra //     ./ex1-volume
2566087c08SValeria Barra //     ./ex1-volume -ceed /cpu/self
2628688798Sjeremylt //     ./ex1-volume -ceed /gpu/cuda
2766087c08SValeria Barra //
2866087c08SValeria Barra // Test in 1D-3D
29fa5adaf5SJeremy L Thompson //TESTARGS(name="1D User QFunction") -ceed {ceed_resource} -d 1 -t
30fa5adaf5SJeremy L Thompson //TESTARGS(name="2D User QFunction") -ceed {ceed_resource} -d 2 -t
31fa5adaf5SJeremy L Thompson //TESTARGS(name="3D User QFunction") -ceed {ceed_resource} -d 3 -t
32fa5adaf5SJeremy L Thompson //TESTARGS(name="1D Gallery QFunction") -ceed {ceed_resource} -d 1 -t -g
33fa5adaf5SJeremy L Thompson //TESTARGS(name="2D Gallery QFunction") -ceed {ceed_resource} -d 2 -t -g
34fa5adaf5SJeremy L Thompson //TESTARGS(name="3D Gallery QFunction") -ceed {ceed_resource} -d 3 -t -g
3566087c08SValeria Barra 
3666087c08SValeria Barra /// @file
3766087c08SValeria Barra /// libCEED example using mass operator to compute volume
3866087c08SValeria Barra 
392b730f8bSJeremy L Thompson #include "ex1-volume.h"
402b730f8bSJeremy L Thompson 
4166087c08SValeria Barra #include <ceed.h>
4266087c08SValeria Barra #include <math.h>
4349aac155SJeremy L Thompson #include <stdio.h>
443d576824SJeremy L Thompson #include <stdlib.h>
4566087c08SValeria Barra #include <string.h>
4666087c08SValeria Barra 
472b730f8bSJeremy L Thompson // Auxiliary functions
482b730f8bSJeremy L Thompson int        GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]);
492b730f8bSJeremy L Thompson int        BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
50d37d859eSJeremy L Thompson                                      CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction);
512b730f8bSJeremy L Thompson int        SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords);
522b730f8bSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords);
5366087c08SValeria Barra 
542b730f8bSJeremy L Thompson // Main example
5566087c08SValeria Barra int main(int argc, const char *argv[]) {
5666087c08SValeria Barra   const char *ceed_spec   = "/cpu/self";
57990fdeb6SJeremy L Thompson   CeedInt     dim         = 3;               // dimension of the mesh
58990fdeb6SJeremy L Thompson   CeedInt     num_comp_x  = 3;               // number of x components
59990fdeb6SJeremy L Thompson   CeedInt     mesh_degree = 4;               // polynomial degree for the mesh
60990fdeb6SJeremy L Thompson   CeedInt     sol_degree  = 4;               // polynomial degree for the solution
61990fdeb6SJeremy L Thompson   CeedInt     num_qpts    = sol_degree + 2;  // number of 1D quadrature points
62990fdeb6SJeremy L Thompson   CeedInt     prob_size   = -1;              // approximate problem size
63*82138112SJeremy L Thompson   CeedInt     help = 0, test = 0, gallery = 0, benchmark = 0;
6466087c08SValeria Barra 
6566087c08SValeria Barra   // Process command line arguments.
6666087c08SValeria Barra   for (int ia = 1; ia < argc; ia++) {
67ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
6866087c08SValeria Barra     int next_arg = ((ia + 1) < argc), parse_error = 0;
6966087c08SValeria Barra     if (!strcmp(argv[ia], "-h")) {
7066087c08SValeria Barra       help = 1;
7166087c08SValeria Barra     } else if (!strcmp(argv[ia], "-c") || !strcmp(argv[ia], "-ceed")) {
7266087c08SValeria Barra       parse_error = next_arg ? ceed_spec = argv[++ia], 0 : 1;
7366087c08SValeria Barra     } else if (!strcmp(argv[ia], "-d")) {
7466087c08SValeria Barra       parse_error = next_arg ? dim = atoi(argv[++ia]), 0 : 1;
75d1d35e2fSjeremylt       num_comp_x                   = dim;
7666087c08SValeria Barra     } else if (!strcmp(argv[ia], "-m")) {
77ded9b81dSJeremy L Thompson       parse_error = next_arg ? mesh_degree = atoi(argv[++ia]), 0 : 1;
78ded9b81dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-p")) {
79ded9b81dSJeremy L Thompson       parse_error = next_arg ? sol_degree = atoi(argv[++ia]), 0 : 1;
8066087c08SValeria Barra     } else if (!strcmp(argv[ia], "-q")) {
8166087c08SValeria Barra       parse_error = next_arg ? num_qpts = atoi(argv[++ia]), 0 : 1;
8266087c08SValeria Barra     } else if (!strcmp(argv[ia], "-s")) {
8366087c08SValeria Barra       parse_error = next_arg ? prob_size = atoi(argv[++ia]), 0 : 1;
84*82138112SJeremy L Thompson     } else if (!strcmp(argv[ia], "-b")) {
85*82138112SJeremy L Thompson       parse_error = next_arg ? benchmark = 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);
107d37d859eSJeremy L Thompson     printf("  Num. 1D quadrature 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 
118ea61e9acSJeremy L Thompson   // Select appropriate backend and logical device based on the (-ceed) command line argument.
11966087c08SValeria Barra   Ceed ceed;
12066087c08SValeria Barra   CeedInit(ceed_spec, &ceed);
12166087c08SValeria Barra 
12266087c08SValeria Barra   // Construct the mesh and solution bases.
12366087c08SValeria Barra   CeedBasis mesh_basis, sol_basis;
1242b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis);
1252b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis);
12666087c08SValeria Barra 
12766087c08SValeria Barra   // Determine the mesh size based on the given approximate problem size.
128990fdeb6SJeremy L Thompson   CeedInt num_xyz[dim];
129d1d35e2fSjeremylt   GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz);
13066087c08SValeria Barra   if (!test) {
131ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
132990fdeb6SJeremy L Thompson     printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]);
1332b730f8bSJeremy L Thompson     if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]);
1342b730f8bSJeremy L Thompson     if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]);
13566087c08SValeria Barra     printf("\n");
136ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
13766087c08SValeria Barra   }
13866087c08SValeria Barra 
139ea61e9acSJeremy L Thompson   // Build CeedElemRestriction objects describing the mesh and solution discrete representations.
14066087c08SValeria Barra   CeedInt             mesh_size, sol_size;
141d37d859eSJeremy L Thompson   CeedElemRestriction mesh_restriction, sol_restriction, q_data_restriction;
142d37d859eSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restriction, NULL);
143d37d859eSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restriction, &q_data_restriction);
14466087c08SValeria Barra   if (!test) {
145e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
146990fdeb6SJeremy L Thompson     printf("Number of mesh nodes     : %" CeedInt_FMT "\n", mesh_size / dim);
147990fdeb6SJeremy L Thompson     printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size);
148e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
14966087c08SValeria Barra   }
15066087c08SValeria Barra 
15166087c08SValeria Barra   // Create a CeedVector with the mesh coordinates.
15266087c08SValeria Barra   CeedVector mesh_coords;
15366087c08SValeria Barra   CeedVectorCreate(ceed, mesh_size, &mesh_coords);
154d1d35e2fSjeremylt   SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords);
15566087c08SValeria Barra 
15666087c08SValeria Barra   // Apply a transformation to the mesh.
157d37d859eSJeremy L Thompson   CeedScalar exact_volume = TransformMeshCoords(dim, mesh_size, mesh_coords);
15866087c08SValeria Barra 
159d37d859eSJeremy L Thompson   // Context data to be passed to the 'build_mass' QFunction.
160777ff853SJeremy L Thompson   CeedQFunctionContext build_ctx;
161777ff853SJeremy L Thompson   struct BuildContext  build_ctx_data;
162777ff853SJeremy L Thompson   build_ctx_data.dim = build_ctx_data.space_dim = dim;
163777ff853SJeremy L Thompson   CeedQFunctionContextCreate(ceed, &build_ctx);
1642b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
16566087c08SValeria Barra 
166ea61e9acSJeremy L Thompson   // Create the QFunction that builds the mass operator (i.e. computes its quadrature data) and set its context data.
167d1d35e2fSjeremylt   CeedQFunction qf_build;
168d37d859eSJeremy L Thompson   if (gallery) {
16966087c08SValeria Barra     // This creates the QFunction via the gallery.
17066087c08SValeria Barra     char name[13] = "";
171990fdeb6SJeremy L Thompson     snprintf(name, sizeof name, "Mass%" CeedInt_FMT "DBuild", dim);
172d1d35e2fSjeremylt     CeedQFunctionCreateInteriorByName(ceed, name, &qf_build);
173d37d859eSJeremy L Thompson   } else {
174d37d859eSJeremy L Thompson     // This creates the QFunction directly.
175d37d859eSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, build_mass, build_mass_loc, &qf_build);
176d37d859eSJeremy L Thompson     CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
177d37d859eSJeremy L Thompson     CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
178d37d859eSJeremy L Thompson     CeedQFunctionAddOutput(qf_build, "qdata", 1, CEED_EVAL_NONE);
179d37d859eSJeremy L Thompson     CeedQFunctionSetContext(qf_build, build_ctx);
18066087c08SValeria Barra   }
18166087c08SValeria Barra 
18266087c08SValeria Barra   // Create the operator that builds the quadrature data for the mass operator.
183d1d35e2fSjeremylt   CeedOperator op_build;
1842b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
185d37d859eSJeremy L Thompson   CeedOperatorSetField(op_build, "dx", mesh_restriction, mesh_basis, CEED_VECTOR_ACTIVE);
1862b730f8bSJeremy L Thompson   CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE);
187356036faSJeremy L Thompson   CeedOperatorSetField(op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
18866087c08SValeria Barra 
18966087c08SValeria Barra   // Compute the quadrature data for the mass operator.
190d1d35e2fSjeremylt   CeedVector q_data;
19166087c08SValeria Barra   CeedInt    elem_qpts = CeedIntPow(num_qpts, dim);
19266087c08SValeria Barra   CeedInt    num_elem  = 1;
1932b730f8bSJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d];
194d1d35e2fSjeremylt   CeedVectorCreate(ceed, num_elem * elem_qpts, &q_data);
1952b730f8bSJeremy L Thompson   CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE);
19666087c08SValeria Barra 
197ded9b81dSJeremy L Thompson   // Create the QFunction that defines the action of the mass operator.
198d1d35e2fSjeremylt   CeedQFunction qf_apply;
199d37d859eSJeremy L Thompson   if (gallery) {
200d37d859eSJeremy L Thompson     // This creates the QFunction via the gallery.
201d37d859eSJeremy L Thompson     CeedQFunctionCreateInteriorByName(ceed, "MassApply", &qf_apply);
202d37d859eSJeremy L Thompson   } else {
20366087c08SValeria Barra     // This creates the QFunction directly.
204d37d859eSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, apply_mass, apply_mass_loc, &qf_apply);
205d1d35e2fSjeremylt     CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP);
206d1d35e2fSjeremylt     CeedQFunctionAddInput(qf_apply, "qdata", 1, CEED_EVAL_NONE);
207d1d35e2fSjeremylt     CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP);
20866087c08SValeria Barra   }
20966087c08SValeria Barra 
21066087c08SValeria Barra   // Create the mass operator.
211d1d35e2fSjeremylt   CeedOperator op_apply;
2122b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
213d37d859eSJeremy L Thompson   CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
214356036faSJeremy L Thompson   CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data);
215d37d859eSJeremy L Thompson   CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
21666087c08SValeria Barra 
21766087c08SValeria Barra   // Create auxiliary solution-size vectors.
21866087c08SValeria Barra   CeedVector u, v;
21966087c08SValeria Barra   CeedVectorCreate(ceed, sol_size, &u);
22066087c08SValeria Barra   CeedVectorCreate(ceed, sol_size, &v);
22166087c08SValeria Barra 
222b8bc0c60SPeter Munch   // Initialize 'u' with ones.
22366087c08SValeria Barra   CeedVectorSetValue(u, 1.0);
22466087c08SValeria Barra 
225d37d859eSJeremy L Thompson   // Compute the mesh volume using the mass operator: volume = 1^T \cdot M \cdot 1
226d1d35e2fSjeremylt   CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
22766087c08SValeria Barra 
228*82138112SJeremy L Thompson   // Benchmark runs
229*82138112SJeremy L Thompson   if (!test && benchmark) {
230*82138112SJeremy L Thompson     // LCOV_EXCL_START
231*82138112SJeremy L Thompson     printf(" Executing %d benchmarking runs...\n", benchmark);
232*82138112SJeremy L Thompson     // LCOV_EXCL_STOP
233*82138112SJeremy L Thompson   }
234*82138112SJeremy L Thompson   for (CeedInt i = 0; i < benchmark; i++) {
235*82138112SJeremy L Thompson     // LCOV_EXCL_START
236*82138112SJeremy L Thompson     CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
237*82138112SJeremy L Thompson     // LCOV_EXCL_STOP
238*82138112SJeremy L Thompson   }
239*82138112SJeremy L Thompson 
24066087c08SValeria Barra   // Compute and print the sum of the entries of 'v' giving the mesh volume.
241d37d859eSJeremy L Thompson   CeedScalar volume = 0.;
242d37d859eSJeremy L Thompson   {
243d1d35e2fSjeremylt     const CeedScalar *v_array;
244d1d35e2fSjeremylt     CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
245d37d859eSJeremy L Thompson     for (CeedInt i = 0; i < sol_size; i++) volume += v_array[i];
246d1d35e2fSjeremylt     CeedVectorRestoreArrayRead(v, &v_array);
247d37d859eSJeremy L Thompson   }
24866087c08SValeria Barra   if (!test) {
249ded9b81dSJeremy L Thompson     // LCOV_EXCL_START
25066087c08SValeria Barra     printf(" done.\n");
251d37d859eSJeremy L Thompson     printf("Exact mesh volume    : % .14g\n", exact_volume);
252d37d859eSJeremy L Thompson     printf("Computed mesh volume : % .14g\n", volume);
253d37d859eSJeremy L Thompson     printf("Volume error         : % .14g\n", volume - exact_volume);
254ded9b81dSJeremy L Thompson     // LCOV_EXCL_STOP
25566087c08SValeria Barra   } else {
256bd882c8aSJames Wright     CeedScalar tol = (dim == 1 ? 200. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5);
257d37d859eSJeremy L Thompson     if (fabs(volume - exact_volume) > tol) printf("Volume error : % .1e\n", volume - exact_volume);
25866087c08SValeria Barra   }
25966087c08SValeria Barra 
26066087c08SValeria Barra   // Free dynamically allocated memory.
26166087c08SValeria Barra   CeedVectorDestroy(&u);
26266087c08SValeria Barra   CeedVectorDestroy(&v);
263d1d35e2fSjeremylt   CeedVectorDestroy(&q_data);
26466087c08SValeria Barra   CeedVectorDestroy(&mesh_coords);
265d1d35e2fSjeremylt   CeedOperatorDestroy(&op_apply);
266d1d35e2fSjeremylt   CeedQFunctionDestroy(&qf_apply);
267777ff853SJeremy L Thompson   CeedQFunctionContextDestroy(&build_ctx);
268d1d35e2fSjeremylt   CeedOperatorDestroy(&op_build);
269d1d35e2fSjeremylt   CeedQFunctionDestroy(&qf_build);
270d37d859eSJeremy L Thompson   CeedElemRestrictionDestroy(&sol_restriction);
271d37d859eSJeremy L Thompson   CeedElemRestrictionDestroy(&mesh_restriction);
272d37d859eSJeremy L Thompson   CeedElemRestrictionDestroy(&q_data_restriction);
27366087c08SValeria Barra   CeedBasisDestroy(&sol_basis);
27466087c08SValeria Barra   CeedBasisDestroy(&mesh_basis);
27566087c08SValeria Barra   CeedDestroy(&ceed);
27666087c08SValeria Barra   return 0;
27766087c08SValeria Barra }
27866087c08SValeria Barra 
2792b730f8bSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) {
28066087c08SValeria Barra   // Use the approximate formula:
281ded9b81dSJeremy L Thompson   //    prob_size ~ num_elem * degree^dim
282ded9b81dSJeremy L Thompson   CeedInt num_elem = prob_size / CeedIntPow(degree, dim);
28366087c08SValeria Barra   CeedInt s        = 0;  // find s: num_elem/2 < 2^s <= num_elem
28466087c08SValeria Barra   while (num_elem > 1) {
28566087c08SValeria Barra     num_elem /= 2;
28666087c08SValeria Barra     s++;
28766087c08SValeria Barra   }
28866087c08SValeria Barra   CeedInt r = s % dim;
289990fdeb6SJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
290990fdeb6SJeremy L Thompson     CeedInt sd = s / dim;
2912b730f8bSJeremy L Thompson     if (r > 0) {
2922b730f8bSJeremy L Thompson       sd++;
2932b730f8bSJeremy L Thompson       r--;
2942b730f8bSJeremy L Thompson     }
295d1d35e2fSjeremylt     num_xyz[d] = 1 << sd;
29666087c08SValeria Barra   }
29766087c08SValeria Barra   return 0;
29866087c08SValeria Barra }
29966087c08SValeria Barra 
3002b730f8bSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
301d37d859eSJeremy L Thompson                               CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction) {
302ded9b81dSJeremy L Thompson   CeedInt p         = degree + 1;
303d1d35e2fSjeremylt   CeedInt num_nodes = CeedIntPow(p, dim);         // number of scalar nodes per element
30466087c08SValeria Barra   CeedInt elem_qpts = CeedIntPow(num_qpts, dim);  // number of qpts per element
30566087c08SValeria Barra   CeedInt nd[3], num_elem = 1, scalar_size = 1;
306990fdeb6SJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
307d1d35e2fSjeremylt     num_elem *= num_xyz[d];
308d1d35e2fSjeremylt     nd[d] = num_xyz[d] * (p - 1) + 1;
30966087c08SValeria Barra     scalar_size *= nd[d];
31066087c08SValeria Barra   }
311d1d35e2fSjeremylt   *size = scalar_size * num_comp;
31266087c08SValeria Barra   // elem:         0             1                 n-1
31366087c08SValeria Barra   //           |---*-...-*---|---*-...-*---|- ... -|--...--|
314d1d35e2fSjeremylt   // num_nodes:   0   1    p-1  p  p+1       2*p             n*p
315d1d35e2fSjeremylt   CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes);
31666087c08SValeria Barra   for (CeedInt e = 0; e < num_elem; e++) {
317d1d35e2fSjeremylt     CeedInt e_xyz[3] = {1, 1, 1}, re = e;
3182b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
3192b730f8bSJeremy L Thompson       e_xyz[d] = re % num_xyz[d];
3202b730f8bSJeremy L Thompson       re /= num_xyz[d];
3212b730f8bSJeremy L Thompson     }
322d37d859eSJeremy L Thompson     CeedInt *local_elem_nodes = elem_nodes + e * num_nodes;
323990fdeb6SJeremy L Thompson     for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) {
324d1d35e2fSjeremylt       CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes;
325990fdeb6SJeremy L Thompson       for (CeedInt d = 0; d < dim; d++) {
326d1d35e2fSjeremylt         g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride;
327d1d35e2fSjeremylt         g_nodes_stride *= nd[d];
328d1d35e2fSjeremylt         r_nodes /= p;
32966087c08SValeria Barra       }
330d37d859eSJeremy L Thompson       local_elem_nodes[l_nodes] = g_nodes;
33166087c08SValeria Barra     }
33266087c08SValeria Barra   }
3332b730f8bSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes,
334d37d859eSJeremy L Thompson                             restriction);
3356c10af5dSJeremy L Thompson   if (q_data_restriction) {
336d37d859eSJeremy L Thompson     CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, q_data_restriction);
3376c10af5dSJeremy L Thompson   }
338d1d35e2fSjeremylt   free(elem_nodes);
33966087c08SValeria Barra   return 0;
34066087c08SValeria Barra }
34166087c08SValeria Barra 
3422b730f8bSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) {
343ded9b81dSJeremy L Thompson   CeedInt p = mesh_degree + 1;
3441d2f0dcbSJeremy L Thompson   CeedInt nd[3], scalar_size = 1;
345990fdeb6SJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
346d1d35e2fSjeremylt     nd[d] = num_xyz[d] * (p - 1) + 1;
34766087c08SValeria Barra     scalar_size *= nd[d];
34866087c08SValeria Barra   }
34966087c08SValeria Barra   CeedScalar *coords;
3509c774eddSJeremy L Thompson   CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords);
351ded9b81dSJeremy L Thompson   CeedScalar *nodes = malloc(sizeof(CeedScalar) * p);
35266087c08SValeria Barra   // The H1 basis uses Lobatto quadrature points as nodes.
353ded9b81dSJeremy L Thompson   CeedLobattoQuadrature(p, nodes, NULL);  // nodes are in [-1,1]
354d37d859eSJeremy L Thompson   for (CeedInt i = 0; i < p; i++) nodes[i] = 0.5 + 0.5 * nodes[i];
355d1d35e2fSjeremylt   for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) {
356d1d35e2fSjeremylt     CeedInt r_nodes = gs_nodes;
357990fdeb6SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
358d1d35e2fSjeremylt       CeedInt d_1d                       = r_nodes % nd[d];
3592b730f8bSJeremy L Thompson       coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d];
360d1d35e2fSjeremylt       r_nodes /= nd[d];
36166087c08SValeria Barra     }
36266087c08SValeria Barra   }
36366087c08SValeria Barra   free(nodes);
36466087c08SValeria Barra   CeedVectorRestoreArray(mesh_coords, &coords);
36566087c08SValeria Barra   return 0;
36666087c08SValeria Barra }
36766087c08SValeria Barra 
36866087c08SValeria Barra #ifndef M_PI
36966087c08SValeria Barra #define M_PI 3.14159265358979323846
37066087c08SValeria Barra #define M_PI_2 1.57079632679489661923
37166087c08SValeria Barra #endif
37266087c08SValeria Barra 
3732b730f8bSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) {
37466087c08SValeria Barra   CeedScalar  exact_volume;
37566087c08SValeria Barra   CeedScalar *coords;
37666087c08SValeria Barra   CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
37766087c08SValeria Barra   if (dim == 1) {
37866087c08SValeria Barra     for (CeedInt i = 0; i < mesh_size; i++) {
37966087c08SValeria Barra       // map [0,1] to [0,1] varying the mesh density
38066087c08SValeria Barra       coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5));
38166087c08SValeria Barra     }
38266087c08SValeria Barra     exact_volume = 1.;
38366087c08SValeria Barra   } else {
38466087c08SValeria Barra     CeedInt num_nodes = mesh_size / dim;
38566087c08SValeria Barra     for (CeedInt i = 0; i < num_nodes; i++) {
38666087c08SValeria Barra       // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar
38766087c08SValeria Barra       // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi
38866087c08SValeria Barra       CeedScalar u = coords[i], v = coords[i + num_nodes];
38966087c08SValeria Barra       u                     = 1. + u;
39066087c08SValeria Barra       v                     = M_PI_2 * v;
39166087c08SValeria Barra       coords[i]             = u * cos(v);
39266087c08SValeria Barra       coords[i + num_nodes] = u * sin(v);
39366087c08SValeria Barra     }
39466087c08SValeria Barra     exact_volume = 3. / 4. * M_PI;
39566087c08SValeria Barra   }
39666087c08SValeria Barra   CeedVectorRestoreArray(mesh_coords, &coords);
39766087c08SValeria Barra   return exact_volume;
39866087c08SValeria Barra }
399