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