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