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 2 966087c08SValeria Barra // 1066087c08SValeria Barra // This example illustrates a simple usage of libCEED to compute the surface 1166087c08SValeria Barra // area of a 3D body using matrix-free application of a diffusion operator. 12ded9b81dSJeremy L Thompson // Arbitrary mesh and solution degrees in 1D, 2D and 3D are supported from the 1366087c08SValeria Barra // same code. 1466087c08SValeria Barra // 1566087c08SValeria Barra // The example has no dependencies, and is designed to be self-contained. For 1666087c08SValeria Barra // additional examples that use external discretization libraries (MFEM, PETSc, 1766087c08SValeria Barra // etc.) see the subdirectories in libceed/examples. 1866087c08SValeria Barra // 1966087c08SValeria Barra // All libCEED objects use a Ceed device object constructed based on a command 2066087c08SValeria Barra // line argument (-ceed). 2166087c08SValeria Barra // 2266087c08SValeria Barra // Build with: 2366087c08SValeria Barra // 2466087c08SValeria Barra // make ex2-surface [CEED_DIR=</path/to/libceed>] 2566087c08SValeria Barra // 2666087c08SValeria Barra // Sample runs: 2766087c08SValeria Barra // 2866087c08SValeria Barra // ./ex2-surface 2966087c08SValeria Barra // ./ex2-surface -ceed /cpu/self 3028688798Sjeremylt // ./ex2-surface -ceed /gpu/cuda 3166087c08SValeria Barra // 3266087c08SValeria Barra // Next line is grep'd from tap.sh to set its arguments 3366087c08SValeria Barra // Test in 1D-3D 34dc8efd83SLeila Ghaffari //TESTARGS(name="1D_user_QFunction") -ceed {ceed_resource} -d 1 -t 35dc8efd83SLeila Ghaffari //TESTARGS(name="2D_user_QFunction") -ceed {ceed_resource} -d 2 -t 36dc8efd83SLeila Ghaffari //TESTARGS(name="3D_user_QFunction") -ceed {ceed_resource} -d 3 -t 37dc8efd83SLeila Ghaffari //TESTARGS(name="1D_Gallery_QFunction") -ceed {ceed_resource} -d 1 -t -g 38dc8efd83SLeila Ghaffari //TESTARGS(name="2D_Gallery_QFunction") -ceed {ceed_resource} -d 2 -t -g 39dc8efd83SLeila Ghaffari //TESTARGS(name="3D_Gallery_QFunction") -ceed {ceed_resource} -d 3 -t -g 4066087c08SValeria Barra 4166087c08SValeria Barra /// @file 4266087c08SValeria Barra /// libCEED example using diffusion operator to compute surface area 4366087c08SValeria Barra 4466087c08SValeria Barra #include <ceed.h> 4566087c08SValeria Barra #include <math.h> 463d576824SJeremy L Thompson #include <stdlib.h> 4766087c08SValeria Barra #include <string.h> 4866087c08SValeria Barra #include "ex2-surface.h" 4966087c08SValeria Barra 5066087c08SValeria Barra // Auxiliary functions. 51*990fdeb6SJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, 52*990fdeb6SJeremy L Thompson CeedInt num_xyz[3]); 53*990fdeb6SJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[3], 54*990fdeb6SJeremy L Thompson CeedInt degree, CeedInt num_comp, CeedInt *size, 55*990fdeb6SJeremy L Thompson CeedInt num_qpts, CeedElemRestriction *restr, 5666087c08SValeria Barra CeedElemRestriction *restr_i); 57*990fdeb6SJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[3], CeedInt mesh_degree, 5866087c08SValeria Barra CeedVector mesh_coords); 59*990fdeb6SJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, 60*990fdeb6SJeremy L Thompson CeedVector mesh_coords); 6166087c08SValeria Barra 6266087c08SValeria Barra int main(int argc, const char *argv[]) { 6366087c08SValeria Barra const char *ceed_spec = "/cpu/self"; 64*990fdeb6SJeremy L Thompson CeedInt dim = 3; // dimension of the mesh 65*990fdeb6SJeremy L Thompson CeedInt num_comp_x = 3; // number of x components 66*990fdeb6SJeremy L Thompson CeedInt mesh_degree = 4; // polynomial degree for the mesh 67*990fdeb6SJeremy L Thompson CeedInt sol_degree = 4; // polynomial degree for the solution 68*990fdeb6SJeremy L Thompson CeedInt num_qpts = sol_degree + 2; // number of 1D quadrature points 69*990fdeb6SJeremy L Thompson CeedInt prob_size = -1; // approximate problem size 70*990fdeb6SJeremy L Thompson CeedInt help = 0, test = 0, gallery = 0; 7166087c08SValeria Barra 7266087c08SValeria Barra // Process command line arguments. 7366087c08SValeria Barra for (int ia = 1; ia < argc; ia++) { 74ded9b81dSJeremy L Thompson // LCOV_EXCL_START 7566087c08SValeria Barra int next_arg = ((ia+1) < argc), parse_error = 0; 7666087c08SValeria Barra if (!strcmp(argv[ia],"-h")) { 7766087c08SValeria Barra help = 1; 7866087c08SValeria Barra } else if (!strcmp(argv[ia],"-c") || !strcmp(argv[ia],"-ceed")) { 7966087c08SValeria Barra parse_error = next_arg ? ceed_spec = argv[++ia], 0 : 1; 8066087c08SValeria Barra } else if (!strcmp(argv[ia],"-d")) { 8166087c08SValeria Barra parse_error = next_arg ? dim = atoi(argv[++ia]), 0 : 1; 82d1d35e2fSjeremylt num_comp_x = dim; 8366087c08SValeria Barra } else if (!strcmp(argv[ia],"-m")) { 84ded9b81dSJeremy L Thompson parse_error = next_arg ? mesh_degree = atoi(argv[++ia]), 0 : 1; 85ded9b81dSJeremy L Thompson } else if (!strcmp(argv[ia],"-p")) { 86ded9b81dSJeremy L Thompson parse_error = next_arg ? sol_degree = atoi(argv[++ia]), 0 : 1; 8766087c08SValeria Barra } else if (!strcmp(argv[ia],"-q")) { 8866087c08SValeria Barra parse_error = next_arg ? num_qpts = atoi(argv[++ia]), 0 : 1; 8966087c08SValeria Barra } else if (!strcmp(argv[ia],"-s")) { 9066087c08SValeria Barra parse_error = next_arg ? prob_size = atoi(argv[++ia]), 0 : 1; 9166087c08SValeria Barra } else if (!strcmp(argv[ia],"-t")) { 9266087c08SValeria Barra test = 1; 9366087c08SValeria Barra } else if (!strcmp(argv[ia],"-g")) { 9466087c08SValeria Barra gallery = 1; 9566087c08SValeria Barra } 9666087c08SValeria Barra if (parse_error) { 9766087c08SValeria Barra printf("Error parsing command line options.\n"); 9866087c08SValeria Barra return 1; 9966087c08SValeria Barra } 100ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 10166087c08SValeria Barra } 10266087c08SValeria Barra if (prob_size < 0) prob_size = test ? 16*16*dim*dim : 256*1024; 10366087c08SValeria Barra 104ded9b81dSJeremy L Thompson // Set mesh_degree = sol_degree. 105ded9b81dSJeremy L Thompson mesh_degree = fmax(mesh_degree, sol_degree); 106ded9b81dSJeremy L Thompson sol_degree = mesh_degree; 10766087c08SValeria Barra 10866087c08SValeria Barra // Print the values of all options: 10966087c08SValeria Barra if (!test || help) { 110ded9b81dSJeremy L Thompson // LCOV_EXCL_START 11166087c08SValeria Barra printf("Selected options: [command line option] : <current value>\n"); 11266087c08SValeria Barra printf(" Ceed specification [-c] : %s\n", ceed_spec); 113*990fdeb6SJeremy L Thompson printf(" Mesh dimension [-d] : %" CeedInt_FMT "\n", dim); 114*990fdeb6SJeremy L Thompson printf(" Mesh degree [-m] : %" CeedInt_FMT "\n", mesh_degree); 115*990fdeb6SJeremy L Thompson printf(" Solution degree [-p] : %" CeedInt_FMT "\n", sol_degree); 116*990fdeb6SJeremy L Thompson printf(" Num. 1D quadr. pts [-q] : %" CeedInt_FMT "\n", num_qpts); 117*990fdeb6SJeremy L Thompson printf(" Approx. # unknowns [-s] : %" CeedInt_FMT "\n", prob_size); 11866087c08SValeria Barra printf(" QFunction source [-g] : %s\n", gallery?"gallery":"header"); 11966087c08SValeria Barra if (help) { 12066087c08SValeria Barra printf("Test/quiet mode is %s\n", (test?"ON":"OFF (use -t to enable)")); 12166087c08SValeria Barra return 0; 12266087c08SValeria Barra } 12366087c08SValeria Barra printf("\n"); 124ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 12566087c08SValeria Barra } 12666087c08SValeria Barra 12766087c08SValeria Barra // Select appropriate backend and logical device based on the <ceed-spec> 12866087c08SValeria Barra // command line argument. 12966087c08SValeria Barra Ceed ceed; 13066087c08SValeria Barra CeedInit(ceed_spec, &ceed); 13166087c08SValeria Barra 13266087c08SValeria Barra // Construct the mesh and solution bases. 13366087c08SValeria Barra CeedBasis mesh_basis, sol_basis; 134d1d35e2fSjeremylt CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, 135d1d35e2fSjeremylt num_qpts, CEED_GAUSS, &mesh_basis); 136ded9b81dSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, 13766087c08SValeria Barra CEED_GAUSS, &sol_basis); 13866087c08SValeria Barra 13966087c08SValeria Barra // Determine the mesh size based on the given approximate problem size. 140*990fdeb6SJeremy L Thompson CeedInt num_xyz[3]; 141d1d35e2fSjeremylt GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz); 14266087c08SValeria Barra 14366087c08SValeria Barra if (!test) { 144ded9b81dSJeremy L Thompson // LCOV_EXCL_START 145*990fdeb6SJeremy L Thompson printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]); 146*990fdeb6SJeremy L Thompson if (dim > 1) { printf(", ny = %" CeedInt_FMT, num_xyz[1]); } 147*990fdeb6SJeremy L Thompson if (dim > 2) { printf(", nz = %" CeedInt_FMT, num_xyz[2]); } 14866087c08SValeria Barra printf("\n"); 149ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 15066087c08SValeria Barra } 15166087c08SValeria Barra 15266087c08SValeria Barra // Build CeedElemRestriction objects describing the mesh and solution discrete 15366087c08SValeria Barra // representations. 15466087c08SValeria Barra CeedInt mesh_size, sol_size; 155d1d35e2fSjeremylt CeedElemRestriction mesh_restr, sol_restr, q_data_restr_i; 156d1d35e2fSjeremylt BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, 157d1d35e2fSjeremylt &mesh_size, num_qpts, &mesh_restr, NULL); 158d1d35e2fSjeremylt BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, dim*(dim+1)/2, 159d1d35e2fSjeremylt &sol_size, num_qpts, NULL, &q_data_restr_i); 160d1d35e2fSjeremylt BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, 161ded9b81dSJeremy L Thompson num_qpts, &sol_restr, NULL); 16266087c08SValeria Barra if (!test) { 163ded9b81dSJeremy L Thompson // LCOV_EXCL_START 164*990fdeb6SJeremy L Thompson printf("Number of mesh nodes : %" CeedInt_FMT "\n", mesh_size/dim); 165*990fdeb6SJeremy L Thompson printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size); 166ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 16766087c08SValeria Barra } 16866087c08SValeria Barra 16966087c08SValeria Barra // Create a CeedVector with the mesh coordinates. 17066087c08SValeria Barra CeedVector mesh_coords; 17166087c08SValeria Barra CeedVectorCreate(ceed, mesh_size, &mesh_coords); 172d1d35e2fSjeremylt SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords); 17366087c08SValeria Barra 17466087c08SValeria Barra // Apply a transformation to the mesh. 17566087c08SValeria Barra CeedScalar exact_sa = TransformMeshCoords(dim, mesh_size, mesh_coords); 17666087c08SValeria Barra 177ded9b81dSJeremy L Thompson // Context data to be passed to the 'f_build_diff' QFunction. 178777ff853SJeremy L Thompson CeedQFunctionContext build_ctx; 179777ff853SJeremy L Thompson struct BuildContext build_ctx_data; 180777ff853SJeremy L Thompson build_ctx_data.dim = build_ctx_data.space_dim = dim; 181777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &build_ctx); 182777ff853SJeremy L Thompson CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, 183777ff853SJeremy L Thompson sizeof(build_ctx_data), &build_ctx_data); 18466087c08SValeria Barra 185ded9b81dSJeremy L Thompson // Create the QFunction that builds the diffusion operator (i.e. computes its 18666087c08SValeria Barra // quadrature data) and set its context data. 187d1d35e2fSjeremylt CeedQFunction qf_build; 18866087c08SValeria Barra switch (gallery) { 18966087c08SValeria Barra case 0: 19066087c08SValeria Barra // This creates the QFunction directly. 19166087c08SValeria Barra CeedQFunctionCreateInterior(ceed, 1, f_build_diff, 192d1d35e2fSjeremylt f_build_diff_loc, &qf_build); 193d1d35e2fSjeremylt CeedQFunctionAddInput(qf_build, "dx", num_comp_x*dim, CEED_EVAL_GRAD); 194d1d35e2fSjeremylt CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT); 195d1d35e2fSjeremylt CeedQFunctionAddOutput(qf_build, "qdata", dim*(dim+1)/2, CEED_EVAL_NONE); 196d1d35e2fSjeremylt CeedQFunctionSetContext(qf_build, build_ctx); 19766087c08SValeria Barra break; 19866087c08SValeria Barra case 1: { 19966087c08SValeria Barra // This creates the QFunction via the gallery. 20066087c08SValeria Barra char name[16] = ""; 201*990fdeb6SJeremy L Thompson snprintf(name, sizeof name, "Poisson%" CeedInt_FMT "DBuild", dim); 202d1d35e2fSjeremylt CeedQFunctionCreateInteriorByName(ceed, name, &qf_build); 20366087c08SValeria Barra break; 20466087c08SValeria Barra } 20566087c08SValeria Barra } 20666087c08SValeria Barra 20766087c08SValeria Barra // Create the operator that builds the quadrature data for the diffusion 20866087c08SValeria Barra // operator. 209d1d35e2fSjeremylt CeedOperator op_build; 210d1d35e2fSjeremylt CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, 211d1d35e2fSjeremylt CEED_QFUNCTION_NONE, &op_build); 212d1d35e2fSjeremylt CeedOperatorSetField(op_build, "dx", mesh_restr, mesh_basis, 213a8d32208Sjeremylt CEED_VECTOR_ACTIVE); 214d1d35e2fSjeremylt CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, 21515910d16Sjeremylt mesh_basis, CEED_VECTOR_NONE); 216d1d35e2fSjeremylt CeedOperatorSetField(op_build, "qdata", q_data_restr_i, 21766087c08SValeria Barra CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 21866087c08SValeria Barra 21966087c08SValeria Barra // Compute the quadrature data for the diffusion operator. 220d1d35e2fSjeremylt CeedVector q_data; 22166087c08SValeria Barra CeedInt elem_qpts = CeedIntPow(num_qpts, dim); 22266087c08SValeria Barra CeedInt num_elem = 1; 223*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) 224d1d35e2fSjeremylt num_elem *= num_xyz[d]; 225d1d35e2fSjeremylt CeedVectorCreate(ceed, num_elem*elem_qpts*dim*(dim+1)/2, &q_data); 226d1d35e2fSjeremylt CeedOperatorApply(op_build, mesh_coords, q_data, 22766087c08SValeria Barra CEED_REQUEST_IMMEDIATE); 22866087c08SValeria Barra 229ded9b81dSJeremy L Thompson // Create the QFunction that defines the action of the diffusion operator. 230d1d35e2fSjeremylt CeedQFunction qf_apply; 23166087c08SValeria Barra switch (gallery) { 23266087c08SValeria Barra case 0: 23366087c08SValeria Barra // This creates the QFunction directly. 23466087c08SValeria Barra CeedQFunctionCreateInterior(ceed, 1, f_apply_diff, 235d1d35e2fSjeremylt f_apply_diff_loc, &qf_apply); 236d1d35e2fSjeremylt CeedQFunctionAddInput(qf_apply, "du", dim, CEED_EVAL_GRAD); 237d1d35e2fSjeremylt CeedQFunctionAddInput(qf_apply, "qdata", dim*(dim+1)/2, CEED_EVAL_NONE); 238d1d35e2fSjeremylt CeedQFunctionAddOutput(qf_apply, "dv", dim, CEED_EVAL_GRAD); 239d1d35e2fSjeremylt CeedQFunctionSetContext(qf_apply, build_ctx); 24066087c08SValeria Barra break; 24166087c08SValeria Barra case 1: { 24266087c08SValeria Barra // This creates the QFunction via the gallery. 24366087c08SValeria Barra char name[16] = ""; 244*990fdeb6SJeremy L Thompson snprintf(name, sizeof name, "Poisson%" CeedInt_FMT "DApply", dim); 245d1d35e2fSjeremylt CeedQFunctionCreateInteriorByName(ceed, name, &qf_apply); 24666087c08SValeria Barra break; 24766087c08SValeria Barra } 24866087c08SValeria Barra } 24966087c08SValeria Barra 25066087c08SValeria Barra // Create the diffusion operator. 251d1d35e2fSjeremylt CeedOperator op_apply; 252d1d35e2fSjeremylt CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, 253d1d35e2fSjeremylt CEED_QFUNCTION_NONE, &op_apply); 254d1d35e2fSjeremylt CeedOperatorSetField(op_apply, "du", sol_restr, sol_basis, CEED_VECTOR_ACTIVE); 255d1d35e2fSjeremylt CeedOperatorSetField(op_apply, "qdata", q_data_restr_i, CEED_BASIS_COLLOCATED, 256d1d35e2fSjeremylt q_data); 257d1d35e2fSjeremylt CeedOperatorSetField(op_apply, "dv", sol_restr, sol_basis, CEED_VECTOR_ACTIVE); 25866087c08SValeria Barra 25966087c08SValeria Barra // Create auxiliary solution-size vectors. 26066087c08SValeria Barra CeedVector u, v; 26166087c08SValeria Barra CeedVectorCreate(ceed, sol_size, &u); 26266087c08SValeria Barra CeedVectorCreate(ceed, sol_size, &v); 26366087c08SValeria Barra 26466087c08SValeria Barra // Initialize 'u' with sum of coordinates, x+y+z. 265d1d35e2fSjeremylt CeedScalar *u_array; 266d1d35e2fSjeremylt const CeedScalar *x_array; 2679c774eddSJeremy L Thompson CeedVectorGetArrayWrite(u, CEED_MEM_HOST, &u_array); 268d1d35e2fSjeremylt CeedVectorGetArrayRead(mesh_coords, CEED_MEM_HOST, &x_array); 26966087c08SValeria Barra for (CeedInt i = 0; i < sol_size; i++) { 270d1d35e2fSjeremylt u_array[i] = 0; 27166087c08SValeria Barra for (CeedInt d = 0; d < dim; d++) 272d1d35e2fSjeremylt u_array[i] += x_array[i+d*sol_size]; 27366087c08SValeria Barra } 274d1d35e2fSjeremylt CeedVectorRestoreArray(u, &u_array); 275d1d35e2fSjeremylt CeedVectorRestoreArrayRead(mesh_coords, &x_array); 27666087c08SValeria Barra 277ded9b81dSJeremy L Thompson // Compute the mesh surface area using the diff operator: 278ded9b81dSJeremy L Thompson // sa = 1^T \cdot abs( K \cdot x). 279d1d35e2fSjeremylt CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE); 28066087c08SValeria Barra 28166087c08SValeria Barra // Compute and print the sum of the entries of 'v' giving the mesh surface area. 282d1d35e2fSjeremylt const CeedScalar *v_array; 283d1d35e2fSjeremylt CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array); 28466087c08SValeria Barra CeedScalar sa = 0.; 28566087c08SValeria Barra for (CeedInt i = 0; i < sol_size; i++) { 286d1d35e2fSjeremylt sa += fabs(v_array[i]); 28766087c08SValeria Barra } 288d1d35e2fSjeremylt CeedVectorRestoreArrayRead(v, &v_array); 28966087c08SValeria Barra if (!test) { 290ded9b81dSJeremy L Thompson // LCOV_EXCL_START 29166087c08SValeria Barra printf(" done.\n"); 29266087c08SValeria Barra printf("Exact mesh surface area : % .14g\n", exact_sa); 29366087c08SValeria Barra printf("Computed mesh surface area : % .14g\n", sa); 29466087c08SValeria Barra printf("Surface area error : % .14g\n", sa-exact_sa); 295ded9b81dSJeremy L Thompson // LCOV_EXCL_STOP 29666087c08SValeria Barra } else { 29780a9ef05SNatalie Beams CeedScalar tol = (dim==1 ? 10000.*CEED_EPSILON : dim==2 ? 1E-1 : 1E-1); 29866087c08SValeria Barra if (fabs(sa-exact_sa)>tol) 299e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 30066087c08SValeria Barra printf("Surface area error : % .14g\n", sa-exact_sa); 301e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 30266087c08SValeria Barra } 30366087c08SValeria Barra 30466087c08SValeria Barra // Free dynamically allocated memory. 30566087c08SValeria Barra CeedVectorDestroy(&u); 30666087c08SValeria Barra CeedVectorDestroy(&v); 307d1d35e2fSjeremylt CeedVectorDestroy(&q_data); 30866087c08SValeria Barra CeedVectorDestroy(&mesh_coords); 309d1d35e2fSjeremylt CeedOperatorDestroy(&op_apply); 310d1d35e2fSjeremylt CeedQFunctionDestroy(&qf_apply); 311777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&build_ctx); 312d1d35e2fSjeremylt CeedOperatorDestroy(&op_build); 313d1d35e2fSjeremylt CeedQFunctionDestroy(&qf_build); 31466087c08SValeria Barra CeedElemRestrictionDestroy(&sol_restr); 31566087c08SValeria Barra CeedElemRestrictionDestroy(&mesh_restr); 316d1d35e2fSjeremylt CeedElemRestrictionDestroy(&q_data_restr_i); 31766087c08SValeria Barra CeedBasisDestroy(&sol_basis); 31866087c08SValeria Barra CeedBasisDestroy(&mesh_basis); 31966087c08SValeria Barra CeedDestroy(&ceed); 32066087c08SValeria Barra return 0; 32166087c08SValeria Barra } 32266087c08SValeria Barra 323*990fdeb6SJeremy L Thompson int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, 324*990fdeb6SJeremy L Thompson CeedInt num_xyz[3]) { 32566087c08SValeria Barra // Use the approximate formula: 326ded9b81dSJeremy L Thompson // prob_size ~ num_elem * degree^dim 327ded9b81dSJeremy L Thompson CeedInt num_elem = prob_size / CeedIntPow(degree, dim); 32866087c08SValeria Barra CeedInt s = 0; // find s: num_elem/2 < 2^s <= num_elem 32966087c08SValeria Barra while (num_elem > 1) { 33066087c08SValeria Barra num_elem /= 2; 33166087c08SValeria Barra s++; 33266087c08SValeria Barra } 33366087c08SValeria Barra CeedInt r = s%dim; 334*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 335*990fdeb6SJeremy L Thompson CeedInt sd = s/dim; 33666087c08SValeria Barra if (r > 0) { sd++; r--; } 337d1d35e2fSjeremylt num_xyz[d] = 1 << sd; 33866087c08SValeria Barra } 33966087c08SValeria Barra return 0; 34066087c08SValeria Barra } 34166087c08SValeria Barra 342*990fdeb6SJeremy L Thompson int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[3], 343*990fdeb6SJeremy L Thompson CeedInt degree, CeedInt num_comp, CeedInt *size, 344*990fdeb6SJeremy L Thompson CeedInt num_qpts, CeedElemRestriction *restr, 34566087c08SValeria Barra CeedElemRestriction *restr_i) { 346ded9b81dSJeremy L Thompson CeedInt p = degree + 1; 347d1d35e2fSjeremylt CeedInt num_nodes = CeedIntPow(p, dim); // number of scalar nodes per element 34866087c08SValeria Barra CeedInt elem_qpts = CeedIntPow(num_qpts, dim); // number of qpts per element 34966087c08SValeria Barra CeedInt nd[3], num_elem = 1, scalar_size = 1; 350*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 351d1d35e2fSjeremylt num_elem *= num_xyz[d]; 352d1d35e2fSjeremylt nd[d] = num_xyz[d] * (p - 1) + 1; 35366087c08SValeria Barra scalar_size *= nd[d]; 35466087c08SValeria Barra } 355d1d35e2fSjeremylt *size = scalar_size*num_comp; 35666087c08SValeria Barra // elem: 0 1 n-1 35766087c08SValeria Barra // |---*-...-*---|---*-...-*---|- ... -|--...--| 358d1d35e2fSjeremylt // num_nodes: 0 1 p-1 p p+1 2*p n*p 359d1d35e2fSjeremylt CeedInt *el_nodes = malloc(sizeof(CeedInt)*num_elem*num_nodes); 36066087c08SValeria Barra for (CeedInt e = 0; e < num_elem; e++) { 361d1d35e2fSjeremylt CeedInt e_xyz[3] = {1, 1, 1}, re = e; 362*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { e_xyz[d] = re%num_xyz[d]; re /= num_xyz[d]; } 363d1d35e2fSjeremylt CeedInt *loc_el_nodes = el_nodes + e*num_nodes; 364*990fdeb6SJeremy L Thompson for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) { 365d1d35e2fSjeremylt CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes; 366*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 367d1d35e2fSjeremylt g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride; 368d1d35e2fSjeremylt g_nodes_stride *= nd[d]; 369d1d35e2fSjeremylt r_nodes /= p; 37066087c08SValeria Barra } 371d1d35e2fSjeremylt loc_el_nodes[l_nodes] = g_nodes; 37266087c08SValeria Barra } 37366087c08SValeria Barra } 37466087c08SValeria Barra if (restr) 375d1d35e2fSjeremylt CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, 376d1d35e2fSjeremylt num_comp * scalar_size, CEED_MEM_HOST, 377d979a051Sjeremylt CEED_COPY_VALUES, el_nodes, restr); 37866087c08SValeria Barra free(el_nodes); 3797509a596Sjeremylt 3807509a596Sjeremylt if (restr_i) { 3817509a596Sjeremylt CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, 382d1d35e2fSjeremylt num_comp, num_comp * elem_qpts * num_elem, 383523b8ea0Sjeremylt CEED_STRIDES_BACKEND, restr_i); 3847509a596Sjeremylt } 3857509a596Sjeremylt 38666087c08SValeria Barra return 0; 38766087c08SValeria Barra } 38866087c08SValeria Barra 389*990fdeb6SJeremy L Thompson int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[3], CeedInt mesh_degree, 39066087c08SValeria Barra CeedVector mesh_coords) { 391ded9b81dSJeremy L Thompson CeedInt p = mesh_degree + 1; 39266087c08SValeria Barra CeedInt nd[3], num_elem = 1, scalar_size = 1; 393*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 394d1d35e2fSjeremylt num_elem *= num_xyz[d]; 395d1d35e2fSjeremylt nd[d] = num_xyz[d] * (p - 1) + 1; 39666087c08SValeria Barra scalar_size *= nd[d]; 39766087c08SValeria Barra } 39866087c08SValeria Barra CeedScalar *coords; 3999c774eddSJeremy L Thompson CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords); 400ded9b81dSJeremy L Thompson CeedScalar *nodes = malloc(sizeof(CeedScalar) * p); 40166087c08SValeria Barra // The H1 basis uses Lobatto quadrature points as nodes. 402ded9b81dSJeremy L Thompson CeedLobattoQuadrature(p, nodes, NULL); // nodes are in [-1,1] 403ded9b81dSJeremy L Thompson for (CeedInt i = 0; i < p; i++) { nodes[i] = 0.5 + 0.5 * nodes[i]; } 404d1d35e2fSjeremylt for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) { 405d1d35e2fSjeremylt CeedInt r_nodes = gs_nodes; 406*990fdeb6SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 407d1d35e2fSjeremylt CeedInt d1d = r_nodes % nd[d]; 408d1d35e2fSjeremylt coords[gs_nodes + scalar_size * d] = ((d1d / (p - 1)) + nodes[d1d % 409d1d35e2fSjeremylt (p - 1)]) / num_xyz[d]; 410d1d35e2fSjeremylt r_nodes /= nd[d]; 41166087c08SValeria Barra } 41266087c08SValeria Barra } 41366087c08SValeria Barra free(nodes); 41466087c08SValeria Barra CeedVectorRestoreArray(mesh_coords, &coords); 41566087c08SValeria Barra return 0; 41666087c08SValeria Barra } 41766087c08SValeria Barra 41866087c08SValeria Barra #ifndef M_PI 41966087c08SValeria Barra #define M_PI 3.14159265358979323846 42066087c08SValeria Barra #endif 42166087c08SValeria Barra 422*990fdeb6SJeremy L Thompson CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, 423*990fdeb6SJeremy L Thompson CeedVector mesh_coords) { 42466087c08SValeria Barra CeedScalar exact_sa = (dim == 1 ? 2 : dim == 2 ? 4 : 6); 42566087c08SValeria Barra CeedScalar *coords; 42666087c08SValeria Barra 42766087c08SValeria Barra CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords); 42866087c08SValeria Barra for (CeedInt i = 0; i < mesh_size; i++) { 42966087c08SValeria Barra // map [0,1] to [0,1] varying the mesh density 43066087c08SValeria Barra coords[i] = 0.5 + 1./sqrt(3.) * sin((2./3.) * M_PI * (coords[i] - 0.5)); 43166087c08SValeria Barra } 43266087c08SValeria Barra CeedVectorRestoreArray(mesh_coords, &coords); 43366087c08SValeria Barra 43466087c08SValeria Barra return exact_sa; 43566087c08SValeria Barra } 436