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