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