13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 366087c08SValeria Barra // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 566087c08SValeria Barra // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 766087c08SValeria Barra 866087c08SValeria Barra // libCEED Example 1 966087c08SValeria Barra // 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 29dc8efd83SLeila Ghaffari //TESTARGS(name="1D_User_QFunction") -ceed {ceed_resource} -d 1 -t 30dc8efd83SLeila Ghaffari //TESTARGS(name="2D_User_QFunction") -ceed {ceed_resource} -d 2 -t 31dc8efd83SLeila Ghaffari //TESTARGS(name="3D_User_QFunction") -ceed {ceed_resource} -d 3 -t 32dc8efd83SLeila Ghaffari //TESTARGS(name="1D_Gallery_QFunction") -ceed {ceed_resource} -d 1 -t -g 33dc8efd83SLeila Ghaffari //TESTARGS(name="2D_Gallery_QFunction") -ceed {ceed_resource} -d 2 -t -g 34dc8efd83SLeila Ghaffari //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> 43*49aac155SJeremy 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 63990fdeb6SJeremy L Thompson CeedInt help = 0, test = 0, gallery = 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; 8466087c08SValeria Barra } else if (!strcmp(argv[ia], "-t")) { 8566087c08SValeria Barra test = 1; 8666087c08SValeria Barra } else if (!strcmp(argv[ia], "-g")) { 8766087c08SValeria Barra gallery = 1; 8866087c08SValeria Barra } 8966087c08SValeria Barra if (parse_error) { 9066087c08SValeria Barra printf("Error parsing command line options.\n"); 9166087c08SValeria Barra return 1; 9266087c08SValeria Barra } 93ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 9466087c08SValeria Barra } 9566087c08SValeria Barra if (prob_size < 0) prob_size = test ? 8 * 16 : 256 * 1024; 9666087c08SValeria Barra 9766087c08SValeria Barra // Print the values of all options: 9866087c08SValeria Barra if (!test || help) { 99ded9b81dSJeremy L Thompson // LCOV_EXCL_START 10066087c08SValeria Barra printf("Selected options: [command line option] : <current value>\n"); 10166087c08SValeria Barra printf(" Ceed specification [-c] : %s\n", ceed_spec); 102990fdeb6SJeremy L Thompson printf(" Mesh dimension [-d] : %" CeedInt_FMT "\n", dim); 103990fdeb6SJeremy L Thompson printf(" Mesh degree [-m] : %" CeedInt_FMT "\n", mesh_degree); 104990fdeb6SJeremy L Thompson printf(" Solution degree [-p] : %" CeedInt_FMT "\n", sol_degree); 105d37d859eSJeremy L Thompson printf(" Num. 1D quadrature pts [-q] : %" CeedInt_FMT "\n", num_qpts); 106990fdeb6SJeremy L Thompson printf(" Approx. # unknowns [-s] : %" CeedInt_FMT "\n", prob_size); 10766087c08SValeria Barra printf(" QFunction source [-g] : %s\n", gallery ? "gallery" : "header"); 10866087c08SValeria Barra if (help) { 10966087c08SValeria Barra printf("Test/quiet mode is %s\n", (test ? "ON" : "OFF (use -t to enable)")); 11066087c08SValeria Barra return 0; 11166087c08SValeria Barra } 11266087c08SValeria Barra printf("\n"); 113ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 11466087c08SValeria Barra } 11566087c08SValeria Barra 116ea61e9acSJeremy L Thompson // Select appropriate backend and logical device based on the (-ceed) command line argument. 11766087c08SValeria Barra Ceed ceed; 11866087c08SValeria Barra CeedInit(ceed_spec, &ceed); 11966087c08SValeria Barra 12066087c08SValeria Barra // Construct the mesh and solution bases. 12166087c08SValeria Barra CeedBasis mesh_basis, sol_basis; 1222b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis); 1232b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis); 12466087c08SValeria Barra 12566087c08SValeria Barra // Determine the mesh size based on the given approximate problem size. 126990fdeb6SJeremy L Thompson CeedInt num_xyz[dim]; 127d1d35e2fSjeremylt GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz); 12866087c08SValeria Barra if (!test) { 129ded9b81dSJeremy L Thompson // LCOV_EXCL_START 130990fdeb6SJeremy L Thompson printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]); 1312b730f8bSJeremy L Thompson if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]); 1322b730f8bSJeremy L Thompson if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]); 13366087c08SValeria Barra printf("\n"); 134ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 13566087c08SValeria Barra } 13666087c08SValeria Barra 137ea61e9acSJeremy L Thompson // Build CeedElemRestriction objects describing the mesh and solution discrete representations. 13866087c08SValeria Barra CeedInt mesh_size, sol_size; 139d37d859eSJeremy L Thompson CeedElemRestriction mesh_restriction, sol_restriction, q_data_restriction; 140d37d859eSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restriction, NULL); 141d37d859eSJeremy L Thompson BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restriction, &q_data_restriction); 14266087c08SValeria Barra if (!test) { 143e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 144990fdeb6SJeremy L Thompson printf("Number of mesh nodes : %" CeedInt_FMT "\n", mesh_size / dim); 145990fdeb6SJeremy L Thompson printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size); 146e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 14766087c08SValeria Barra } 14866087c08SValeria Barra 14966087c08SValeria Barra // Create a CeedVector with the mesh coordinates. 15066087c08SValeria Barra CeedVector mesh_coords; 15166087c08SValeria Barra CeedVectorCreate(ceed, mesh_size, &mesh_coords); 152d1d35e2fSjeremylt SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords); 15366087c08SValeria Barra 15466087c08SValeria Barra // Apply a transformation to the mesh. 155d37d859eSJeremy L Thompson CeedScalar exact_volume = TransformMeshCoords(dim, mesh_size, mesh_coords); 15666087c08SValeria Barra 157d37d859eSJeremy L Thompson // Context data to be passed to the 'build_mass' QFunction. 158777ff853SJeremy L Thompson CeedQFunctionContext build_ctx; 159777ff853SJeremy L Thompson struct BuildContext build_ctx_data; 160777ff853SJeremy L Thompson build_ctx_data.dim = build_ctx_data.space_dim = dim; 161777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &build_ctx); 1622b730f8bSJeremy L Thompson CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data); 16366087c08SValeria Barra 164ea61e9acSJeremy L Thompson // Create the QFunction that builds the mass operator (i.e. computes its quadrature data) and set its context data. 165d1d35e2fSjeremylt CeedQFunction qf_build; 166d37d859eSJeremy L Thompson if (gallery) { 16766087c08SValeria Barra // This creates the QFunction via the gallery. 16866087c08SValeria Barra char name[13] = ""; 169990fdeb6SJeremy L Thompson snprintf(name, sizeof name, "Mass%" CeedInt_FMT "DBuild", dim); 170d1d35e2fSjeremylt CeedQFunctionCreateInteriorByName(ceed, name, &qf_build); 171d37d859eSJeremy L Thompson } else { 172d37d859eSJeremy L Thompson // This creates the QFunction directly. 173d37d859eSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, build_mass, build_mass_loc, &qf_build); 174d37d859eSJeremy L Thompson CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 175d37d859eSJeremy L Thompson CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT); 176d37d859eSJeremy L Thompson CeedQFunctionAddOutput(qf_build, "qdata", 1, CEED_EVAL_NONE); 177d37d859eSJeremy L Thompson CeedQFunctionSetContext(qf_build, build_ctx); 17866087c08SValeria Barra } 17966087c08SValeria Barra 18066087c08SValeria Barra // Create the operator that builds the quadrature data for the mass operator. 181d1d35e2fSjeremylt CeedOperator op_build; 1822b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build); 183d37d859eSJeremy L Thompson CeedOperatorSetField(op_build, "dx", mesh_restriction, mesh_basis, CEED_VECTOR_ACTIVE); 1842b730f8bSJeremy L Thompson CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE); 185d37d859eSJeremy L Thompson CeedOperatorSetField(op_build, "qdata", q_data_restriction, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 18666087c08SValeria Barra 18766087c08SValeria Barra // Compute the quadrature data for the mass operator. 188d1d35e2fSjeremylt CeedVector q_data; 18966087c08SValeria Barra CeedInt elem_qpts = CeedIntPow(num_qpts, dim); 19066087c08SValeria Barra CeedInt num_elem = 1; 1912b730f8bSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d]; 192d1d35e2fSjeremylt CeedVectorCreate(ceed, num_elem * elem_qpts, &q_data); 1932b730f8bSJeremy L Thompson CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE); 19466087c08SValeria Barra 195ded9b81dSJeremy L Thompson // Create the QFunction that defines the action of the mass operator. 196d1d35e2fSjeremylt CeedQFunction qf_apply; 197d37d859eSJeremy L Thompson if (gallery) { 198d37d859eSJeremy L Thompson // This creates the QFunction via the gallery. 199d37d859eSJeremy L Thompson CeedQFunctionCreateInteriorByName(ceed, "MassApply", &qf_apply); 200d37d859eSJeremy L Thompson } else { 20166087c08SValeria Barra // This creates the QFunction directly. 202d37d859eSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, apply_mass, apply_mass_loc, &qf_apply); 203d1d35e2fSjeremylt CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP); 204d1d35e2fSjeremylt CeedQFunctionAddInput(qf_apply, "qdata", 1, CEED_EVAL_NONE); 205d1d35e2fSjeremylt CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP); 20666087c08SValeria Barra } 20766087c08SValeria Barra 20866087c08SValeria Barra // Create the mass operator. 209d1d35e2fSjeremylt CeedOperator op_apply; 2102b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply); 211d37d859eSJeremy L Thompson CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 212d37d859eSJeremy L Thompson CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_COLLOCATED, q_data); 213d37d859eSJeremy L Thompson CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 21466087c08SValeria Barra 21566087c08SValeria Barra // Create auxiliary solution-size vectors. 21666087c08SValeria Barra CeedVector u, v; 21766087c08SValeria Barra CeedVectorCreate(ceed, sol_size, &u); 21866087c08SValeria Barra CeedVectorCreate(ceed, sol_size, &v); 21966087c08SValeria Barra 220b8bc0c60SPeter Munch // Initialize 'u' with ones. 22166087c08SValeria Barra CeedVectorSetValue(u, 1.0); 22266087c08SValeria Barra 223d37d859eSJeremy L Thompson // Compute the mesh volume using the mass operator: volume = 1^T \cdot M \cdot 1 224d1d35e2fSjeremylt CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE); 22566087c08SValeria Barra 22666087c08SValeria Barra // Compute and print the sum of the entries of 'v' giving the mesh volume. 227d37d859eSJeremy L Thompson CeedScalar volume = 0.; 228d37d859eSJeremy L Thompson { 229d1d35e2fSjeremylt const CeedScalar *v_array; 230d1d35e2fSjeremylt CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array); 231d37d859eSJeremy L Thompson for (CeedInt i = 0; i < sol_size; i++) volume += v_array[i]; 232d1d35e2fSjeremylt CeedVectorRestoreArrayRead(v, &v_array); 233d37d859eSJeremy L Thompson } 23466087c08SValeria Barra if (!test) { 235ded9b81dSJeremy L Thompson // LCOV_EXCL_START 23666087c08SValeria Barra printf(" done.\n"); 237d37d859eSJeremy L Thompson printf("Exact mesh volume : % .14g\n", exact_volume); 238d37d859eSJeremy L Thompson printf("Computed mesh volume : % .14g\n", volume); 239d37d859eSJeremy L Thompson printf("Volume error : % .14g\n", volume - exact_volume); 240ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 24166087c08SValeria Barra } else { 24280a9ef05SNatalie Beams CeedScalar tol = (dim == 1 ? 100. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5); 243d37d859eSJeremy L Thompson if (fabs(volume - exact_volume) > tol) printf("Volume error : % .1e\n", volume - exact_volume); 24466087c08SValeria Barra } 24566087c08SValeria Barra 24666087c08SValeria Barra // Free dynamically allocated memory. 24766087c08SValeria Barra CeedVectorDestroy(&u); 24866087c08SValeria Barra CeedVectorDestroy(&v); 249d1d35e2fSjeremylt CeedVectorDestroy(&q_data); 25066087c08SValeria Barra CeedVectorDestroy(&mesh_coords); 251d1d35e2fSjeremylt CeedOperatorDestroy(&op_apply); 252d1d35e2fSjeremylt CeedQFunctionDestroy(&qf_apply); 253777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&build_ctx); 254d1d35e2fSjeremylt CeedOperatorDestroy(&op_build); 255d1d35e2fSjeremylt CeedQFunctionDestroy(&qf_build); 256d37d859eSJeremy L Thompson CeedElemRestrictionDestroy(&sol_restriction); 257d37d859eSJeremy L Thompson CeedElemRestrictionDestroy(&mesh_restriction); 258d37d859eSJeremy L Thompson CeedElemRestrictionDestroy(&q_data_restriction); 25966087c08SValeria Barra CeedBasisDestroy(&sol_basis); 26066087c08SValeria Barra CeedBasisDestroy(&mesh_basis); 26166087c08SValeria Barra CeedDestroy(&ceed); 26266087c08SValeria Barra return 0; 26366087c08SValeria Barra } 26466087c08SValeria Barra 2652b730f8bSJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) { 26666087c08SValeria Barra // Use the approximate formula: 267ded9b81dSJeremy L Thompson // prob_size ~ num_elem * degree^dim 268ded9b81dSJeremy L Thompson CeedInt num_elem = prob_size / CeedIntPow(degree, dim); 26966087c08SValeria Barra CeedInt s = 0; // find s: num_elem/2 < 2^s <= num_elem 27066087c08SValeria Barra while (num_elem > 1) { 27166087c08SValeria Barra num_elem /= 2; 27266087c08SValeria Barra s++; 27366087c08SValeria Barra } 27466087c08SValeria Barra CeedInt r = s % dim; 275990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 276990fdeb6SJeremy L Thompson CeedInt sd = s / dim; 2772b730f8bSJeremy L Thompson if (r > 0) { 2782b730f8bSJeremy L Thompson sd++; 2792b730f8bSJeremy L Thompson r--; 2802b730f8bSJeremy L Thompson } 281d1d35e2fSjeremylt num_xyz[d] = 1 << sd; 28266087c08SValeria Barra } 28366087c08SValeria Barra return 0; 28466087c08SValeria Barra } 28566087c08SValeria Barra 2862b730f8bSJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts, 287d37d859eSJeremy L Thompson CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction) { 288ded9b81dSJeremy L Thompson CeedInt p = degree + 1; 289d1d35e2fSjeremylt CeedInt num_nodes = CeedIntPow(p, dim); // number of scalar nodes per element 29066087c08SValeria Barra CeedInt elem_qpts = CeedIntPow(num_qpts, dim); // number of qpts per element 29166087c08SValeria Barra CeedInt nd[3], num_elem = 1, scalar_size = 1; 292990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 293d1d35e2fSjeremylt num_elem *= num_xyz[d]; 294d1d35e2fSjeremylt nd[d] = num_xyz[d] * (p - 1) + 1; 29566087c08SValeria Barra scalar_size *= nd[d]; 29666087c08SValeria Barra } 297d1d35e2fSjeremylt *size = scalar_size * num_comp; 29866087c08SValeria Barra // elem: 0 1 n-1 29966087c08SValeria Barra // |---*-...-*---|---*-...-*---|- ... -|--...--| 300d1d35e2fSjeremylt // num_nodes: 0 1 p-1 p p+1 2*p n*p 301d1d35e2fSjeremylt CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes); 30266087c08SValeria Barra for (CeedInt e = 0; e < num_elem; e++) { 303d1d35e2fSjeremylt CeedInt e_xyz[3] = {1, 1, 1}, re = e; 3042b730f8bSJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 3052b730f8bSJeremy L Thompson e_xyz[d] = re % num_xyz[d]; 3062b730f8bSJeremy L Thompson re /= num_xyz[d]; 3072b730f8bSJeremy L Thompson } 308d37d859eSJeremy L Thompson CeedInt *local_elem_nodes = elem_nodes + e * num_nodes; 309990fdeb6SJeremy L Thompson for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) { 310d1d35e2fSjeremylt CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes; 311990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 312d1d35e2fSjeremylt g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride; 313d1d35e2fSjeremylt g_nodes_stride *= nd[d]; 314d1d35e2fSjeremylt r_nodes /= p; 31566087c08SValeria Barra } 316d37d859eSJeremy L Thompson local_elem_nodes[l_nodes] = g_nodes; 31766087c08SValeria Barra } 31866087c08SValeria Barra } 3192b730f8bSJeremy L Thompson CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes, 320d37d859eSJeremy L Thompson restriction); 321d37d859eSJeremy L Thompson if (q_data_restriction) 322d37d859eSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, q_data_restriction); 323d1d35e2fSjeremylt free(elem_nodes); 32466087c08SValeria Barra return 0; 32566087c08SValeria Barra } 32666087c08SValeria Barra 3272b730f8bSJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) { 328ded9b81dSJeremy L Thompson CeedInt p = mesh_degree + 1; 32966087c08SValeria Barra CeedInt nd[3], num_elem = 1, scalar_size = 1; 330990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 331d1d35e2fSjeremylt num_elem *= num_xyz[d]; 332d1d35e2fSjeremylt nd[d] = num_xyz[d] * (p - 1) + 1; 33366087c08SValeria Barra scalar_size *= nd[d]; 33466087c08SValeria Barra } 33566087c08SValeria Barra CeedScalar *coords; 3369c774eddSJeremy L Thompson CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords); 337ded9b81dSJeremy L Thompson CeedScalar *nodes = malloc(sizeof(CeedScalar) * p); 33866087c08SValeria Barra // The H1 basis uses Lobatto quadrature points as nodes. 339ded9b81dSJeremy L Thompson CeedLobattoQuadrature(p, nodes, NULL); // nodes are in [-1,1] 340d37d859eSJeremy L Thompson for (CeedInt i = 0; i < p; i++) nodes[i] = 0.5 + 0.5 * nodes[i]; 341d1d35e2fSjeremylt for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) { 342d1d35e2fSjeremylt CeedInt r_nodes = gs_nodes; 343990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 344d1d35e2fSjeremylt CeedInt d_1d = r_nodes % nd[d]; 3452b730f8bSJeremy L Thompson coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d]; 346d1d35e2fSjeremylt r_nodes /= nd[d]; 34766087c08SValeria Barra } 34866087c08SValeria Barra } 34966087c08SValeria Barra free(nodes); 35066087c08SValeria Barra CeedVectorRestoreArray(mesh_coords, &coords); 35166087c08SValeria Barra return 0; 35266087c08SValeria Barra } 35366087c08SValeria Barra 35466087c08SValeria Barra #ifndef M_PI 35566087c08SValeria Barra #define M_PI 3.14159265358979323846 35666087c08SValeria Barra #define M_PI_2 1.57079632679489661923 35766087c08SValeria Barra #endif 35866087c08SValeria Barra 3592b730f8bSJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) { 36066087c08SValeria Barra CeedScalar exact_volume; 36166087c08SValeria Barra CeedScalar *coords; 36266087c08SValeria Barra CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords); 36366087c08SValeria Barra if (dim == 1) { 36466087c08SValeria Barra for (CeedInt i = 0; i < mesh_size; i++) { 36566087c08SValeria Barra // map [0,1] to [0,1] varying the mesh density 36666087c08SValeria Barra coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5)); 36766087c08SValeria Barra } 36866087c08SValeria Barra exact_volume = 1.; 36966087c08SValeria Barra } else { 37066087c08SValeria Barra CeedInt num_nodes = mesh_size / dim; 37166087c08SValeria Barra for (CeedInt i = 0; i < num_nodes; i++) { 37266087c08SValeria Barra // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar 37366087c08SValeria Barra // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi 37466087c08SValeria Barra CeedScalar u = coords[i], v = coords[i + num_nodes]; 37566087c08SValeria Barra u = 1. + u; 37666087c08SValeria Barra v = M_PI_2 * v; 37766087c08SValeria Barra coords[i] = u * cos(v); 37866087c08SValeria Barra coords[i + num_nodes] = u * sin(v); 37966087c08SValeria Barra } 38066087c08SValeria Barra exact_volume = 3. / 4. * M_PI; 38166087c08SValeria Barra } 38266087c08SValeria Barra CeedVectorRestoreArray(mesh_coords, &coords); 38366087c08SValeria Barra return exact_volume; 38466087c08SValeria Barra } 385