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