xref: /libCEED/examples/ceed/ex1-volume.c (revision 59fa3f92dbea5a7e09ccfccd0c5c76c4fcbd3373)
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 
121   CeedInit(ceed_spec, &ceed);
122 
123   // Construct the mesh and solution bases.
124   CeedBasis mesh_basis, sol_basis;
125 
126   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, mesh_degree + 1, num_qpts, CEED_GAUSS, &mesh_basis);
127   CeedBasisCreateTensorH1Lagrange(ceed, dim, 1, sol_degree + 1, num_qpts, CEED_GAUSS, &sol_basis);
128 
129   // Determine the mesh size based on the given approximate problem size.
130   CeedInt num_xyz[dim];
131 
132   GetCartesianMeshSize(dim, sol_degree, prob_size, num_xyz);
133   if (!test) {
134     // LCOV_EXCL_START
135     printf("Mesh size: nx = %" CeedInt_FMT, num_xyz[0]);
136     if (dim > 1) printf(", ny = %" CeedInt_FMT, num_xyz[1]);
137     if (dim > 2) printf(", nz = %" CeedInt_FMT, num_xyz[2]);
138     printf("\n");
139     // LCOV_EXCL_STOP
140   }
141 
142   // Build CeedElemRestriction objects describing the mesh and solution discrete representations.
143   CeedInt             mesh_size, sol_size;
144   CeedElemRestriction mesh_restriction, sol_restriction, q_data_restriction;
145 
146   BuildCartesianRestriction(ceed, dim, num_xyz, mesh_degree, num_comp_x, &mesh_size, num_qpts, &mesh_restriction, NULL);
147   BuildCartesianRestriction(ceed, dim, num_xyz, sol_degree, 1, &sol_size, num_qpts, &sol_restriction, &q_data_restriction);
148   if (!test) {
149     // LCOV_EXCL_START
150     printf("Number of mesh nodes     : %" CeedInt_FMT "\n", mesh_size / dim);
151     printf("Number of solution nodes : %" CeedInt_FMT "\n", sol_size);
152     // LCOV_EXCL_STOP
153   }
154 
155   // Create a CeedVector with the mesh coordinates.
156   CeedVector mesh_coords;
157 
158   CeedVectorCreate(ceed, mesh_size, &mesh_coords);
159   SetCartesianMeshCoords(dim, num_xyz, mesh_degree, mesh_coords);
160 
161   // Apply a transformation to the mesh.
162   CeedScalar exact_volume = TransformMeshCoords(dim, mesh_size, mesh_coords);
163 
164   // Context data to be passed to the 'build_mass' QFunction.
165   CeedQFunctionContext build_ctx;
166   struct BuildContext  build_ctx_data;
167 
168   build_ctx_data.dim = build_ctx_data.space_dim = dim;
169   CeedQFunctionContextCreate(ceed, &build_ctx);
170   CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
171 
172   // Create the QFunction that builds the mass operator (i.e. computes its quadrature data) and set its context data.
173   CeedQFunction qf_build;
174 
175   if (gallery) {
176     // This creates the QFunction via the gallery.
177     char name[13] = "";
178     snprintf(name, sizeof name, "Mass%" CeedInt_FMT "DBuild", dim);
179     CeedQFunctionCreateInteriorByName(ceed, name, &qf_build);
180   } else {
181     // This creates the QFunction directly.
182     CeedQFunctionCreateInterior(ceed, 1, build_mass, build_mass_loc, &qf_build);
183     CeedQFunctionAddInput(qf_build, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
184     CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
185     CeedQFunctionAddOutput(qf_build, "qdata", 1, CEED_EVAL_NONE);
186     CeedQFunctionSetContext(qf_build, build_ctx);
187   }
188 
189   // Create the operator that builds the quadrature data for the mass operator.
190   CeedOperator op_build;
191 
192   CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
193   CeedOperatorSetField(op_build, "dx", mesh_restriction, mesh_basis, CEED_VECTOR_ACTIVE);
194   CeedOperatorSetField(op_build, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE);
195   CeedOperatorSetField(op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
196 
197   // Compute the quadrature data for the mass operator.
198   CeedVector q_data;
199   CeedInt    elem_qpts = CeedIntPow(num_qpts, dim);
200   CeedInt    num_elem  = 1;
201 
202   for (CeedInt d = 0; d < dim; d++) num_elem *= num_xyz[d];
203   CeedVectorCreate(ceed, num_elem * elem_qpts, &q_data);
204   CeedOperatorApply(op_build, mesh_coords, q_data, CEED_REQUEST_IMMEDIATE);
205 
206   // Create the QFunction that defines the action of the mass operator.
207   CeedQFunction qf_apply;
208 
209   if (gallery) {
210     // This creates the QFunction via the gallery.
211     CeedQFunctionCreateInteriorByName(ceed, "MassApply", &qf_apply);
212   } else {
213     // This creates the QFunction directly.
214     CeedQFunctionCreateInterior(ceed, 1, apply_mass, apply_mass_loc, &qf_apply);
215     CeedQFunctionAddInput(qf_apply, "u", 1, CEED_EVAL_INTERP);
216     CeedQFunctionAddInput(qf_apply, "qdata", 1, CEED_EVAL_NONE);
217     CeedQFunctionAddOutput(qf_apply, "v", 1, CEED_EVAL_INTERP);
218   }
219 
220   // Create the mass operator.
221   CeedOperator op_apply;
222 
223   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
224   CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
225   CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data);
226   CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
227 
228   // Create auxiliary solution-size vectors.
229   CeedVector u, v;
230 
231   CeedVectorCreate(ceed, sol_size, &u);
232   CeedVectorCreate(ceed, sol_size, &v);
233 
234   // Initialize 'u' with ones.
235   CeedVectorSetValue(u, 1.0);
236 
237   // Compute the mesh volume using the mass operator: volume = 1^T \cdot M \cdot 1
238   CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
239 
240   // Benchmark runs
241   if (!test && benchmark) {
242     // LCOV_EXCL_START
243     printf(" Executing %d benchmarking runs...\n", benchmark);
244     // LCOV_EXCL_STOP
245   }
246   for (CeedInt i = 0; i < benchmark; i++) {
247     // LCOV_EXCL_START
248     CeedOperatorApply(op_apply, u, v, CEED_REQUEST_IMMEDIATE);
249     // LCOV_EXCL_STOP
250   }
251 
252   // Compute and print the sum of the entries of 'v' giving the mesh volume.
253   CeedScalar volume = 0.;
254 
255   {
256     const CeedScalar *v_array;
257 
258     CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
259     for (CeedInt i = 0; i < sol_size; i++) volume += v_array[i];
260     CeedVectorRestoreArrayRead(v, &v_array);
261   }
262   if (!test) {
263     // LCOV_EXCL_START
264     printf(" done.\n");
265     printf("Exact mesh volume    : % .14g\n", exact_volume);
266     printf("Computed mesh volume : % .14g\n", volume);
267     printf("Volume error         : % .14g\n", volume - exact_volume);
268     // LCOV_EXCL_STOP
269   } else {
270     CeedScalar tol = (dim == 1 ? 200. * CEED_EPSILON : dim == 2 ? 1E-5 : 1E-5);
271 
272     if (fabs(volume - exact_volume) > tol) printf("Volume error : % .1e\n", volume - exact_volume);
273   }
274 
275   // Free dynamically allocated memory.
276   CeedVectorDestroy(&u);
277   CeedVectorDestroy(&v);
278   CeedVectorDestroy(&q_data);
279   CeedVectorDestroy(&mesh_coords);
280   CeedOperatorDestroy(&op_apply);
281   CeedQFunctionDestroy(&qf_apply);
282   CeedQFunctionContextDestroy(&build_ctx);
283   CeedOperatorDestroy(&op_build);
284   CeedQFunctionDestroy(&qf_build);
285   CeedElemRestrictionDestroy(&sol_restriction);
286   CeedElemRestrictionDestroy(&mesh_restriction);
287   CeedElemRestrictionDestroy(&q_data_restriction);
288   CeedBasisDestroy(&sol_basis);
289   CeedBasisDestroy(&mesh_basis);
290   CeedDestroy(&ceed);
291   return 0;
292 }
293 
294 int GetCartesianMeshSize(CeedInt dim, CeedInt degree, CeedInt prob_size, CeedInt num_xyz[dim]) {
295   // Use the approximate formula:
296   //    prob_size ~ num_elem * degree^dim
297   CeedInt num_elem = prob_size / CeedIntPow(degree, dim);
298   CeedInt s        = 0;  // find s: num_elem/2 < 2^s <= num_elem
299 
300   while (num_elem > 1) {
301     num_elem /= 2;
302     s++;
303   }
304   CeedInt r = s % dim;
305 
306   for (CeedInt d = 0; d < dim; d++) {
307     CeedInt sd = s / dim;
308 
309     if (r > 0) {
310       sd++;
311       r--;
312     }
313     num_xyz[d] = 1 << sd;
314   }
315   return 0;
316 }
317 
318 int BuildCartesianRestriction(Ceed ceed, CeedInt dim, CeedInt num_xyz[dim], CeedInt degree, CeedInt num_comp, CeedInt *size, CeedInt num_qpts,
319                               CeedElemRestriction *restriction, CeedElemRestriction *q_data_restriction) {
320   CeedInt p         = degree + 1;
321   CeedInt num_nodes = CeedIntPow(p, dim);         // number of scalar nodes per element
322   CeedInt elem_qpts = CeedIntPow(num_qpts, dim);  // number of qpts per element
323   CeedInt nd[3], num_elem = 1, scalar_size = 1;
324 
325   for (CeedInt d = 0; d < dim; d++) {
326     num_elem *= num_xyz[d];
327     nd[d] = num_xyz[d] * (p - 1) + 1;
328     scalar_size *= nd[d];
329   }
330   *size = scalar_size * num_comp;
331   // elem:         0             1                 n-1
332   //           |---*-...-*---|---*-...-*---|- ... -|--...--|
333   // num_nodes:   0   1    p-1  p  p+1       2*p             n*p
334   CeedInt *elem_nodes = malloc(sizeof(CeedInt) * num_elem * num_nodes);
335 
336   for (CeedInt e = 0; e < num_elem; e++) {
337     CeedInt e_xyz[3] = {1, 1, 1}, re = e;
338 
339     for (CeedInt d = 0; d < dim; d++) {
340       e_xyz[d] = re % num_xyz[d];
341       re /= num_xyz[d];
342     }
343     CeedInt *local_elem_nodes = elem_nodes + e * num_nodes;
344 
345     for (CeedInt l_nodes = 0; l_nodes < num_nodes; l_nodes++) {
346       CeedInt g_nodes = 0, g_nodes_stride = 1, r_nodes = l_nodes;
347 
348       for (CeedInt d = 0; d < dim; d++) {
349         g_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * g_nodes_stride;
350         g_nodes_stride *= nd[d];
351         r_nodes /= p;
352       }
353       local_elem_nodes[l_nodes] = g_nodes;
354     }
355   }
356   CeedElemRestrictionCreate(ceed, num_elem, num_nodes, num_comp, scalar_size, num_comp * scalar_size, CEED_MEM_HOST, CEED_COPY_VALUES, elem_nodes,
357                             restriction);
358   if (q_data_restriction) {
359     CeedElemRestrictionCreateStrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem, CEED_STRIDES_BACKEND, q_data_restriction);
360   }
361   free(elem_nodes);
362   return 0;
363 }
364 
365 int SetCartesianMeshCoords(CeedInt dim, CeedInt num_xyz[dim], CeedInt mesh_degree, CeedVector mesh_coords) {
366   CeedInt p = mesh_degree + 1;
367   CeedInt nd[3], scalar_size = 1;
368 
369   for (CeedInt d = 0; d < dim; d++) {
370     nd[d] = num_xyz[d] * (p - 1) + 1;
371     scalar_size *= nd[d];
372   }
373   CeedScalar *coords;
374 
375   CeedVectorGetArrayWrite(mesh_coords, CEED_MEM_HOST, &coords);
376   CeedScalar *nodes = malloc(sizeof(CeedScalar) * p);
377 
378   // The H1 basis uses Lobatto quadrature points as nodes.
379   CeedLobattoQuadrature(p, nodes, NULL);  // nodes are in [-1,1]
380   for (CeedInt i = 0; i < p; i++) nodes[i] = 0.5 + 0.5 * nodes[i];
381   for (CeedInt gs_nodes = 0; gs_nodes < scalar_size; gs_nodes++) {
382     CeedInt r_nodes = gs_nodes;
383 
384     for (CeedInt d = 0; d < dim; d++) {
385       CeedInt d_1d = r_nodes % nd[d];
386 
387       coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) + nodes[d_1d % (p - 1)]) / num_xyz[d];
388       r_nodes /= nd[d];
389     }
390   }
391   free(nodes);
392   CeedVectorRestoreArray(mesh_coords, &coords);
393   return 0;
394 }
395 
396 #ifndef M_PI
397 #define M_PI 3.14159265358979323846
398 #define M_PI_2 1.57079632679489661923
399 #endif
400 
401 CeedScalar TransformMeshCoords(CeedInt dim, CeedInt mesh_size, CeedVector mesh_coords) {
402   CeedScalar  exact_volume;
403   CeedScalar *coords;
404 
405   CeedVectorGetArray(mesh_coords, CEED_MEM_HOST, &coords);
406   if (dim == 1) {
407     for (CeedInt i = 0; i < mesh_size; i++) {
408       // map [0,1] to [0,1] varying the mesh density
409       coords[i] = 0.5 + 1. / sqrt(3.) * sin((2. / 3.) * M_PI * (coords[i] - 0.5));
410     }
411     exact_volume = 1.;
412   } else {
413     CeedInt num_nodes = mesh_size / dim;
414 
415     for (CeedInt i = 0; i < num_nodes; i++) {
416       // map (x,y) from [0,1]x[0,1] to the quarter annulus with polar
417       // coordinates, (r,phi) in [1,2]x[0,pi/2] with area = 3/4*pi
418       CeedScalar u = coords[i], v = coords[i + num_nodes];
419 
420       u                     = 1. + u;
421       v                     = M_PI_2 * v;
422       coords[i]             = u * cos(v);
423       coords[i + num_nodes] = u * sin(v);
424     }
425     exact_volume = 3. / 4. * M_PI;
426   }
427   CeedVectorRestoreArray(mesh_coords, &coords);
428   return exact_volume;
429 }
430