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