1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors. 20b96b02dSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 30b96b02dSJeremy L Thompson // 40b96b02dSJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 50b96b02dSJeremy L Thompson // 60b96b02dSJeremy L Thompson // This file is part of CEED: http://github.com/ceed 70b96b02dSJeremy L Thompson 80b96b02dSJeremy L Thompson // libCEED Example 1 90b96b02dSJeremy L Thompson // 100b96b02dSJeremy 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. 110b96b02dSJeremy L Thompson // This example also uses a diffusion operator, which provides zero contribution to the computed volume but demonstrates libCEED's ability 120b96b02dSJeremy L Thompson // to handle multiple basis evaluation modes for the same input and output vectors. 130b96b02dSJeremy L Thompson // Arbitrary mesh and solution degrees in 1D, 2D and 3D are supported from the same code. 140b96b02dSJeremy L Thompson // 150b96b02dSJeremy L Thompson // The example has no dependencies, and is designed to be self-contained. 160b96b02dSJeremy L Thompson // For additional examples that use external discretization libraries (MFEM, PETSc, etc.) see the subdirectories in libceed/examples. 170b96b02dSJeremy L Thompson // 180b96b02dSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command line argument (-ceed). 190b96b02dSJeremy L Thompson // 200b96b02dSJeremy L Thompson // Build with: 210b96b02dSJeremy L Thompson // 220806d17aSJeremy L Thompson // make ex3-volume [CEED_DIR=</path/to/libceed>] 230b96b02dSJeremy L Thompson // 240b96b02dSJeremy L Thompson // Sample runs: 250b96b02dSJeremy L Thompson // 260b96b02dSJeremy L Thompson // ./ex3-volume 270b96b02dSJeremy L Thompson // ./ex3-volume -ceed /cpu/self 280b96b02dSJeremy L Thompson // ./ex3-volume -ceed /gpu/cuda 290b96b02dSJeremy L Thompson // 300b96b02dSJeremy L Thompson // Test in 1D-3D 310b96b02dSJeremy L Thompson //TESTARGS(name="1D User QFunction") -ceed {ceed_resource} -d 1 -t 320b96b02dSJeremy L Thompson //TESTARGS(name="2D User QFunction") -ceed {ceed_resource} -d 2 -t 330b96b02dSJeremy L Thompson //TESTARGS(name="3D User QFunction") -ceed {ceed_resource} -d 3 -t 340b96b02dSJeremy L Thompson 350b96b02dSJeremy L Thompson /// @file 360b96b02dSJeremy L Thompson /// libCEED example using mass operator to compute volume 370b96b02dSJeremy L Thompson 380b96b02dSJeremy L Thompson #include "ex3-volume.h" 390b96b02dSJeremy L Thompson 400b96b02dSJeremy L Thompson #include <ceed.h> 410b96b02dSJeremy L Thompson #include <math.h> 420b96b02dSJeremy L Thompson #include <stdio.h> 430b96b02dSJeremy L Thompson #include <stdlib.h> 440b96b02dSJeremy L Thompson #include <string.h> 450b96b02dSJeremy L Thompson 460b96b02dSJeremy L Thompson // Auxiliary functions 470b96b02dSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]); 480b96b02dSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts, 490b96b02dSJeremy L Thompson CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction); 500b96b02dSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords); 510b96b02dSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords); 520b96b02dSJeremy L Thompson 530b96b02dSJeremy L Thompson // Main example 540b96b02dSJeremy L Thompson int main(int argc, const char *argv[]) { 550b96b02dSJeremy L Thompson const char *ceed_spec = "/cpu/self"; 560b96b02dSJeremy L Thompson CeedInt dim = 3; // dimension of the mesh 570b96b02dSJeremy L Thompson CeedInt num_comp_x = 3; // number of x components 580b96b02dSJeremy L Thompson CeedInt mesh_degree = 4; // polynomial degree for the mesh 590b96b02dSJeremy L Thompson CeedInt sol_degree = 4; // polynomial degree for the solution 600b96b02dSJeremy L Thompson CeedInt num_qpts = sol_degree + 2; // number of 1D quadrature points 610b96b02dSJeremy L Thompson CeedInt prob_size = -1; // approximate problem size 620b96b02dSJeremy L Thompson CeedInt help = 0, test = 0, benchmark = 0; 630b96b02dSJeremy L Thompson 640b96b02dSJeremy L Thompson // Process command line arguments. 650b96b02dSJeremy L Thompson for (int ia = 1; ia < argc; ia++) { 660b96b02dSJeremy L Thompson // LCOV_EXCL_START 670b96b02dSJeremy L Thompson int next_arg = ((ia + 1) < argc), parse_error = 0; 680b96b02dSJeremy L Thompson if (!strcmp(argv[ia], "-h")) { 690b96b02dSJeremy L Thompson help = 1; 700b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-c") || !strcmp(argv[ia], "-ceed")) { 710b96b02dSJeremy L Thompson parse_error = next_arg ? ceed_spec = argv[++ia], 0 : 1; 720b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-d")) { 730b96b02dSJeremy L Thompson parse_error = next_arg ? dim = atoi(argv[++ia]), 0 : 1; 740b96b02dSJeremy L Thompson num_comp_x = dim; 750b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-m")) { 760b96b02dSJeremy L Thompson parse_error = next_arg ? mesh_degree = atoi(argv[++ia]), 0 : 1; 770b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-p")) { 780b96b02dSJeremy L Thompson parse_error = next_arg ? sol_degree = atoi(argv[++ia]), 0 : 1; 790b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-q")) { 800b96b02dSJeremy L Thompson parse_error = next_arg ? num_qpts = atoi(argv[++ia]), 0 : 1; 810b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-s")) { 820b96b02dSJeremy L Thompson parse_error = next_arg ? prob_size = atoi(argv[++ia]), 0 : 1; 830b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-b")) { 840b96b02dSJeremy L Thompson parse_error = next_arg ? benchmark = atoi(argv[++ia]), 0 : 1; 850b96b02dSJeremy L Thompson } else if (!strcmp(argv[ia], "-t")) { 860b96b02dSJeremy L Thompson test = 1; 870b96b02dSJeremy L Thompson } 880b96b02dSJeremy L Thompson if (parse_error) { 890b96b02dSJeremy L Thompson printf("Error parsing command line options.\n"); 900b96b02dSJeremy L Thompson return 1; 910b96b02dSJeremy L Thompson } 920b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 930b96b02dSJeremy L Thompson } 940b96b02dSJeremy L Thompson if (prob_size < 0) prob_size = test ? 8 * 16 : 256 * 1024; 950b96b02dSJeremy L Thompson 960b96b02dSJeremy L Thompson // Print the values of all options: 970b96b02dSJeremy L Thompson if (!test || help) { 980b96b02dSJeremy L Thompson // LCOV_EXCL_START 990b96b02dSJeremy L Thompson printf("Selected options: [command line option] : <current value>\n"); 1000b96b02dSJeremy L Thompson printf(" Ceed specification [-c] : %s\n", ceed_spec); 1010b96b02dSJeremy L Thompson printf(" Mesh dimension [-d] : %" CeedInt_FMT "\n", dim); 1020b96b02dSJeremy L Thompson printf(" Mesh degree [-m] : %" CeedInt_FMT "\n", mesh_degree); 1030b96b02dSJeremy L Thompson printf(" Solution degree [-p] : %" CeedInt_FMT "\n", sol_degree); 1040b96b02dSJeremy L Thompson printf(" Num. 1D quadrature pts [-q] : %" CeedInt_FMT "\n", num_qpts); 1050b96b02dSJeremy L Thompson printf(" Approx. # unknowns [-s] : %" CeedInt_FMT "\n", prob_size); 1060b96b02dSJeremy L Thompson printf(" QFunction source : header"); 1070b96b02dSJeremy L Thompson if (help) { 1080b96b02dSJeremy L Thompson printf("Test/quiet mode is %s\n", (test ? "ON" : "OFF (use -t to enable)")); 1090b96b02dSJeremy L Thompson return 0; 1100b96b02dSJeremy L Thompson } 1110b96b02dSJeremy L Thompson printf("\n"); 1120b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 1130b96b02dSJeremy L Thompson } 1140b96b02dSJeremy L Thompson 1150b96b02dSJeremy L Thompson // Select appropriate backend and logical device based on the (-ceed) command line argument. 1160b96b02dSJeremy L Thompson Ceed ceed; 1170b96b02dSJeremy L Thompson 1180b96b02dSJeremy L Thompson CeedInit(ceed_spec, &ceed); 1190b96b02dSJeremy L Thompson 1200b96b02dSJeremy L Thompson // Construct the mesh and solution bases. 1210b96b02dSJeremy L Thompson CeedBasis mesh_basis, sol_basis; 1220b96b02dSJeremy L Thompson 1230b96b02dSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis); 1240b96b02dSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis); 1250b96b02dSJeremy L Thompson 1260b96b02dSJeremy L Thompson // Determine the mesh size based on the given approximate problem size. 1270b96b02dSJeremy L Thompson CeedInt num_xyz[dim]; 1280b96b02dSJeremy L Thompson 1290b96b02dSJeremy L Thompson GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz); 1300b96b02dSJeremy L Thompson if (!test) { 1310b96b02dSJeremy L Thompson // LCOV_EXCL_START 1320b96b02dSJeremy L Thompson printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]); 1330b96b02dSJeremy L Thompson if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]); 1340b96b02dSJeremy L Thompson if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]); 1350b96b02dSJeremy L Thompson printf("\n"); 1360b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 1370b96b02dSJeremy L Thompson } 1380b96b02dSJeremy L Thompson 1390b96b02dSJeremy L Thompson // Build CeedElemRestriction objects describing the mesh and solution discrete representations. 1400b96b02dSJeremy L Thompson CeedInt mesh_size, sol_size; 1410b96b02dSJeremy L Thompson CeedElemRestriction mesh_restriction, sol_restriction, q_data_restriction; 1420b96b02dSJeremy L Thompson 1430b96b02dSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restriction, NULL); 1440b96b02dSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1 + dim * (dim + 1) / 2, &sol_size, num_qpts, NULL, &q_data_restriction); 1450b96b02dSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restriction, NULL); 1460b96b02dSJeremy L Thompson if (!test) { 1470b96b02dSJeremy L Thompson // LCOV_EXCL_START 1480b96b02dSJeremy L Thompson printf("Number of mesh nodes : %" CeedInt_FMT "\n", mesh_size / dim); 1490b96b02dSJeremy L Thompson printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size); 1500b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 1510b96b02dSJeremy L Thompson } 1520b96b02dSJeremy L Thompson 1530b96b02dSJeremy L Thompson // Create a CeedVector with the mesh coordinates. 1540b96b02dSJeremy L Thompson CeedVector mesh_coords; 1550b96b02dSJeremy L Thompson 1560b96b02dSJeremy L Thompson CeedVectorCreate(ceed, mesh_size, &mesh_coords); 1570b96b02dSJeremy L Thompson SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords); 1580b96b02dSJeremy L Thompson 1590b96b02dSJeremy L Thompson // Apply a transformation to the mesh. 1600b96b02dSJeremy L Thompson CeedScalar exact_volume = TransformMeshCoords(dim, mesh_size, mesh_coords); 1610b96b02dSJeremy L Thompson 1620b96b02dSJeremy L Thompson // Context data to be passed to the 'build_mass_diff' QFunction. 1630b96b02dSJeremy L Thompson CeedQFunctionContext build_ctx; 1640b96b02dSJeremy L Thompson struct BuildContext build_ctx_data; 1650b96b02dSJeremy L Thompson 1660b96b02dSJeremy L Thompson build_ctx_data.dim = build_ctx_data.space_dim = dim; 1670b96b02dSJeremy L Thompson CeedQFunctionContextCreate(ceed, &build_ctx); 1680b96b02dSJeremy L Thompson CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data); 1690b96b02dSJeremy L Thompson 1700b96b02dSJeremy L Thompson // Create the QFunction that builds the mass + diffusion operator (i.e. computes its quadrature data) and set its context data. 1710b96b02dSJeremy L Thompson CeedQFunction qf_build; 1720b96b02dSJeremy L Thompson 1730b96b02dSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, build_mass_diff, build_mass_diff_loc, &qf_build); 1740b96b02dSJeremy L Thompson CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 1750b96b02dSJeremy L Thompson CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT); 1760b96b02dSJeremy L Thompson CeedQFunctionAddOutput(qf_build, "qdata", 1 + dim * (dim + 1) / 2, CEED_EVAL_NONE); 1770b96b02dSJeremy L Thompson CeedQFunctionSetContext(qf_build, build_ctx); 1780b96b02dSJeremy L Thompson 1790b96b02dSJeremy L Thompson // Create the operator that builds the quadrature data for the mass + diffusion operator. 1800b96b02dSJeremy L Thompson CeedOperator op_build; 1810b96b02dSJeremy L Thompson 1820b96b02dSJeremy L Thompson CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build); 1830b96b02dSJeremy L Thompson CeedOperatorSetField(op_build, "dx", mesh_restriction, mesh_basis, CEED_VECTOR_ACTIVE); 1840b96b02dSJeremy L Thompson CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE); 1850b96b02dSJeremy L Thompson CeedOperatorSetField(op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 1860b96b02dSJeremy L Thompson 1870b96b02dSJeremy L Thompson // Compute the quadrature data for the mass + diffusion operator. 1880b96b02dSJeremy L Thompson CeedVector q_data; 1890b96b02dSJeremy L Thompson CeedInt elem_qpts = CeedIntPow(num_qpts, dim); 1900b96b02dSJeremy L Thompson CeedInt num_elem = 1; 1910b96b02dSJeremy L Thompson 1920b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d]; 1930b96b02dSJeremy L Thompson CeedVectorCreate(ceed, num_elem * elem_qpts * (1 + dim * (dim + 1) / 2), &q_data); 1940b96b02dSJeremy L Thompson CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE); 1950b96b02dSJeremy L Thompson 1960b96b02dSJeremy L Thompson // Create the QFunction that defines the action of the mass + diffusion operator. 1970b96b02dSJeremy L Thompson CeedQFunction qf_apply; 1980b96b02dSJeremy L Thompson 1990b96b02dSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, apply_mass_diff, apply_mass_diff_loc, &qf_apply); 2000b96b02dSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP); 2010b96b02dSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "du", dim, CEED_EVAL_GRAD); 2020b96b02dSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "qdata", 1 + dim * (dim + 1) / 2, CEED_EVAL_NONE); 2030b96b02dSJeremy L Thompson CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP); 2040b96b02dSJeremy L Thompson CeedQFunctionAddOutput(qf_apply, "dv", dim, CEED_EVAL_GRAD); 2050b96b02dSJeremy L Thompson CeedQFunctionSetContext(qf_apply, build_ctx); 2060b96b02dSJeremy L Thompson 2070b96b02dSJeremy L Thompson // Create the mass + diffusion operator. 2080b96b02dSJeremy L Thompson CeedOperator op_apply; 2090b96b02dSJeremy L Thompson 2100b96b02dSJeremy L Thompson CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply); 2110b96b02dSJeremy L Thompson CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 2120b96b02dSJeremy L Thompson CeedOperatorSetField(op_apply, "du", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 2130b96b02dSJeremy L Thompson CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data); 2140b96b02dSJeremy L Thompson CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 2150b96b02dSJeremy L Thompson CeedOperatorSetField(op_apply, "dv", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 2160b96b02dSJeremy L Thompson 2170b96b02dSJeremy L Thompson // Create auxiliary solution-size vectors. 2180b96b02dSJeremy L Thompson CeedVector u, v; 2190b96b02dSJeremy L Thompson 2200b96b02dSJeremy L Thompson CeedVectorCreate(ceed, sol_size, &u); 2210b96b02dSJeremy L Thompson CeedVectorCreate(ceed, sol_size, &v); 2220b96b02dSJeremy L Thompson 2230b96b02dSJeremy L Thompson // Initialize 'u' with ones. 2240b96b02dSJeremy L Thompson CeedVectorSetValue(u, 1.0); 2250b96b02dSJeremy L Thompson 2260b96b02dSJeremy L Thompson // Compute the mesh volume using the mass + diffusion operator: volume = 1^T \cdot M \cdot 1 2270b96b02dSJeremy L Thompson CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE); 2280b96b02dSJeremy L Thompson 2290b96b02dSJeremy L Thompson // Benchmark runs 2300b96b02dSJeremy L Thompson if (!test && benchmark) { 2310b96b02dSJeremy L Thompson // LCOV_EXCL_START 2320b96b02dSJeremy L Thompson printf(" Executing %d benchmarking runs...\n", benchmark); 2330b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 2340b96b02dSJeremy L Thompson } 2350b96b02dSJeremy L Thompson for (CeedInt i = 0; i < benchmark; i++) { 2360b96b02dSJeremy L Thompson // LCOV_EXCL_START 2370b96b02dSJeremy L Thompson CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE); 2380b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 2390b96b02dSJeremy L Thompson } 2400b96b02dSJeremy L Thompson 2410b96b02dSJeremy L Thompson // Compute and print the sum of the entries of 'v' giving the mesh volume. 2420b96b02dSJeremy L Thompson CeedScalar volume = 0.; 2430b96b02dSJeremy L Thompson 2440b96b02dSJeremy L Thompson { 2450b96b02dSJeremy L Thompson const CeedScalar *v_array; 2460b96b02dSJeremy L Thompson 2470b96b02dSJeremy L Thompson CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array); 2480b96b02dSJeremy L Thompson for (CeedInt i = 0; i < sol_size; i++) volume += v_array[i]; 2490b96b02dSJeremy L Thompson CeedVectorRestoreArrayRead(v, &v_array); 2500b96b02dSJeremy L Thompson } 2510b96b02dSJeremy L Thompson if (!test) { 2520b96b02dSJeremy L Thompson // LCOV_EXCL_START 2530b96b02dSJeremy L Thompson printf(" done.\n"); 2540b96b02dSJeremy L Thompson printf("Exact mesh volume : % .14g\n", exact_volume); 2550b96b02dSJeremy L Thompson printf("Computed mesh volume : % .14g\n", volume); 2560b96b02dSJeremy L Thompson printf("Volume error : % .14g\n", volume - exact_volume); 2570b96b02dSJeremy L Thompson // LCOV_EXCL_STOP 2580b96b02dSJeremy L Thompson } else { 2590b96b02dSJeremy L Thompson CeedScalar tol = (dim == 1 ? 200. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5); 2600b96b02dSJeremy L Thompson 2610b96b02dSJeremy L Thompson if (fabs(volume - exact_volume) > tol) printf("Volume error : % .1e\n", volume - exact_volume); 2620b96b02dSJeremy L Thompson } 2630b96b02dSJeremy L Thompson 2640b96b02dSJeremy L Thompson // Free dynamically allocated memory. 2650b96b02dSJeremy L Thompson CeedVectorDestroy(&u); 2660b96b02dSJeremy L Thompson CeedVectorDestroy(&v); 2670b96b02dSJeremy L Thompson CeedVectorDestroy(&q_data); 2680b96b02dSJeremy L Thompson CeedVectorDestroy(&mesh_coords); 2690b96b02dSJeremy L Thompson CeedOperatorDestroy(&op_apply); 2700b96b02dSJeremy L Thompson CeedQFunctionDestroy(&qf_apply); 2710b96b02dSJeremy L Thompson CeedQFunctionContextDestroy(&build_ctx); 2720b96b02dSJeremy L Thompson CeedOperatorDestroy(&op_build); 2730b96b02dSJeremy L Thompson CeedQFunctionDestroy(&qf_build); 2740b96b02dSJeremy L Thompson CeedElemRestrictionDestroy(&sol_restriction); 2750b96b02dSJeremy L Thompson CeedElemRestrictionDestroy(&mesh_restriction); 2760b96b02dSJeremy L Thompson CeedElemRestrictionDestroy(&q_data_restriction); 2770b96b02dSJeremy L Thompson CeedBasisDestroy(&sol_basis); 2780b96b02dSJeremy L Thompson CeedBasisDestroy(&mesh_basis); 2790b96b02dSJeremy L Thompson CeedDestroy(&ceed); 2800b96b02dSJeremy L Thompson return 0; 2810b96b02dSJeremy L Thompson } 2820b96b02dSJeremy L Thompson 2830b96b02dSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) { 2840b96b02dSJeremy L Thompson // Use the approximate formula: 2850b96b02dSJeremy L Thompson // prob_size ~ num_elem * degree^dim 2860b96b02dSJeremy L Thompson CeedInt num_elem = prob_size / CeedIntPow(degree, dim); 2870b96b02dSJeremy L Thompson CeedInt s = 0; // find s: num_elem/2 < 2^s <= num_elem 2880b96b02dSJeremy L Thompson 2890b96b02dSJeremy L Thompson while (num_elem > 1) { 2900b96b02dSJeremy L Thompson num_elem /= 2; 2910b96b02dSJeremy L Thompson s++; 2920b96b02dSJeremy L Thompson } 2930b96b02dSJeremy L Thompson CeedInt r = s % dim; 2940b96b02dSJeremy L Thompson 2950b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 2960b96b02dSJeremy L Thompson CeedInt sd = s / dim; 2970b96b02dSJeremy L Thompson 2980b96b02dSJeremy L Thompson if (r > 0) { 2990b96b02dSJeremy L Thompson sd++; 3000b96b02dSJeremy L Thompson r--; 3010b96b02dSJeremy L Thompson } 3020b96b02dSJeremy L Thompson num_xyz[d] = 1 << sd; 3030b96b02dSJeremy L Thompson } 3040b96b02dSJeremy L Thompson return 0; 3050b96b02dSJeremy L Thompson } 3060b96b02dSJeremy L Thompson 3070b96b02dSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts, 3080b96b02dSJeremy L Thompson CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction) { 3090b96b02dSJeremy L Thompson CeedInt p = degree + 1; 3100b96b02dSJeremy L Thompson CeedInt num_nodes = CeedIntPow(p, dim); // number of scalar nodes per element 3110b96b02dSJeremy L Thompson CeedInt elem_qpts = CeedIntPow(num_qpts, dim); // number of qpts per element 3120b96b02dSJeremy L Thompson CeedInt nd[3], num_elem = 1, scalar_size = 1; 3130b96b02dSJeremy L Thompson 3140b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 3150b96b02dSJeremy L Thompson num_elem *= num_xyz[d]; 3160b96b02dSJeremy L Thompson nd[d] = num_xyz[d] * (p - 1) + 1; 3170b96b02dSJeremy L Thompson scalar_size *= nd[d]; 3180b96b02dSJeremy L Thompson } 3190b96b02dSJeremy L Thompson *size = scalar_size * num_comp; 3200b96b02dSJeremy L Thompson // elem: 0 1 n-1 3210b96b02dSJeremy L Thompson // |---*-...-*---|---*-...-*---|- ... -|--...--| 3220b96b02dSJeremy L Thompson // num_nodes: 0 1 p-1 p p+1 2*p n*p 3230b96b02dSJeremy L Thompson CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes); 3240b96b02dSJeremy L Thompson 3250b96b02dSJeremy L Thompson for (CeedInt e = 0; e < num_elem; e++) { 3260b96b02dSJeremy L Thompson CeedInt e_xyz[3] = {1, 1, 1}, re = e; 3270b96b02dSJeremy L Thompson 3280b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 3290b96b02dSJeremy L Thompson e_xyz[d] = re % num_xyz[d]; 3300b96b02dSJeremy L Thompson re /= num_xyz[d]; 3310b96b02dSJeremy L Thompson } 3320b96b02dSJeremy L Thompson CeedInt *local_elem_nodes = elem_nodes + e * num_nodes; 3330b96b02dSJeremy L Thompson 3340b96b02dSJeremy L Thompson for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) { 3350b96b02dSJeremy L Thompson CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes; 3360b96b02dSJeremy L Thompson 3370b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 3380b96b02dSJeremy L Thompson g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride; 3390b96b02dSJeremy L Thompson g_nodes_stride *= nd[d]; 3400b96b02dSJeremy L Thompson r_nodes /= p; 3410b96b02dSJeremy L Thompson } 3420b96b02dSJeremy L Thompson local_elem_nodes[l_nodes] = g_nodes; 3430b96b02dSJeremy L Thompson } 3440b96b02dSJeremy L Thompson } 3450b96b02dSJeremy L Thompson if (restriction) { 3460b96b02dSJeremy L Thompson CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes, 3470b96b02dSJeremy L Thompson restriction); 3480b96b02dSJeremy L Thompson } 3490b96b02dSJeremy L Thompson if (q_data_restriction) { 3500b96b02dSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, q_data_restriction); 3510b96b02dSJeremy L Thompson } 3520b96b02dSJeremy L Thompson free(elem_nodes); 3530b96b02dSJeremy L Thompson return 0; 3540b96b02dSJeremy L Thompson } 3550b96b02dSJeremy L Thompson 3560b96b02dSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) { 3570b96b02dSJeremy L Thompson CeedInt p = mesh_degree + 1; 3580b96b02dSJeremy L Thompson CeedInt nd[3], scalar_size = 1; 3590b96b02dSJeremy L Thompson 3600b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 3610b96b02dSJeremy L Thompson nd[d] = num_xyz[d] * (p - 1) + 1; 3620b96b02dSJeremy L Thompson scalar_size *= nd[d]; 3630b96b02dSJeremy L Thompson } 3640b96b02dSJeremy L Thompson CeedScalar *coords; 3650b96b02dSJeremy L Thompson 3660b96b02dSJeremy L Thompson CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords); 3670b96b02dSJeremy L Thompson CeedScalar *nodes = malloc(sizeof(CeedScalar) * p); 3680b96b02dSJeremy L Thompson 3690b96b02dSJeremy L Thompson // The H1 basis uses Lobatto quadrature points as nodes. 3700b96b02dSJeremy L Thompson CeedLobattoQuadrature(p, nodes, NULL); // nodes are in [-1,1] 3710b96b02dSJeremy L Thompson for (CeedInt i = 0; i < p; i++) nodes[i] = 0.5 + 0.5 * nodes[i]; 3720b96b02dSJeremy L Thompson for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) { 3730b96b02dSJeremy L Thompson CeedInt r_nodes = gs_nodes; 3740b96b02dSJeremy L Thompson 3750b96b02dSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 3760b96b02dSJeremy L Thompson CeedInt d_1d = r_nodes % nd[d]; 3770b96b02dSJeremy L Thompson coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d]; 3780b96b02dSJeremy L Thompson r_nodes /= nd[d]; 3790b96b02dSJeremy L Thompson } 3800b96b02dSJeremy L Thompson } 3810b96b02dSJeremy L Thompson free(nodes); 3820b96b02dSJeremy L Thompson CeedVectorRestoreArray(mesh_coords, &coords); 3830b96b02dSJeremy L Thompson return 0; 3840b96b02dSJeremy L Thompson } 3850b96b02dSJeremy L Thompson 3860b96b02dSJeremy L Thompson #ifndef M_PI 3870b96b02dSJeremy L Thompson #define M_PI 3.14159265358979323846 3880b96b02dSJeremy L Thompson #define M_PI_2 1.57079632679489661923 3890b96b02dSJeremy L Thompson #endif 3900b96b02dSJeremy L Thompson 3910b96b02dSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) { 3920b96b02dSJeremy L Thompson CeedScalar exact_volume; 3930b96b02dSJeremy L Thompson CeedScalar *coords; 3940b96b02dSJeremy L Thompson 3950b96b02dSJeremy L Thompson CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords); 3960b96b02dSJeremy L Thompson if (dim == 1) { 3970b96b02dSJeremy L Thompson for (CeedInt i = 0; i < mesh_size; i++) { 3980b96b02dSJeremy L Thompson // map [0,1] to [0,1] varying the mesh density 3990b96b02dSJeremy L Thompson coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5)); 4000b96b02dSJeremy L Thompson } 4010b96b02dSJeremy L Thompson exact_volume = 1.; 4020b96b02dSJeremy L Thompson } else { 4030b96b02dSJeremy L Thompson CeedInt num_nodes = mesh_size / dim; 4040b96b02dSJeremy L Thompson for (CeedInt i = 0; i < num_nodes; i++) { 4050b96b02dSJeremy L Thompson // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar 4060b96b02dSJeremy L Thompson // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi 4070b96b02dSJeremy L Thompson CeedScalar u = coords[i], v = coords[i + num_nodes]; 4080b96b02dSJeremy L Thompson 4090b96b02dSJeremy L Thompson u = 1. + u; 4100b96b02dSJeremy L Thompson v = M_PI_2 * v; 4110b96b02dSJeremy L Thompson coords[i] = u * cos(v); 4120b96b02dSJeremy L Thompson coords[i + num_nodes] = u * sin(v); 4130b96b02dSJeremy L Thompson } 4140b96b02dSJeremy L Thompson exact_volume = 3. / 4. * M_PI; 4150b96b02dSJeremy L Thompson } 4160b96b02dSJeremy L Thompson CeedVectorRestoreArray(mesh_coords, &coords); 4170b96b02dSJeremy L Thompson return exact_volume; 4180b96b02dSJeremy L Thompson } 419