xref: /libCEED/examples/ceed/ex3-volume.c (revision 0b96b02dcf565e79672f1c0ce41352c1d7cf48d3)
1*0b96b02dSJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2*0b96b02dSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*0b96b02dSJeremy L Thompson //
4*0b96b02dSJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5*0b96b02dSJeremy L Thompson //
6*0b96b02dSJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7*0b96b02dSJeremy L Thompson 
8*0b96b02dSJeremy L Thompson //                             libCEED Example 1
9*0b96b02dSJeremy L Thompson //
10*0b96b02dSJeremy 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.
11*0b96b02dSJeremy L Thompson // This example also uses a diffusion operator, which provides zero contribution to the computed volume but demonstrates libCEED's ability
12*0b96b02dSJeremy L Thompson // to handle multiple basis evaluation modes for the same input and output vectors.
13*0b96b02dSJeremy L Thompson // Arbitrary mesh and solution degrees in 1D, 2D and 3D are supported from the same code.
14*0b96b02dSJeremy L Thompson //
15*0b96b02dSJeremy L Thompson // The example has no dependencies, and is designed to be self-contained.
16*0b96b02dSJeremy L Thompson // For additional examples that use external discretization libraries (MFEM, PETSc, etc.) see the subdirectories in libceed/examples.
17*0b96b02dSJeremy L Thompson //
18*0b96b02dSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command line argument (-ceed).
19*0b96b02dSJeremy L Thompson //
20*0b96b02dSJeremy L Thompson // Build with:
21*0b96b02dSJeremy L Thompson //
22*0b96b02dSJeremy L Thompson //     make ex1-volume [CEED_DIR=</path/to/libceed>]
23*0b96b02dSJeremy L Thompson //
24*0b96b02dSJeremy L Thompson // Sample runs:
25*0b96b02dSJeremy L Thompson //
26*0b96b02dSJeremy L Thompson //     ./ex3-volume
27*0b96b02dSJeremy L Thompson //     ./ex3-volume -ceed /cpu/self
28*0b96b02dSJeremy L Thompson //     ./ex3-volume -ceed /gpu/cuda
29*0b96b02dSJeremy L Thompson //
30*0b96b02dSJeremy L Thompson // Test in 1D-3D
31*0b96b02dSJeremy L Thompson //TESTARGS(name="1D User QFunction") -ceed {ceed_resource} -d 1 -t
32*0b96b02dSJeremy L Thompson //TESTARGS(name="2D User QFunction") -ceed {ceed_resource} -d 2 -t
33*0b96b02dSJeremy L Thompson //TESTARGS(name="3D User QFunction") -ceed {ceed_resource} -d 3 -t
34*0b96b02dSJeremy L Thompson 
35*0b96b02dSJeremy L Thompson /// @file
36*0b96b02dSJeremy L Thompson /// libCEED example using mass operator to compute volume
37*0b96b02dSJeremy L Thompson 
38*0b96b02dSJeremy L Thompson #include "ex3-volume.h"
39*0b96b02dSJeremy L Thompson 
40*0b96b02dSJeremy L Thompson #include <ceed.h>
41*0b96b02dSJeremy L Thompson #include <math.h>
42*0b96b02dSJeremy L Thompson #include <stdio.h>
43*0b96b02dSJeremy L Thompson #include <stdlib.h>
44*0b96b02dSJeremy L Thompson #include <string.h>
45*0b96b02dSJeremy L Thompson 
46*0b96b02dSJeremy L Thompson // Auxiliary functions
47*0b96b02dSJeremy L Thompson int        GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]);
48*0b96b02dSJeremy L Thompson int        BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
49*0b96b02dSJeremy L Thompson                                      CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction);
50*0b96b02dSJeremy L Thompson int        SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords);
51*0b96b02dSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords);
52*0b96b02dSJeremy L Thompson 
53*0b96b02dSJeremy L Thompson // Main example
54*0b96b02dSJeremy L Thompson int main(int argc, const char *argv[]) {
55*0b96b02dSJeremy L Thompson   const char *ceed_spec   = "/cpu/self";
56*0b96b02dSJeremy L Thompson   CeedInt     dim         = 3;               // dimension of the mesh
57*0b96b02dSJeremy L Thompson   CeedInt     num_comp_x  = 3;               // number of x components
58*0b96b02dSJeremy L Thompson   CeedInt     mesh_degree = 4;               // polynomial degree for the mesh
59*0b96b02dSJeremy L Thompson   CeedInt     sol_degree  = 4;               // polynomial degree for the solution
60*0b96b02dSJeremy L Thompson   CeedInt     num_qpts    = sol_degree + 2;  // number of 1D quadrature points
61*0b96b02dSJeremy L Thompson   CeedInt     prob_size   = -1;              // approximate problem size
62*0b96b02dSJeremy L Thompson   CeedInt     help = 0, test = 0, benchmark = 0;
63*0b96b02dSJeremy L Thompson 
64*0b96b02dSJeremy L Thompson   // Process command line arguments.
65*0b96b02dSJeremy L Thompson   for (int ia = 1; ia < argc; ia++) {
66*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
67*0b96b02dSJeremy L Thompson     int next_arg = ((ia + 1) < argc), parse_error = 0;
68*0b96b02dSJeremy L Thompson     if (!strcmp(argv[ia], "-h")) {
69*0b96b02dSJeremy L Thompson       help = 1;
70*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-c") || !strcmp(argv[ia], "-ceed")) {
71*0b96b02dSJeremy L Thompson       parse_error = next_arg ? ceed_spec = argv[++ia], 0 : 1;
72*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-d")) {
73*0b96b02dSJeremy L Thompson       parse_error = next_arg ? dim = atoi(argv[++ia]), 0 : 1;
74*0b96b02dSJeremy L Thompson       num_comp_x                   = dim;
75*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-m")) {
76*0b96b02dSJeremy L Thompson       parse_error = next_arg ? mesh_degree = atoi(argv[++ia]), 0 : 1;
77*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-p")) {
78*0b96b02dSJeremy L Thompson       parse_error = next_arg ? sol_degree = atoi(argv[++ia]), 0 : 1;
79*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-q")) {
80*0b96b02dSJeremy L Thompson       parse_error = next_arg ? num_qpts = atoi(argv[++ia]), 0 : 1;
81*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-s")) {
82*0b96b02dSJeremy L Thompson       parse_error = next_arg ? prob_size = atoi(argv[++ia]), 0 : 1;
83*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-b")) {
84*0b96b02dSJeremy L Thompson       parse_error = next_arg ? benchmark = atoi(argv[++ia]), 0 : 1;
85*0b96b02dSJeremy L Thompson     } else if (!strcmp(argv[ia], "-t")) {
86*0b96b02dSJeremy L Thompson       test = 1;
87*0b96b02dSJeremy L Thompson     }
88*0b96b02dSJeremy L Thompson     if (parse_error) {
89*0b96b02dSJeremy L Thompson       printf("Error parsing command line options.\n");
90*0b96b02dSJeremy L Thompson       return 1;
91*0b96b02dSJeremy L Thompson     }
92*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
93*0b96b02dSJeremy L Thompson   }
94*0b96b02dSJeremy L Thompson   if (prob_size < 0) prob_size = test ? 8 * 16 : 256 * 1024;
95*0b96b02dSJeremy L Thompson 
96*0b96b02dSJeremy L Thompson   // Print the values of all options:
97*0b96b02dSJeremy L Thompson   if (!test || help) {
98*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
99*0b96b02dSJeremy L Thompson     printf("Selected options: [command line option] : <current value>\n");
100*0b96b02dSJeremy L Thompson     printf("  Ceed specification     [-c] : %s\n", ceed_spec);
101*0b96b02dSJeremy L Thompson     printf("  Mesh dimension         [-d] : %" CeedInt_FMT "\n", dim);
102*0b96b02dSJeremy L Thompson     printf("  Mesh degree            [-m] : %" CeedInt_FMT "\n", mesh_degree);
103*0b96b02dSJeremy L Thompson     printf("  Solution degree        [-p] : %" CeedInt_FMT "\n", sol_degree);
104*0b96b02dSJeremy L Thompson     printf("  Num. 1D quadrature pts [-q] : %" CeedInt_FMT "\n", num_qpts);
105*0b96b02dSJeremy L Thompson     printf("  Approx. # unknowns     [-s] : %" CeedInt_FMT "\n", prob_size);
106*0b96b02dSJeremy L Thompson     printf("  QFunction source            : header");
107*0b96b02dSJeremy L Thompson     if (help) {
108*0b96b02dSJeremy L Thompson       printf("Test/quiet mode is %s\n", (test ? "ON" : "OFF (use -t to enable)"));
109*0b96b02dSJeremy L Thompson       return 0;
110*0b96b02dSJeremy L Thompson     }
111*0b96b02dSJeremy L Thompson     printf("\n");
112*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
113*0b96b02dSJeremy L Thompson   }
114*0b96b02dSJeremy L Thompson 
115*0b96b02dSJeremy L Thompson   // Select appropriate backend and logical device based on the (-ceed) command line argument.
116*0b96b02dSJeremy L Thompson   Ceed ceed;
117*0b96b02dSJeremy L Thompson 
118*0b96b02dSJeremy L Thompson   CeedInit(ceed_spec, &ceed);
119*0b96b02dSJeremy L Thompson 
120*0b96b02dSJeremy L Thompson   // Construct the mesh and solution bases.
121*0b96b02dSJeremy L Thompson   CeedBasis mesh_basis, sol_basis;
122*0b96b02dSJeremy L Thompson 
123*0b96b02dSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis);
124*0b96b02dSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis);
125*0b96b02dSJeremy L Thompson 
126*0b96b02dSJeremy L Thompson   // Determine the mesh size based on the given approximate problem size.
127*0b96b02dSJeremy L Thompson   CeedInt num_xyz[dim];
128*0b96b02dSJeremy L Thompson 
129*0b96b02dSJeremy L Thompson   GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz);
130*0b96b02dSJeremy L Thompson   if (!test) {
131*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
132*0b96b02dSJeremy L Thompson     printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]);
133*0b96b02dSJeremy L Thompson     if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]);
134*0b96b02dSJeremy L Thompson     if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]);
135*0b96b02dSJeremy L Thompson     printf("\n");
136*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
137*0b96b02dSJeremy L Thompson   }
138*0b96b02dSJeremy L Thompson 
139*0b96b02dSJeremy L Thompson   // Build CeedElemRestriction objects describing the mesh and solution discrete representations.
140*0b96b02dSJeremy L Thompson   CeedInt             mesh_size, sol_size;
141*0b96b02dSJeremy L Thompson   CeedElemRestriction mesh_restriction, sol_restriction, q_data_restriction;
142*0b96b02dSJeremy L Thompson 
143*0b96b02dSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restriction, NULL);
144*0b96b02dSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1 + dim * (dim + 1) / 2, &sol_size, num_qpts, NULL, &q_data_restriction);
145*0b96b02dSJeremy L Thompson   BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restriction, NULL);
146*0b96b02dSJeremy L Thompson   if (!test) {
147*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
148*0b96b02dSJeremy L Thompson     printf("Number of mesh nodes     : %" CeedInt_FMT "\n", mesh_size / dim);
149*0b96b02dSJeremy L Thompson     printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size);
150*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
151*0b96b02dSJeremy L Thompson   }
152*0b96b02dSJeremy L Thompson 
153*0b96b02dSJeremy L Thompson   // Create a CeedVector with the mesh coordinates.
154*0b96b02dSJeremy L Thompson   CeedVector mesh_coords;
155*0b96b02dSJeremy L Thompson 
156*0b96b02dSJeremy L Thompson   CeedVectorCreate(ceed, mesh_size, &mesh_coords);
157*0b96b02dSJeremy L Thompson   SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords);
158*0b96b02dSJeremy L Thompson 
159*0b96b02dSJeremy L Thompson   // Apply a transformation to the mesh.
160*0b96b02dSJeremy L Thompson   CeedScalar exact_volume = TransformMeshCoords(dim, mesh_size, mesh_coords);
161*0b96b02dSJeremy L Thompson 
162*0b96b02dSJeremy L Thompson   // Context data to be passed to the 'build_mass_diff' QFunction.
163*0b96b02dSJeremy L Thompson   CeedQFunctionContext build_ctx;
164*0b96b02dSJeremy L Thompson   struct BuildContext  build_ctx_data;
165*0b96b02dSJeremy L Thompson 
166*0b96b02dSJeremy L Thompson   build_ctx_data.dim = build_ctx_data.space_dim = dim;
167*0b96b02dSJeremy L Thompson   CeedQFunctionContextCreate(ceed, &build_ctx);
168*0b96b02dSJeremy L Thompson   CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
169*0b96b02dSJeremy L Thompson 
170*0b96b02dSJeremy L Thompson   // Create the QFunction that builds the mass + diffusion operator (i.e. computes its quadrature data) and set its context data.
171*0b96b02dSJeremy L Thompson   CeedQFunction qf_build;
172*0b96b02dSJeremy L Thompson 
173*0b96b02dSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, build_mass_diff, build_mass_diff_loc, &qf_build);
174*0b96b02dSJeremy L Thompson   CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
175*0b96b02dSJeremy L Thompson   CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
176*0b96b02dSJeremy L Thompson   CeedQFunctionAddOutput(qf_build, "qdata", 1 + dim * (dim + 1) / 2, CEED_EVAL_NONE);
177*0b96b02dSJeremy L Thompson   CeedQFunctionSetContext(qf_build, build_ctx);
178*0b96b02dSJeremy L Thompson 
179*0b96b02dSJeremy L Thompson   // Create the operator that builds the quadrature data for the mass + diffusion operator.
180*0b96b02dSJeremy L Thompson   CeedOperator op_build;
181*0b96b02dSJeremy L Thompson 
182*0b96b02dSJeremy L Thompson   CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
183*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_build, "dx", mesh_restriction, mesh_basis, CEED_VECTOR_ACTIVE);
184*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE);
185*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
186*0b96b02dSJeremy L Thompson 
187*0b96b02dSJeremy L Thompson   // Compute the quadrature data for the mass + diffusion operator.
188*0b96b02dSJeremy L Thompson   CeedVector q_data;
189*0b96b02dSJeremy L Thompson   CeedInt    elem_qpts = CeedIntPow(num_qpts, dim);
190*0b96b02dSJeremy L Thompson   CeedInt    num_elem  = 1;
191*0b96b02dSJeremy L Thompson 
192*0b96b02dSJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d];
193*0b96b02dSJeremy L Thompson   CeedVectorCreate(ceed, num_elem * elem_qpts * (1 + dim * (dim + 1) / 2), &q_data);
194*0b96b02dSJeremy L Thompson   CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE);
195*0b96b02dSJeremy L Thompson 
196*0b96b02dSJeremy L Thompson   // Create the QFunction that defines the action of the mass + diffusion operator.
197*0b96b02dSJeremy L Thompson   CeedQFunction qf_apply;
198*0b96b02dSJeremy L Thompson 
199*0b96b02dSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, apply_mass_diff, apply_mass_diff_loc, &qf_apply);
200*0b96b02dSJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP);
201*0b96b02dSJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "du", dim, CEED_EVAL_GRAD);
202*0b96b02dSJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "qdata", 1 + dim * (dim + 1) / 2, CEED_EVAL_NONE);
203*0b96b02dSJeremy L Thompson   CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP);
204*0b96b02dSJeremy L Thompson   CeedQFunctionAddOutput(qf_apply, "dv", dim, CEED_EVAL_GRAD);
205*0b96b02dSJeremy L Thompson   CeedQFunctionSetContext(qf_apply, build_ctx);
206*0b96b02dSJeremy L Thompson 
207*0b96b02dSJeremy L Thompson   // Create the mass +diffusion operator.
208*0b96b02dSJeremy L Thompson   CeedOperator op_apply;
209*0b96b02dSJeremy L Thompson 
210*0b96b02dSJeremy L Thompson   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
211*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
212*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_apply, "du", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
213*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data);
214*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
215*0b96b02dSJeremy L Thompson   CeedOperatorSetField(op_apply, "dv", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
216*0b96b02dSJeremy L Thompson 
217*0b96b02dSJeremy L Thompson   // Create auxiliary solution-size vectors.
218*0b96b02dSJeremy L Thompson   CeedVector u, v;
219*0b96b02dSJeremy L Thompson 
220*0b96b02dSJeremy L Thompson   CeedVectorCreate(ceed, sol_size, &u);
221*0b96b02dSJeremy L Thompson   CeedVectorCreate(ceed, sol_size, &v);
222*0b96b02dSJeremy L Thompson 
223*0b96b02dSJeremy L Thompson   // Initialize 'u' with ones.
224*0b96b02dSJeremy L Thompson   CeedVectorSetValue(u, 1.0);
225*0b96b02dSJeremy L Thompson 
226*0b96b02dSJeremy L Thompson   // Compute the mesh volume using the mass + diffusion operator: volume = 1^T \cdot M \cdot 1
227*0b96b02dSJeremy L Thompson   CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
228*0b96b02dSJeremy L Thompson 
229*0b96b02dSJeremy L Thompson   // Benchmark runs
230*0b96b02dSJeremy L Thompson   if (!test && benchmark) {
231*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
232*0b96b02dSJeremy L Thompson     printf(" Executing %d benchmarking runs...\n", benchmark);
233*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
234*0b96b02dSJeremy L Thompson   }
235*0b96b02dSJeremy L Thompson   for (CeedInt i = 0; i < benchmark; i++) {
236*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
237*0b96b02dSJeremy L Thompson     CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
238*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
239*0b96b02dSJeremy L Thompson   }
240*0b96b02dSJeremy L Thompson 
241*0b96b02dSJeremy L Thompson   // Compute and print the sum of the entries of 'v' giving the mesh volume.
242*0b96b02dSJeremy L Thompson   CeedScalar volume = 0.;
243*0b96b02dSJeremy L Thompson 
244*0b96b02dSJeremy L Thompson   {
245*0b96b02dSJeremy L Thompson     const CeedScalar *v_array;
246*0b96b02dSJeremy L Thompson 
247*0b96b02dSJeremy L Thompson     CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
248*0b96b02dSJeremy L Thompson     for (CeedInt i = 0; i < sol_size; i++) volume += v_array[i];
249*0b96b02dSJeremy L Thompson     CeedVectorRestoreArrayRead(v, &v_array);
250*0b96b02dSJeremy L Thompson   }
251*0b96b02dSJeremy L Thompson   if (!test) {
252*0b96b02dSJeremy L Thompson     // LCOV_EXCL_START
253*0b96b02dSJeremy L Thompson     printf(" done.\n");
254*0b96b02dSJeremy L Thompson     printf("Exact mesh volume    : % .14g\n", exact_volume);
255*0b96b02dSJeremy L Thompson     printf("Computed mesh volume : % .14g\n", volume);
256*0b96b02dSJeremy L Thompson     printf("Volume error         : % .14g\n", volume - exact_volume);
257*0b96b02dSJeremy L Thompson     // LCOV_EXCL_STOP
258*0b96b02dSJeremy L Thompson   } else {
259*0b96b02dSJeremy L Thompson     CeedScalar tol = (dim == 1 ? 200. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5);
260*0b96b02dSJeremy L Thompson 
261*0b96b02dSJeremy L Thompson     if (fabs(volume - exact_volume) > tol) printf("Volume error : % .1e\n", volume - exact_volume);
262*0b96b02dSJeremy L Thompson   }
263*0b96b02dSJeremy L Thompson 
264*0b96b02dSJeremy L Thompson   // Free dynamically allocated memory.
265*0b96b02dSJeremy L Thompson   CeedVectorDestroy(&u);
266*0b96b02dSJeremy L Thompson   CeedVectorDestroy(&v);
267*0b96b02dSJeremy L Thompson   CeedVectorDestroy(&q_data);
268*0b96b02dSJeremy L Thompson   CeedVectorDestroy(&mesh_coords);
269*0b96b02dSJeremy L Thompson   CeedOperatorDestroy(&op_apply);
270*0b96b02dSJeremy L Thompson   CeedQFunctionDestroy(&qf_apply);
271*0b96b02dSJeremy L Thompson   CeedQFunctionContextDestroy(&build_ctx);
272*0b96b02dSJeremy L Thompson   CeedOperatorDestroy(&op_build);
273*0b96b02dSJeremy L Thompson   CeedQFunctionDestroy(&qf_build);
274*0b96b02dSJeremy L Thompson   CeedElemRestrictionDestroy(&sol_restriction);
275*0b96b02dSJeremy L Thompson   CeedElemRestrictionDestroy(&mesh_restriction);
276*0b96b02dSJeremy L Thompson   CeedElemRestrictionDestroy(&q_data_restriction);
277*0b96b02dSJeremy L Thompson   CeedBasisDestroy(&sol_basis);
278*0b96b02dSJeremy L Thompson   CeedBasisDestroy(&mesh_basis);
279*0b96b02dSJeremy L Thompson   CeedDestroy(&ceed);
280*0b96b02dSJeremy L Thompson   return 0;
281*0b96b02dSJeremy L Thompson }
282*0b96b02dSJeremy L Thompson 
283*0b96b02dSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) {
284*0b96b02dSJeremy L Thompson   // Use the approximate formula:
285*0b96b02dSJeremy L Thompson   //    prob_size ~ num_elem * degree^dim
286*0b96b02dSJeremy L Thompson   CeedInt num_elem = prob_size / CeedIntPow(degree, dim);
287*0b96b02dSJeremy L Thompson   CeedInt s        = 0;  // find s: num_elem/2 < 2^s <= num_elem
288*0b96b02dSJeremy L Thompson 
289*0b96b02dSJeremy L Thompson   while (num_elem > 1) {
290*0b96b02dSJeremy L Thompson     num_elem /= 2;
291*0b96b02dSJeremy L Thompson     s++;
292*0b96b02dSJeremy L Thompson   }
293*0b96b02dSJeremy L Thompson   CeedInt r = s % dim;
294*0b96b02dSJeremy L Thompson 
295*0b96b02dSJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
296*0b96b02dSJeremy L Thompson     CeedInt sd = s / dim;
297*0b96b02dSJeremy L Thompson 
298*0b96b02dSJeremy L Thompson     if (r > 0) {
299*0b96b02dSJeremy L Thompson       sd++;
300*0b96b02dSJeremy L Thompson       r--;
301*0b96b02dSJeremy L Thompson     }
302*0b96b02dSJeremy L Thompson     num_xyz[d] = 1 << sd;
303*0b96b02dSJeremy L Thompson   }
304*0b96b02dSJeremy L Thompson   return 0;
305*0b96b02dSJeremy L Thompson }
306*0b96b02dSJeremy L Thompson 
307*0b96b02dSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
308*0b96b02dSJeremy L Thompson                               CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction) {
309*0b96b02dSJeremy L Thompson   CeedInt p         = degree + 1;
310*0b96b02dSJeremy L Thompson   CeedInt num_nodes = CeedIntPow(p, dim);         // number of scalar nodes per element
311*0b96b02dSJeremy L Thompson   CeedInt elem_qpts = CeedIntPow(num_qpts, dim);  // number of qpts per element
312*0b96b02dSJeremy L Thompson   CeedInt nd[3], num_elem = 1, scalar_size = 1;
313*0b96b02dSJeremy L Thompson 
314*0b96b02dSJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
315*0b96b02dSJeremy L Thompson     num_elem *= num_xyz[d];
316*0b96b02dSJeremy L Thompson     nd[d] = num_xyz[d] * (p - 1) + 1;
317*0b96b02dSJeremy L Thompson     scalar_size *= nd[d];
318*0b96b02dSJeremy L Thompson   }
319*0b96b02dSJeremy L Thompson   *size = scalar_size * num_comp;
320*0b96b02dSJeremy L Thompson   // elem:         0             1                 n-1
321*0b96b02dSJeremy L Thompson   //           |---*-...-*---|---*-...-*---|- ... -|--...--|
322*0b96b02dSJeremy L Thompson   // num_nodes:   0   1    p-1  p  p+1       2*p             n*p
323*0b96b02dSJeremy L Thompson   CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes);
324*0b96b02dSJeremy L Thompson 
325*0b96b02dSJeremy L Thompson   for (CeedInt e = 0; e < num_elem; e++) {
326*0b96b02dSJeremy L Thompson     CeedInt e_xyz[3] = {1, 1, 1}, re = e;
327*0b96b02dSJeremy L Thompson 
328*0b96b02dSJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
329*0b96b02dSJeremy L Thompson       e_xyz[d] = re % num_xyz[d];
330*0b96b02dSJeremy L Thompson       re /= num_xyz[d];
331*0b96b02dSJeremy L Thompson     }
332*0b96b02dSJeremy L Thompson     CeedInt *local_elem_nodes = elem_nodes + e * num_nodes;
333*0b96b02dSJeremy L Thompson 
334*0b96b02dSJeremy L Thompson     for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) {
335*0b96b02dSJeremy L Thompson       CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes;
336*0b96b02dSJeremy L Thompson 
337*0b96b02dSJeremy L Thompson       for (CeedInt d = 0; d < dim; d++) {
338*0b96b02dSJeremy L Thompson         g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride;
339*0b96b02dSJeremy L Thompson         g_nodes_stride *= nd[d];
340*0b96b02dSJeremy L Thompson         r_nodes /= p;
341*0b96b02dSJeremy L Thompson       }
342*0b96b02dSJeremy L Thompson       local_elem_nodes[l_nodes] = g_nodes;
343*0b96b02dSJeremy L Thompson     }
344*0b96b02dSJeremy L Thompson   }
345*0b96b02dSJeremy L Thompson   if (restriction) {
346*0b96b02dSJeremy L Thompson     CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes,
347*0b96b02dSJeremy L Thompson                               restriction);
348*0b96b02dSJeremy L Thompson   }
349*0b96b02dSJeremy L Thompson   if (q_data_restriction) {
350*0b96b02dSJeremy L Thompson     CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, q_data_restriction);
351*0b96b02dSJeremy L Thompson   }
352*0b96b02dSJeremy L Thompson   free(elem_nodes);
353*0b96b02dSJeremy L Thompson   return 0;
354*0b96b02dSJeremy L Thompson }
355*0b96b02dSJeremy L Thompson 
356*0b96b02dSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) {
357*0b96b02dSJeremy L Thompson   CeedInt p = mesh_degree + 1;
358*0b96b02dSJeremy L Thompson   CeedInt nd[3], scalar_size = 1;
359*0b96b02dSJeremy L Thompson 
360*0b96b02dSJeremy L Thompson   for (CeedInt d = 0; d < dim; d++) {
361*0b96b02dSJeremy L Thompson     nd[d] = num_xyz[d] * (p - 1) + 1;
362*0b96b02dSJeremy L Thompson     scalar_size *= nd[d];
363*0b96b02dSJeremy L Thompson   }
364*0b96b02dSJeremy L Thompson   CeedScalar *coords;
365*0b96b02dSJeremy L Thompson 
366*0b96b02dSJeremy L Thompson   CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords);
367*0b96b02dSJeremy L Thompson   CeedScalar *nodes = malloc(sizeof(CeedScalar) * p);
368*0b96b02dSJeremy L Thompson 
369*0b96b02dSJeremy L Thompson   // The H1 basis uses Lobatto quadrature points as nodes.
370*0b96b02dSJeremy L Thompson   CeedLobattoQuadrature(p, nodes, NULL);  // nodes are in [-1,1]
371*0b96b02dSJeremy L Thompson   for (CeedInt i = 0; i < p; i++) nodes[i] = 0.5 + 0.5 * nodes[i];
372*0b96b02dSJeremy L Thompson   for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) {
373*0b96b02dSJeremy L Thompson     CeedInt r_nodes = gs_nodes;
374*0b96b02dSJeremy L Thompson 
375*0b96b02dSJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
376*0b96b02dSJeremy L Thompson       CeedInt d_1d                       = r_nodes % nd[d];
377*0b96b02dSJeremy L Thompson       coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d];
378*0b96b02dSJeremy L Thompson       r_nodes /= nd[d];
379*0b96b02dSJeremy L Thompson     }
380*0b96b02dSJeremy L Thompson   }
381*0b96b02dSJeremy L Thompson   free(nodes);
382*0b96b02dSJeremy L Thompson   CeedVectorRestoreArray(mesh_coords, &coords);
383*0b96b02dSJeremy L Thompson   return 0;
384*0b96b02dSJeremy L Thompson }
385*0b96b02dSJeremy L Thompson 
386*0b96b02dSJeremy L Thompson #ifndef M_PI
387*0b96b02dSJeremy L Thompson #define M_PI 3.14159265358979323846
388*0b96b02dSJeremy L Thompson #define M_PI_2 1.57079632679489661923
389*0b96b02dSJeremy L Thompson #endif
390*0b96b02dSJeremy L Thompson 
391*0b96b02dSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) {
392*0b96b02dSJeremy L Thompson   CeedScalar  exact_volume;
393*0b96b02dSJeremy L Thompson   CeedScalar *coords;
394*0b96b02dSJeremy L Thompson 
395*0b96b02dSJeremy L Thompson   CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
396*0b96b02dSJeremy L Thompson   if (dim == 1) {
397*0b96b02dSJeremy L Thompson     for (CeedInt i = 0; i < mesh_size; i++) {
398*0b96b02dSJeremy L Thompson       // map [0,1] to [0,1] varying the mesh density
399*0b96b02dSJeremy L Thompson       coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5));
400*0b96b02dSJeremy L Thompson     }
401*0b96b02dSJeremy L Thompson     exact_volume = 1.;
402*0b96b02dSJeremy L Thompson   } else {
403*0b96b02dSJeremy L Thompson     CeedInt num_nodes = mesh_size / dim;
404*0b96b02dSJeremy L Thompson     for (CeedInt i = 0; i < num_nodes; i++) {
405*0b96b02dSJeremy L Thompson       // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar
406*0b96b02dSJeremy L Thompson       // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi
407*0b96b02dSJeremy L Thompson       CeedScalar u = coords[i], v = coords[i + num_nodes];
408*0b96b02dSJeremy L Thompson 
409*0b96b02dSJeremy L Thompson       u                     = 1. + u;
410*0b96b02dSJeremy L Thompson       v                     = M_PI_2 * v;
411*0b96b02dSJeremy L Thompson       coords[i]             = u * cos(v);
412*0b96b02dSJeremy L Thompson       coords[i + num_nodes] = u * sin(v);
413*0b96b02dSJeremy L Thompson     }
414*0b96b02dSJeremy L Thompson     exact_volume = 3. / 4. * M_PI;
415*0b96b02dSJeremy L Thompson   }
416*0b96b02dSJeremy L Thompson   CeedVectorRestoreArray(mesh_coords, &coords);
417*0b96b02dSJeremy L Thompson   return exact_volume;
418*0b96b02dSJeremy L Thompson }
419