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