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