1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
main(int argc,const char * argv[])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
6382138112SJeremy 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;
8482138112SJeremy L Thompson } else if (!strcmp(argv[ia], "-b")) {
8582138112SJeremy 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;
1200b96b02dSJeremy L Thompson
12166087c08SValeria Barra CeedInit(ceed_spec, &ceed);
12266087c08SValeria Barra
12366087c08SValeria Barra // Construct the mesh and solution bases.
12466087c08SValeria Barra CeedBasis mesh_basis, sol_basis;
1250b96b02dSJeremy L Thompson
1262b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis);
1272b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis);
12866087c08SValeria Barra
12966087c08SValeria Barra // Determine the mesh size based on the given approximate problem size.
130990fdeb6SJeremy L Thompson CeedInt num_xyz[dim];
1310b96b02dSJeremy L Thompson
132d1d35e2fSjeremylt GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz);
13366087c08SValeria Barra if (!test) {
134ded9b81dSJeremy L Thompson // LCOV_EXCL_START
135990fdeb6SJeremy L Thompson printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]);
1362b730f8bSJeremy L Thompson if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]);
1372b730f8bSJeremy L Thompson if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]);
13866087c08SValeria Barra printf("\n");
139ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP
14066087c08SValeria Barra }
14166087c08SValeria Barra
142ea61e9acSJeremy L Thompson // Build CeedElemRestriction objects describing the mesh and solution discrete representations.
14366087c08SValeria Barra CeedInt mesh_size, sol_size;
144d37d859eSJeremy L Thompson CeedElemRestriction mesh_restriction, sol_restriction, q_data_restriction;
1450b96b02dSJeremy L Thompson
146d37d859eSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restriction, NULL);
147d37d859eSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restriction, &q_data_restriction);
14866087c08SValeria Barra if (!test) {
149e15f9bd0SJeremy L Thompson // LCOV_EXCL_START
150990fdeb6SJeremy L Thompson printf("Number of mesh nodes : %" CeedInt_FMT "\n", mesh_size / dim);
151990fdeb6SJeremy L Thompson printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size);
152e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP
15366087c08SValeria Barra }
15466087c08SValeria Barra
15566087c08SValeria Barra // Create a CeedVector with the mesh coordinates.
15666087c08SValeria Barra CeedVector mesh_coords;
1570b96b02dSJeremy L Thompson
15866087c08SValeria Barra CeedVectorCreate(ceed, mesh_size, &mesh_coords);
159d1d35e2fSjeremylt SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords);
16066087c08SValeria Barra
16166087c08SValeria Barra // Apply a transformation to the mesh.
162d37d859eSJeremy L Thompson CeedScalar exact_volume = TransformMeshCoords(dim, mesh_size, mesh_coords);
16366087c08SValeria Barra
164d37d859eSJeremy L Thompson // Context data to be passed to the 'build_mass' QFunction.
165777ff853SJeremy L Thompson CeedQFunctionContext build_ctx;
166777ff853SJeremy L Thompson struct BuildContext build_ctx_data;
1670b96b02dSJeremy L Thompson
168777ff853SJeremy L Thompson build_ctx_data.dim = build_ctx_data.space_dim = dim;
169777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &build_ctx);
1702b730f8bSJeremy L Thompson CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
17166087c08SValeria Barra
172ea61e9acSJeremy L Thompson // Create the QFunction that builds the mass operator (i.e. computes its quadrature data) and set its context data.
173d1d35e2fSjeremylt CeedQFunction qf_build;
1740b96b02dSJeremy L Thompson
175d37d859eSJeremy L Thompson if (gallery) {
17666087c08SValeria Barra // This creates the QFunction via the gallery.
17766087c08SValeria Barra char name[13] = "";
178990fdeb6SJeremy L Thompson snprintf(name, sizeof name, "Mass%" CeedInt_FMT "DBuild", dim);
179d1d35e2fSjeremylt CeedQFunctionCreateInteriorByName(ceed, name, &qf_build);
180d37d859eSJeremy L Thompson } else {
181d37d859eSJeremy L Thompson // This creates the QFunction directly.
182d37d859eSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, build_mass, build_mass_loc, &qf_build);
183d37d859eSJeremy L Thompson CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
184d37d859eSJeremy L Thompson CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
185d37d859eSJeremy L Thompson CeedQFunctionAddOutput(qf_build, "qdata", 1, CEED_EVAL_NONE);
186d37d859eSJeremy L Thompson CeedQFunctionSetContext(qf_build, build_ctx);
18766087c08SValeria Barra }
18866087c08SValeria Barra
18966087c08SValeria Barra // Create the operator that builds the quadrature data for the mass operator.
190d1d35e2fSjeremylt CeedOperator op_build;
1910b96b02dSJeremy L Thompson
1922b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
193d37d859eSJeremy L Thompson CeedOperatorSetField(op_build, "dx", mesh_restriction, mesh_basis, CEED_VECTOR_ACTIVE);
1942b730f8bSJeremy L Thompson CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE);
195356036faSJeremy L Thompson CeedOperatorSetField(op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
19666087c08SValeria Barra
19766087c08SValeria Barra // Compute the quadrature data for the mass operator.
198d1d35e2fSjeremylt CeedVector q_data;
19966087c08SValeria Barra CeedInt elem_qpts = CeedIntPow(num_qpts, dim);
20066087c08SValeria Barra CeedInt num_elem = 1;
2010b96b02dSJeremy L Thompson
2022b730f8bSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d];
203d1d35e2fSjeremylt CeedVectorCreate(ceed, num_elem * elem_qpts, &q_data);
2042b730f8bSJeremy L Thompson CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE);
20566087c08SValeria Barra
206ded9b81dSJeremy L Thompson // Create the QFunction that defines the action of the mass operator.
207d1d35e2fSjeremylt CeedQFunction qf_apply;
2080b96b02dSJeremy L Thompson
209d37d859eSJeremy L Thompson if (gallery) {
210d37d859eSJeremy L Thompson // This creates the QFunction via the gallery.
211d37d859eSJeremy L Thompson CeedQFunctionCreateInteriorByName(ceed, "MassApply", &qf_apply);
212d37d859eSJeremy L Thompson } else {
21366087c08SValeria Barra // This creates the QFunction directly.
214d37d859eSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, apply_mass, apply_mass_loc, &qf_apply);
215d1d35e2fSjeremylt CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP);
216d1d35e2fSjeremylt CeedQFunctionAddInput(qf_apply, "qdata", 1, CEED_EVAL_NONE);
217d1d35e2fSjeremylt CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP);
21866087c08SValeria Barra }
21966087c08SValeria Barra
22066087c08SValeria Barra // Create the mass operator.
221d1d35e2fSjeremylt CeedOperator op_apply;
2220b96b02dSJeremy L Thompson
2232b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
224d37d859eSJeremy L Thompson CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
225356036faSJeremy L Thompson CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data);
226d37d859eSJeremy L Thompson CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
22766087c08SValeria Barra
22866087c08SValeria Barra // Create auxiliary solution-size vectors.
22966087c08SValeria Barra CeedVector u, v;
2300b96b02dSJeremy L Thompson
23166087c08SValeria Barra CeedVectorCreate(ceed, sol_size, &u);
23266087c08SValeria Barra CeedVectorCreate(ceed, sol_size, &v);
23366087c08SValeria Barra
234b8bc0c60SPeter Munch // Initialize 'u' with ones.
23566087c08SValeria Barra CeedVectorSetValue(u, 1.0);
23666087c08SValeria Barra
237d37d859eSJeremy L Thompson // Compute the mesh volume using the mass operator: volume = 1^T \cdot M \cdot 1
238d1d35e2fSjeremylt CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
23966087c08SValeria Barra
24082138112SJeremy L Thompson // Benchmark runs
24182138112SJeremy L Thompson if (!test && benchmark) {
24282138112SJeremy L Thompson // LCOV_EXCL_START
24382138112SJeremy L Thompson printf(" Executing %d benchmarking runs...\n", benchmark);
24482138112SJeremy L Thompson // LCOV_EXCL_STOP
24582138112SJeremy L Thompson }
24682138112SJeremy L Thompson for (CeedInt i = 0; i < benchmark; i++) {
24782138112SJeremy L Thompson // LCOV_EXCL_START
24882138112SJeremy L Thompson CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
24982138112SJeremy L Thompson // LCOV_EXCL_STOP
25082138112SJeremy L Thompson }
25182138112SJeremy L Thompson
25266087c08SValeria Barra // Compute and print the sum of the entries of 'v' giving the mesh volume.
253d37d859eSJeremy L Thompson CeedScalar volume = 0.;
2540b96b02dSJeremy L Thompson
255d37d859eSJeremy L Thompson {
256d1d35e2fSjeremylt const CeedScalar *v_array;
2570b96b02dSJeremy L Thompson
258d1d35e2fSjeremylt CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
259d37d859eSJeremy L Thompson for (CeedInt i = 0; i < sol_size; i++) volume += v_array[i];
260d1d35e2fSjeremylt CeedVectorRestoreArrayRead(v, &v_array);
261d37d859eSJeremy L Thompson }
26266087c08SValeria Barra if (!test) {
263ded9b81dSJeremy L Thompson // LCOV_EXCL_START
26466087c08SValeria Barra printf(" done.\n");
265d37d859eSJeremy L Thompson printf("Exact mesh volume : % .14g\n", exact_volume);
266d37d859eSJeremy L Thompson printf("Computed mesh volume : % .14g\n", volume);
267d37d859eSJeremy L Thompson printf("Volume error : % .14g\n", volume - exact_volume);
268ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP
26966087c08SValeria Barra } else {
270bd882c8aSJames Wright CeedScalar tol = (dim == 1 ? 200. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5);
2710b96b02dSJeremy L Thompson
272d37d859eSJeremy L Thompson if (fabs(volume - exact_volume) > tol) printf("Volume error : % .1e\n", volume - exact_volume);
27366087c08SValeria Barra }
27466087c08SValeria Barra
27566087c08SValeria Barra // Free dynamically allocated memory.
27666087c08SValeria Barra CeedVectorDestroy(&u);
27766087c08SValeria Barra CeedVectorDestroy(&v);
278d1d35e2fSjeremylt CeedVectorDestroy(&q_data);
27966087c08SValeria Barra CeedVectorDestroy(&mesh_coords);
280d1d35e2fSjeremylt CeedOperatorDestroy(&op_apply);
281d1d35e2fSjeremylt CeedQFunctionDestroy(&qf_apply);
282777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&build_ctx);
283d1d35e2fSjeremylt CeedOperatorDestroy(&op_build);
284d1d35e2fSjeremylt CeedQFunctionDestroy(&qf_build);
285d37d859eSJeremy L Thompson CeedElemRestrictionDestroy(&sol_restriction);
286d37d859eSJeremy L Thompson CeedElemRestrictionDestroy(&mesh_restriction);
287d37d859eSJeremy L Thompson CeedElemRestrictionDestroy(&q_data_restriction);
28866087c08SValeria Barra CeedBasisDestroy(&sol_basis);
28966087c08SValeria Barra CeedBasisDestroy(&mesh_basis);
29066087c08SValeria Barra CeedDestroy(&ceed);
29166087c08SValeria Barra return 0;
29266087c08SValeria Barra }
29366087c08SValeria Barra
GetCartesianMeshSize(CeedInt dim,CeedInt degree,CeedInt prob_size,CeedInt num_xyz[dim])2942b730f8bSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) {
29566087c08SValeria Barra // Use the approximate formula:
296ded9b81dSJeremy L Thompson // prob_size ~ num_elem * degree^dim
297ded9b81dSJeremy L Thompson CeedInt num_elem = prob_size / CeedIntPow(degree, dim);
29866087c08SValeria Barra CeedInt s = 0; // find s: num_elem/2 < 2^s <= num_elem
2990b96b02dSJeremy L Thompson
30066087c08SValeria Barra while (num_elem > 1) {
30166087c08SValeria Barra num_elem /= 2;
30266087c08SValeria Barra s++;
30366087c08SValeria Barra }
30466087c08SValeria Barra CeedInt r = s % dim;
3050b96b02dSJeremy L Thompson
306990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) {
307990fdeb6SJeremy L Thompson CeedInt sd = s / dim;
3080b96b02dSJeremy L Thompson
3092b730f8bSJeremy L Thompson if (r > 0) {
3102b730f8bSJeremy L Thompson sd++;
3112b730f8bSJeremy L Thompson r--;
3122b730f8bSJeremy L Thompson }
313d1d35e2fSjeremylt num_xyz[d] = 1 << sd;
31466087c08SValeria Barra }
31566087c08SValeria Barra return 0;
31666087c08SValeria Barra }
31766087c08SValeria Barra
BuildCartesianRestriction(Ceed ceed,CeedInt dim,CeedInt num_xyz[dim],CeedInt degree,CeedInt num_comp,CeedInt * size,CeedInt num_qpts,CeedElemRestriction * restriction,CeedElemRestriction * q_data_restriction)3182b730f8bSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
319d37d859eSJeremy L Thompson CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction) {
320ded9b81dSJeremy L Thompson CeedInt p = degree + 1;
321d1d35e2fSjeremylt CeedInt num_nodes = CeedIntPow(p, dim); // number of scalar nodes per element
32266087c08SValeria Barra CeedInt elem_qpts = CeedIntPow(num_qpts, dim); // number of qpts per element
32366087c08SValeria Barra CeedInt nd[3], num_elem = 1, scalar_size = 1;
3240b96b02dSJeremy L Thompson
325990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) {
326d1d35e2fSjeremylt num_elem *= num_xyz[d];
327d1d35e2fSjeremylt nd[d] = num_xyz[d] * (p - 1) + 1;
32866087c08SValeria Barra scalar_size *= nd[d];
32966087c08SValeria Barra }
330d1d35e2fSjeremylt *size = scalar_size * num_comp;
33166087c08SValeria Barra // elem: 0 1 n-1
33266087c08SValeria Barra // |---*-...-*---|---*-...-*---|- ... -|--...--|
333d1d35e2fSjeremylt // num_nodes: 0 1 p-1 p p+1 2*p n*p
334d1d35e2fSjeremylt CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes);
3350b96b02dSJeremy L Thompson
33666087c08SValeria Barra for (CeedInt e = 0; e < num_elem; e++) {
337d1d35e2fSjeremylt CeedInt e_xyz[3] = {1, 1, 1}, re = e;
3380b96b02dSJeremy L Thompson
3392b730f8bSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) {
3402b730f8bSJeremy L Thompson e_xyz[d] = re % num_xyz[d];
3412b730f8bSJeremy L Thompson re /= num_xyz[d];
3422b730f8bSJeremy L Thompson }
343d37d859eSJeremy L Thompson CeedInt *local_elem_nodes = elem_nodes + e * num_nodes;
3440b96b02dSJeremy L Thompson
345990fdeb6SJeremy L Thompson for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) {
346d1d35e2fSjeremylt CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes;
3470b96b02dSJeremy L Thompson
348990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) {
349d1d35e2fSjeremylt g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride;
350d1d35e2fSjeremylt g_nodes_stride *= nd[d];
351d1d35e2fSjeremylt r_nodes /= p;
35266087c08SValeria Barra }
353d37d859eSJeremy L Thompson local_elem_nodes[l_nodes] = g_nodes;
35466087c08SValeria Barra }
35566087c08SValeria Barra }
3562b730f8bSJeremy L Thompson CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes,
357d37d859eSJeremy L Thompson restriction);
3586c10af5dSJeremy L Thompson if (q_data_restriction) {
359d37d859eSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, q_data_restriction);
3606c10af5dSJeremy L Thompson }
361d1d35e2fSjeremylt free(elem_nodes);
36266087c08SValeria Barra return 0;
36366087c08SValeria Barra }
36466087c08SValeria Barra
SetCartesianMeshCoords(CeedInt dim,CeedInt num_xyz[dim],CeedInt mesh_degree,CeedVector mesh_coords)3652b730f8bSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) {
366ded9b81dSJeremy L Thompson CeedInt p = mesh_degree + 1;
3671d2f0dcbSJeremy L Thompson CeedInt nd[3], scalar_size = 1;
3680b96b02dSJeremy L Thompson
369990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) {
370d1d35e2fSjeremylt nd[d] = num_xyz[d] * (p - 1) + 1;
37166087c08SValeria Barra scalar_size *= nd[d];
37266087c08SValeria Barra }
37366087c08SValeria Barra CeedScalar *coords;
3740b96b02dSJeremy L Thompson
3759c774eddSJeremy L Thompson CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords);
376ded9b81dSJeremy L Thompson CeedScalar *nodes = malloc(sizeof(CeedScalar) * p);
3770b96b02dSJeremy L Thompson
37866087c08SValeria Barra // The H1 basis uses Lobatto quadrature points as nodes.
379ded9b81dSJeremy L Thompson CeedLobattoQuadrature(p, nodes, NULL); // nodes are in [-1,1]
380d37d859eSJeremy L Thompson for (CeedInt i = 0; i < p; i++) nodes[i] = 0.5 + 0.5 * nodes[i];
381d1d35e2fSjeremylt for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) {
382d1d35e2fSjeremylt CeedInt r_nodes = gs_nodes;
3830b96b02dSJeremy L Thompson
384990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) {
385d1d35e2fSjeremylt CeedInt d_1d = r_nodes % nd[d];
3860b96b02dSJeremy L Thompson
3872b730f8bSJeremy L Thompson coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d];
388d1d35e2fSjeremylt 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
TransformMeshCoords(CeedInt dim,CeedInt mesh_size,CeedVector mesh_coords)4012b730f8bSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) {
40266087c08SValeria Barra CeedScalar exact_volume;
40366087c08SValeria Barra CeedScalar *coords;
4040b96b02dSJeremy L Thompson
40566087c08SValeria Barra CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
40666087c08SValeria Barra if (dim == 1) {
40766087c08SValeria Barra for (CeedInt i = 0; i < mesh_size; i++) {
40866087c08SValeria Barra // map [0,1] to [0,1] varying the mesh density
40966087c08SValeria Barra coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5));
41066087c08SValeria Barra }
41166087c08SValeria Barra exact_volume = 1.;
41266087c08SValeria Barra } else {
41366087c08SValeria Barra CeedInt num_nodes = mesh_size / dim;
4140b96b02dSJeremy L Thompson
41566087c08SValeria Barra for (CeedInt i = 0; i < num_nodes; i++) {
41666087c08SValeria Barra // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar
41766087c08SValeria Barra // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi
41866087c08SValeria Barra CeedScalar u = coords[i], v = coords[i + num_nodes];
4190b96b02dSJeremy L Thompson
42066087c08SValeria Barra u = 1. + u;
42166087c08SValeria Barra v = M_PI_2 * v;
42266087c08SValeria Barra coords[i] = u * cos(v);
42366087c08SValeria Barra coords[i + num_nodes] = u * sin(v);
42466087c08SValeria Barra }
42566087c08SValeria Barra exact_volume = 3. / 4. * M_PI;
42666087c08SValeria Barra }
42766087c08SValeria Barra CeedVectorRestoreArray(mesh_coords, &coords);
42866087c08SValeria Barra return exact_volume;
42966087c08SValeria Barra }
430