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