xref: /petsc/src/dm/impls/plex/transform/impls/extrude/plextrextrude.c (revision f78ce678aed4d3e450e422d3a2725e4faeb03dcb)
1 #include <petsc/private/dmplextransformimpl.h> /*I "petscdmplextransform.h" I*/
2 
3 static PetscErrorCode DMPlexTransformView_Extrude(DMPlexTransform tr, PetscViewer viewer)
4 {
5   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
6   PetscBool                isascii;
7 
8   PetscFunctionBegin;
9   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
10   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
11   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii));
12   if (isascii) {
13     const char *name;
14 
15     PetscCall(PetscObjectGetName((PetscObject) tr, &name));
16     PetscCall(PetscViewerASCIIPrintf(viewer, "Extrusion transformation %s\n", name ? name : ""));
17     PetscCall(PetscViewerASCIIPrintf(viewer, "  number of layers: %" PetscInt_FMT "\n", ex->layers));
18     PetscCall(PetscViewerASCIIPrintf(viewer, "  create tensor cells: %s\n", ex->useTensor ? "YES" : "NO"));
19   } else {
20     SETERRQ(PetscObjectComm((PetscObject) tr), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlexTransform writing", ((PetscObject) viewer)->type_name);
21   }
22   PetscFunctionReturn(0);
23 }
24 
25 static PetscErrorCode DMPlexTransformSetFromOptions_Extrude(PetscOptionItems *PetscOptionsObject, DMPlexTransform tr)
26 {
27   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
28   PetscReal                th, normal[3], *thicknesses;
29   PetscInt                 nl, Nc;
30   PetscBool                tensor, sym, flg;
31   char                     funcname[PETSC_MAX_PATH_LEN];
32 
33   PetscFunctionBegin;
34   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 2);
35   PetscOptionsHeadBegin(PetscOptionsObject, "DMPlexTransform Extrusion Options");
36   PetscCall(PetscOptionsBoundedInt("-dm_plex_transform_extrude_layers", "Number of layers to extrude", "", ex->layers, &nl, &flg, 1));
37   if (flg) PetscCall(DMPlexTransformExtrudeSetLayers(tr, nl));
38   PetscCall(PetscOptionsReal("-dm_plex_transform_extrude_thickness", "Total thickness of extruded layers", "", ex->thickness, &th, &flg));
39   if (flg) PetscCall(DMPlexTransformExtrudeSetThickness(tr, th));
40   PetscCall(PetscOptionsBool("-dm_plex_transform_extrude_use_tensor", "Create tensor cells", "", ex->useTensor, &tensor, &flg));
41   if (flg) PetscCall(DMPlexTransformExtrudeSetTensor(tr, tensor));
42   PetscCall(PetscOptionsBool("-dm_plex_transform_extrude_symmetric", "Extrude layers symmetrically about the surface", "", ex->symmetric, &sym, &flg));
43   if (flg) PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, sym));
44   Nc = 3;
45   PetscCall(PetscOptionsRealArray("-dm_plex_transform_extrude_normal", "Input normal vector for extrusion", "DMPlexTransformExtrudeSetNormal", normal, &Nc, &flg));
46   if (flg) {
47     // Extrusion dimension might not yet be determined
48     PetscCheck(!ex->cdimEx || Nc == ex->cdimEx,PetscObjectComm((PetscObject) tr), PETSC_ERR_ARG_SIZ, "Input normal has size %" PetscInt_FMT " != %" PetscInt_FMT " extruded coordinate dimension", Nc, ex->cdimEx);
49     PetscCall(DMPlexTransformExtrudeSetNormal(tr, normal));
50   }
51   PetscCall(PetscOptionsString("-dm_plex_transform_extrude_normal_function", "Function to determine normal vector", "DMPlexTransformExtrudeSetNormalFunction", funcname, funcname, sizeof(funcname), &flg));
52   if (flg) {
53     PetscSimplePointFunc normalFunc;
54 
55     PetscCall(PetscDLSym(NULL, funcname, (void **) &normalFunc));
56     PetscCall(DMPlexTransformExtrudeSetNormalFunction(tr, normalFunc));
57   }
58   nl   = ex->layers;
59   PetscCall(PetscMalloc1(nl, &thicknesses));
60   PetscCall(PetscOptionsRealArray("-dm_plex_transform_extrude_thicknesses", "Thickness of each individual extruded layer", "", thicknesses, &nl, &flg));
61   if (flg) {
62     PetscCheck(nl,PetscObjectComm((PetscObject) tr), PETSC_ERR_ARG_OUTOFRANGE, "Must give at least one thickness for -dm_plex_transform_extrude_thicknesses");
63     PetscCall(DMPlexTransformExtrudeSetThicknesses(tr, nl, thicknesses));
64   }
65   PetscCall(PetscFree(thicknesses));
66   PetscOptionsHeadEnd();
67   PetscFunctionReturn(0);
68 }
69 
70 /* Determine the implicit dimension pre-extrusion (either the implicit dimension of the DM or of a point in the active set for the transform).
71    If that dimension is the same as the current coordinate dimension (ex->dim), the extruded mesh will have a coordinate dimension one greater;
72    Otherwise the coordinate dimension will be kept. */
73 static PetscErrorCode DMPlexTransformExtrudeComputeExtrusionDim(DMPlexTransform tr)
74 {
75   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
76   DM                       dm;
77   DMLabel                  active;
78   PetscInt                 dim, dimOrig;
79 
80   PetscFunctionBegin;
81   PetscCall(DMPlexTransformGetDM(tr, &dm));
82   PetscCall(DMGetDimension(dm, &dim));
83   dimOrig = dim;
84   PetscCall(DMPlexTransformGetActive(tr, &active));
85   if (active) {
86     PetscInt pStart, pEnd, p;
87 
88     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
89     for (p = pStart; p < pEnd; ++p) {
90       DMPolytopeType ct;
91       PetscInt       val;
92 
93       PetscCall(DMLabelGetValue(active, p, &val));
94       if (val < 0) continue;
95       PetscCall(DMPlexGetCellType(dm, p, &ct));
96       dim = DMPolytopeTypeGetDim(ct);
97       break;
98     }
99   }
100   ex->dimEx  = PetscMax(dimOrig, dim+1);
101   ex->cdimEx = ex->cdim == dim ? ex->cdim+1 : ex->cdim;
102   PetscFunctionReturn(0);
103 }
104 
105 static PetscErrorCode DMPlexTransformSetDimensions_Extrude(DMPlexTransform tr, DM dm, DM tdm)
106 {
107   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
108   PetscInt                 dim;
109 
110   PetscFunctionBegin;
111   PetscCall(DMGetDimension(dm, &dim));
112   PetscCall(DMSetDimension(tdm, ex->dimEx));
113   PetscCall(DMSetCoordinateDim(tdm, ex->cdimEx));
114   PetscFunctionReturn(0);
115 }
116 
117 /*
118   The refine types for extrusion are:
119 
120   ct:       For any normally extruded point
121   ct + 100: For any point which should just return itself
122 */
123 static PetscErrorCode DMPlexTransformSetUp_Extrude(DMPlexTransform tr)
124 {
125   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
126   DM                       dm;
127   DMLabel                  active;
128   DMPolytopeType           ct;
129   PetscInt                 Nl = ex->layers, l, i, ict, Nc, No, coff, ooff;
130 
131   PetscFunctionBegin;
132   PetscCall(DMPlexTransformExtrudeComputeExtrusionDim(tr));
133   PetscCall(DMPlexTransformGetDM(tr, &dm));
134   PetscCall(DMPlexTransformGetActive(tr, &active));
135   if (active) {
136     DMLabel  celltype;
137     PetscInt pStart, pEnd, p;
138 
139     PetscCall(DMPlexGetCellTypeLabel(dm, &celltype));
140     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Refine Type", &tr->trType));
141     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
142     for (p = pStart; p < pEnd; ++p) {
143       PetscInt ct, val;
144 
145       PetscCall(DMLabelGetValue(celltype, p, &ct));
146       PetscCall(DMLabelGetValue(active, p, &val));
147       if (val < 0) {PetscCall(DMLabelSetValue(tr->trType, p, ct + 100));}
148       else         {PetscCall(DMLabelSetValue(tr->trType, p, ct));}
149     }
150   }
151   PetscCall(PetscMalloc5(DM_NUM_POLYTOPES, &ex->Nt, DM_NUM_POLYTOPES, &ex->target, DM_NUM_POLYTOPES, &ex->size, DM_NUM_POLYTOPES, &ex->cone, DM_NUM_POLYTOPES, &ex->ornt));
152   for (ict = 0; ict < DM_NUM_POLYTOPES; ++ict) {
153     ex->Nt[ict]     = -1;
154     ex->target[ict] = NULL;
155     ex->size[ict]   = NULL;
156     ex->cone[ict]   = NULL;
157     ex->ornt[ict]   = NULL;
158   }
159   /* DM_POLYTOPE_POINT produces Nl+1 points and Nl segments/tensor segments */
160   ct = DM_POLYTOPE_POINT;
161   ex->Nt[ct] = 2;
162   Nc = 6*Nl;
163   No = 2*Nl;
164   PetscCall(PetscMalloc4(ex->Nt[ct], &ex->target[ct], ex->Nt[ct], &ex->size[ct], Nc, &ex->cone[ct], No, &ex->ornt[ct]));
165   ex->target[ct][0] = DM_POLYTOPE_POINT;
166   ex->target[ct][1] = ex->useTensor ? DM_POLYTOPE_POINT_PRISM_TENSOR : DM_POLYTOPE_SEGMENT;
167   ex->size[ct][0]   = Nl+1;
168   ex->size[ct][1]   = Nl;
169   /*   cones for segments/tensor segments */
170   for (i = 0; i < Nl; ++i) {
171     ex->cone[ct][6*i+0] = DM_POLYTOPE_POINT;
172     ex->cone[ct][6*i+1] = 0;
173     ex->cone[ct][6*i+2] = i;
174     ex->cone[ct][6*i+3] = DM_POLYTOPE_POINT;
175     ex->cone[ct][6*i+4] = 0;
176     ex->cone[ct][6*i+5] = i+1;
177   }
178   for (i = 0; i < No; ++i) ex->ornt[ct][i] = 0;
179   /* DM_POLYTOPE_SEGMENT produces Nl+1 segments and Nl quads/tensor quads */
180   ct = DM_POLYTOPE_SEGMENT;
181   ex->Nt[ct] = 2;
182   Nc = 8*(Nl+1) + 14*Nl;
183   No = 2*(Nl+1) + 4*Nl;
184   PetscCall(PetscMalloc4(ex->Nt[ct], &ex->target[ct], ex->Nt[ct], &ex->size[ct], Nc, &ex->cone[ct], No, &ex->ornt[ct]));
185   ex->target[ct][0] = DM_POLYTOPE_SEGMENT;
186   ex->target[ct][1] = ex->useTensor ? DM_POLYTOPE_SEG_PRISM_TENSOR : DM_POLYTOPE_QUADRILATERAL;
187   ex->size[ct][0]   = Nl+1;
188   ex->size[ct][1]   = Nl;
189   /*   cones for segments */
190   for (i = 0; i < Nl+1; ++i) {
191     ex->cone[ct][8*i+0] = DM_POLYTOPE_POINT;
192     ex->cone[ct][8*i+1] = 1;
193     ex->cone[ct][8*i+2] = 0;
194     ex->cone[ct][8*i+3] = i;
195     ex->cone[ct][8*i+4] = DM_POLYTOPE_POINT;
196     ex->cone[ct][8*i+5] = 1;
197     ex->cone[ct][8*i+6] = 1;
198     ex->cone[ct][8*i+7] = i;
199   }
200   for (i = 0; i < 2*(Nl+1); ++i) ex->ornt[ct][i] = 0;
201   /*   cones for quads/tensor quads */
202   coff = 8*(Nl+1);
203   ooff = 2*(Nl+1);
204   for (i = 0; i < Nl; ++i) {
205     if (ex->useTensor) {
206       ex->cone[ct][coff+14*i+0]  = DM_POLYTOPE_SEGMENT;
207       ex->cone[ct][coff+14*i+1]  = 0;
208       ex->cone[ct][coff+14*i+2]  = i;
209       ex->cone[ct][coff+14*i+3]  = DM_POLYTOPE_SEGMENT;
210       ex->cone[ct][coff+14*i+4]  = 0;
211       ex->cone[ct][coff+14*i+5]  = i+1;
212       ex->cone[ct][coff+14*i+6]  = DM_POLYTOPE_POINT_PRISM_TENSOR;
213       ex->cone[ct][coff+14*i+7]  = 1;
214       ex->cone[ct][coff+14*i+8]  = 0;
215       ex->cone[ct][coff+14*i+9]  = i;
216       ex->cone[ct][coff+14*i+10] = DM_POLYTOPE_POINT_PRISM_TENSOR;
217       ex->cone[ct][coff+14*i+11] = 1;
218       ex->cone[ct][coff+14*i+12] = 1;
219       ex->cone[ct][coff+14*i+13] = i;
220       ex->ornt[ct][ooff+4*i+0] = 0;
221       ex->ornt[ct][ooff+4*i+1] = 0;
222       ex->ornt[ct][ooff+4*i+2] = 0;
223       ex->ornt[ct][ooff+4*i+3] = 0;
224     } else {
225       ex->cone[ct][coff+14*i+0]  = DM_POLYTOPE_SEGMENT;
226       ex->cone[ct][coff+14*i+1]  = 0;
227       ex->cone[ct][coff+14*i+2]  = i;
228       ex->cone[ct][coff+14*i+3]  = DM_POLYTOPE_SEGMENT;
229       ex->cone[ct][coff+14*i+4]  = 1;
230       ex->cone[ct][coff+14*i+5]  = 1;
231       ex->cone[ct][coff+14*i+6]  = i;
232       ex->cone[ct][coff+14*i+7]  = DM_POLYTOPE_SEGMENT;
233       ex->cone[ct][coff+14*i+8]  = 0;
234       ex->cone[ct][coff+14*i+9]  = i+1;
235       ex->cone[ct][coff+14*i+10] = DM_POLYTOPE_SEGMENT;
236       ex->cone[ct][coff+14*i+11] = 1;
237       ex->cone[ct][coff+14*i+12] = 0;
238       ex->cone[ct][coff+14*i+13] = i;
239       ex->ornt[ct][ooff+4*i+0] =  0;
240       ex->ornt[ct][ooff+4*i+1] =  0;
241       ex->ornt[ct][ooff+4*i+2] = -1;
242       ex->ornt[ct][ooff+4*i+3] = -1;
243     }
244   }
245   /* DM_POLYTOPE_TRIANGLE produces Nl+1 triangles and Nl triangular prisms/tensor triangular prisms */
246   ct = DM_POLYTOPE_TRIANGLE;
247   ex->Nt[ct] = 2;
248   Nc = 12*(Nl+1) + 18*Nl;
249   No =  3*(Nl+1) +  5*Nl;
250   PetscCall(PetscMalloc4(ex->Nt[ct], &ex->target[ct], ex->Nt[ct], &ex->size[ct], Nc, &ex->cone[ct], No, &ex->ornt[ct]));
251   ex->target[ct][0] = DM_POLYTOPE_TRIANGLE;
252   ex->target[ct][1] = ex->useTensor ? DM_POLYTOPE_TRI_PRISM_TENSOR : DM_POLYTOPE_TRI_PRISM;
253   ex->size[ct][0]   = Nl+1;
254   ex->size[ct][1]   = Nl;
255   /*   cones for triangles */
256   for (i = 0; i < Nl+1; ++i) {
257     ex->cone[ct][12*i+0]  = DM_POLYTOPE_SEGMENT;
258     ex->cone[ct][12*i+1]  = 1;
259     ex->cone[ct][12*i+2]  = 0;
260     ex->cone[ct][12*i+3]  = i;
261     ex->cone[ct][12*i+4]  = DM_POLYTOPE_SEGMENT;
262     ex->cone[ct][12*i+5]  = 1;
263     ex->cone[ct][12*i+6]  = 1;
264     ex->cone[ct][12*i+7]  = i;
265     ex->cone[ct][12*i+8]  = DM_POLYTOPE_SEGMENT;
266     ex->cone[ct][12*i+9]  = 1;
267     ex->cone[ct][12*i+10] = 2;
268     ex->cone[ct][12*i+11] = i;
269   }
270   for (i = 0; i < 3*(Nl+1); ++i) ex->ornt[ct][i] = 0;
271   /*   cones for triangular prisms/tensor triangular prisms */
272   coff = 12*(Nl+1);
273   ooff = 3*(Nl+1);
274   for (i = 0; i < Nl; ++i) {
275     if (ex->useTensor) {
276       ex->cone[ct][coff+18*i+0]  = DM_POLYTOPE_TRIANGLE;
277       ex->cone[ct][coff+18*i+1]  = 0;
278       ex->cone[ct][coff+18*i+2]  = i;
279       ex->cone[ct][coff+18*i+3]  = DM_POLYTOPE_TRIANGLE;
280       ex->cone[ct][coff+18*i+4]  = 0;
281       ex->cone[ct][coff+18*i+5]  = i+1;
282       ex->cone[ct][coff+18*i+6]  = DM_POLYTOPE_SEG_PRISM_TENSOR;
283       ex->cone[ct][coff+18*i+7]  = 1;
284       ex->cone[ct][coff+18*i+8]  = 0;
285       ex->cone[ct][coff+18*i+9]  = i;
286       ex->cone[ct][coff+18*i+10] = DM_POLYTOPE_SEG_PRISM_TENSOR;
287       ex->cone[ct][coff+18*i+11] = 1;
288       ex->cone[ct][coff+18*i+12] = 1;
289       ex->cone[ct][coff+18*i+13] = i;
290       ex->cone[ct][coff+18*i+14] = DM_POLYTOPE_SEG_PRISM_TENSOR;
291       ex->cone[ct][coff+18*i+15] = 1;
292       ex->cone[ct][coff+18*i+16] = 2;
293       ex->cone[ct][coff+18*i+17] = i;
294       ex->ornt[ct][ooff+5*i+0] = 0;
295       ex->ornt[ct][ooff+5*i+1] = 0;
296       ex->ornt[ct][ooff+5*i+2] = 0;
297       ex->ornt[ct][ooff+5*i+3] = 0;
298       ex->ornt[ct][ooff+5*i+4] = 0;
299     } else {
300       ex->cone[ct][coff+18*i+0]  = DM_POLYTOPE_TRIANGLE;
301       ex->cone[ct][coff+18*i+1]  = 0;
302       ex->cone[ct][coff+18*i+2]  = i;
303       ex->cone[ct][coff+18*i+3]  = DM_POLYTOPE_TRIANGLE;
304       ex->cone[ct][coff+18*i+4]  = 0;
305       ex->cone[ct][coff+18*i+5]  = i+1;
306       ex->cone[ct][coff+18*i+6]  = DM_POLYTOPE_QUADRILATERAL;
307       ex->cone[ct][coff+18*i+7]  = 1;
308       ex->cone[ct][coff+18*i+8]  = 0;
309       ex->cone[ct][coff+18*i+9]  = i;
310       ex->cone[ct][coff+18*i+10] = DM_POLYTOPE_QUADRILATERAL;
311       ex->cone[ct][coff+18*i+11] = 1;
312       ex->cone[ct][coff+18*i+12] = 1;
313       ex->cone[ct][coff+18*i+13] = i;
314       ex->cone[ct][coff+18*i+14] = DM_POLYTOPE_QUADRILATERAL;
315       ex->cone[ct][coff+18*i+15] = 1;
316       ex->cone[ct][coff+18*i+16] = 2;
317       ex->cone[ct][coff+18*i+17] = i;
318       ex->ornt[ct][ooff+5*i+0] = -2;
319       ex->ornt[ct][ooff+5*i+1] =  0;
320       ex->ornt[ct][ooff+5*i+2] =  0;
321       ex->ornt[ct][ooff+5*i+3] =  0;
322       ex->ornt[ct][ooff+5*i+4] =  0;
323     }
324   }
325   /* DM_POLYTOPE_QUADRILATERAL produces Nl+1 quads and Nl hexes/tensor hexes */
326   ct = DM_POLYTOPE_QUADRILATERAL;
327   ex->Nt[ct] = 2;
328   Nc = 16*(Nl+1) + 22*Nl;
329   No =  4*(Nl+1) +  6*Nl;
330   PetscCall(PetscMalloc4(ex->Nt[ct], &ex->target[ct], ex->Nt[ct], &ex->size[ct], Nc, &ex->cone[ct], No, &ex->ornt[ct]));
331   ex->target[ct][0] = DM_POLYTOPE_QUADRILATERAL;
332   ex->target[ct][1] = ex->useTensor ? DM_POLYTOPE_QUAD_PRISM_TENSOR : DM_POLYTOPE_HEXAHEDRON;
333   ex->size[ct][0]   = Nl+1;
334   ex->size[ct][1]   = Nl;
335   /*   cones for quads */
336   for (i = 0; i < Nl+1; ++i) {
337     ex->cone[ct][16*i+0]  = DM_POLYTOPE_SEGMENT;
338     ex->cone[ct][16*i+1]  = 1;
339     ex->cone[ct][16*i+2]  = 0;
340     ex->cone[ct][16*i+3]  = i;
341     ex->cone[ct][16*i+4]  = DM_POLYTOPE_SEGMENT;
342     ex->cone[ct][16*i+5]  = 1;
343     ex->cone[ct][16*i+6]  = 1;
344     ex->cone[ct][16*i+7]  = i;
345     ex->cone[ct][16*i+8]  = DM_POLYTOPE_SEGMENT;
346     ex->cone[ct][16*i+9]  = 1;
347     ex->cone[ct][16*i+10] = 2;
348     ex->cone[ct][16*i+11] = i;
349     ex->cone[ct][16*i+12] = DM_POLYTOPE_SEGMENT;
350     ex->cone[ct][16*i+13] = 1;
351     ex->cone[ct][16*i+14] = 3;
352     ex->cone[ct][16*i+15] = i;
353   }
354   for (i = 0; i < 4*(Nl+1); ++i) ex->ornt[ct][i] = 0;
355   /*   cones for hexes/tensor hexes */
356   coff = 16*(Nl+1);
357   ooff = 4*(Nl+1);
358   for (i = 0; i < Nl; ++i) {
359     if (ex->useTensor) {
360       ex->cone[ct][coff+22*i+0]  = DM_POLYTOPE_QUADRILATERAL;
361       ex->cone[ct][coff+22*i+1]  = 0;
362       ex->cone[ct][coff+22*i+2]  = i;
363       ex->cone[ct][coff+22*i+3]  = DM_POLYTOPE_QUADRILATERAL;
364       ex->cone[ct][coff+22*i+4]  = 0;
365       ex->cone[ct][coff+22*i+5]  = i+1;
366       ex->cone[ct][coff+22*i+6]  = DM_POLYTOPE_SEG_PRISM_TENSOR;
367       ex->cone[ct][coff+22*i+7]  = 1;
368       ex->cone[ct][coff+22*i+8]  = 0;
369       ex->cone[ct][coff+22*i+9]  = i;
370       ex->cone[ct][coff+22*i+10] = DM_POLYTOPE_SEG_PRISM_TENSOR;
371       ex->cone[ct][coff+22*i+11] = 1;
372       ex->cone[ct][coff+22*i+12] = 1;
373       ex->cone[ct][coff+22*i+13] = i;
374       ex->cone[ct][coff+22*i+14] = DM_POLYTOPE_SEG_PRISM_TENSOR;
375       ex->cone[ct][coff+22*i+15] = 1;
376       ex->cone[ct][coff+22*i+16] = 2;
377       ex->cone[ct][coff+22*i+17] = i;
378       ex->cone[ct][coff+22*i+18] = DM_POLYTOPE_SEG_PRISM_TENSOR;
379       ex->cone[ct][coff+22*i+19] = 1;
380       ex->cone[ct][coff+22*i+20] = 3;
381       ex->cone[ct][coff+22*i+21] = i;
382       ex->ornt[ct][ooff+6*i+0] = 0;
383       ex->ornt[ct][ooff+6*i+1] = 0;
384       ex->ornt[ct][ooff+6*i+2] = 0;
385       ex->ornt[ct][ooff+6*i+3] = 0;
386       ex->ornt[ct][ooff+6*i+4] = 0;
387       ex->ornt[ct][ooff+6*i+5] = 0;
388     } else {
389       ex->cone[ct][coff+22*i+0]  = DM_POLYTOPE_QUADRILATERAL;
390       ex->cone[ct][coff+22*i+1]  = 0;
391       ex->cone[ct][coff+22*i+2]  = i;
392       ex->cone[ct][coff+22*i+3]  = DM_POLYTOPE_QUADRILATERAL;
393       ex->cone[ct][coff+22*i+4]  = 0;
394       ex->cone[ct][coff+22*i+5]  = i+1;
395       ex->cone[ct][coff+22*i+6]  = DM_POLYTOPE_QUADRILATERAL;
396       ex->cone[ct][coff+22*i+7]  = 1;
397       ex->cone[ct][coff+22*i+8]  = 0;
398       ex->cone[ct][coff+22*i+9]  = i;
399       ex->cone[ct][coff+22*i+10] = DM_POLYTOPE_QUADRILATERAL;
400       ex->cone[ct][coff+22*i+11] = 1;
401       ex->cone[ct][coff+22*i+12] = 2;
402       ex->cone[ct][coff+22*i+13] = i;
403       ex->cone[ct][coff+22*i+14] = DM_POLYTOPE_QUADRILATERAL;
404       ex->cone[ct][coff+22*i+15] = 1;
405       ex->cone[ct][coff+22*i+16] = 1;
406       ex->cone[ct][coff+22*i+17] = i;
407       ex->cone[ct][coff+22*i+18] = DM_POLYTOPE_QUADRILATERAL;
408       ex->cone[ct][coff+22*i+19] = 1;
409       ex->cone[ct][coff+22*i+20] = 3;
410       ex->cone[ct][coff+22*i+21] = i;
411       ex->ornt[ct][ooff+6*i+0] = -2;
412       ex->ornt[ct][ooff+6*i+1] =  0;
413       ex->ornt[ct][ooff+6*i+2] =  0;
414       ex->ornt[ct][ooff+6*i+3] =  0;
415       ex->ornt[ct][ooff+6*i+4] =  0;
416       ex->ornt[ct][ooff+6*i+5] =  1;
417     }
418   }
419   /* Layers positions */
420   if (!ex->Nth) {
421     if (ex->symmetric) for (l = 0; l <= ex->layers; ++l) ex->layerPos[l] = (ex->thickness*l)/ex->layers - ex->thickness/2;
422     else               for (l = 0; l <= ex->layers; ++l) ex->layerPos[l] = (ex->thickness*l)/ex->layers;
423   } else {
424     ex->thickness   = 0.;
425     ex->layerPos[0] = 0.;
426     for (l = 0; l < ex->layers; ++l) {
427       const PetscReal t = ex->thicknesses[PetscMin(l, ex->Nth-1)];
428       ex->layerPos[l+1] = ex->layerPos[l] + t;
429       ex->thickness    += t;
430     }
431     if (ex->symmetric) for (l = 0; l <= ex->layers; ++l) ex->layerPos[l] -= ex->thickness/2.;
432   }
433   PetscFunctionReturn(0);
434 }
435 
436 static PetscErrorCode DMPlexTransformDestroy_Extrude(DMPlexTransform tr)
437 {
438   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
439   PetscInt                 ct;
440 
441   PetscFunctionBegin;
442   for (ct = 0; ct < DM_NUM_POLYTOPES; ++ct) {
443     PetscCall(PetscFree4(ex->target[ct], ex->size[ct], ex->cone[ct], ex->ornt[ct]));
444   }
445   PetscCall(PetscFree5(ex->Nt, ex->target, ex->size, ex->cone, ex->ornt));
446   PetscCall(PetscFree(ex->layerPos));
447   PetscCall(PetscFree(ex));
448   PetscFunctionReturn(0);
449 }
450 
451 static PetscErrorCode DMPlexTransformGetSubcellOrientation_Extrude(DMPlexTransform tr, DMPolytopeType sct, PetscInt sp, PetscInt so, DMPolytopeType tct, PetscInt r, PetscInt o, PetscInt *rnew, PetscInt *onew)
452 {
453   DMPlexTransform_Extrude *ex     = (DMPlexTransform_Extrude *) tr->data;
454   DMLabel                  trType = tr->trType;
455   PetscInt                 rt;
456 
457   PetscFunctionBeginHot;
458   *rnew = r;
459   *onew = DMPolytopeTypeComposeOrientation(tct, o, so);
460   if (!so) PetscFunctionReturn(0);
461   if (trType) {
462     PetscCall(DMLabelGetValue(tr->trType, sp, &rt));
463     if (rt >= 100) PetscFunctionReturn(0);
464   }
465   if (ex->useTensor) {
466     switch (sct) {
467       case DM_POLYTOPE_POINT: break;
468       case DM_POLYTOPE_SEGMENT:
469       switch (tct) {
470         case DM_POLYTOPE_SEGMENT: break;
471         case DM_POLYTOPE_SEG_PRISM_TENSOR:
472           *onew = DMPolytopeTypeComposeOrientation(tct, o, so ? -1 : 0);
473           break;
474         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
475       }
476       break;
477       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell type %s", DMPolytopeTypes[sct]);
478     }
479   } else {
480     switch (sct) {
481       case DM_POLYTOPE_POINT: break;
482       case DM_POLYTOPE_SEGMENT:
483       switch (tct) {
484         case DM_POLYTOPE_SEGMENT: break;
485         case DM_POLYTOPE_QUADRILATERAL:
486           *onew = DMPolytopeTypeComposeOrientation(tct, o, so ? -3 : 0);
487           break;
488         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
489       }
490       break;
491       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell type %s", DMPolytopeTypes[sct]);
492     }
493   }
494   PetscFunctionReturn(0);
495 }
496 
497 static PetscErrorCode DMPlexTransformCellTransform_Extrude(DMPlexTransform tr, DMPolytopeType source, PetscInt p, PetscInt *rt, PetscInt *Nt, DMPolytopeType *target[], PetscInt *size[], PetscInt *cone[], PetscInt *ornt[])
498 {
499   DMPlexTransform_Extrude *ex     = (DMPlexTransform_Extrude *) tr->data;
500   DMLabel                  trType = tr->trType;
501   PetscBool                ignore = PETSC_FALSE, identity = PETSC_FALSE;
502   PetscInt                 val    = 0;
503 
504   PetscFunctionBegin;
505   if (trType) {
506     PetscCall(DMLabelGetValue(trType, p, &val));
507     identity = val >= 100 ? PETSC_TRUE : PETSC_FALSE;
508   } else {
509     ignore = ex->Nt[source] < 0 ? PETSC_TRUE : PETSC_FALSE;
510   }
511   if (rt) *rt = val;
512   if (ignore) {
513     /* Ignore cells that cannot be extruded */
514     *Nt     = 0;
515     *target = NULL;
516     *size   = NULL;
517     *cone   = NULL;
518     *ornt   = NULL;
519   } else if (identity) {
520     PetscCall(DMPlexTransformCellTransformIdentity(tr, source, p, NULL, Nt, target, size, cone, ornt));
521   } else {
522     *Nt     = ex->Nt[source];
523     *target = ex->target[source];
524     *size   = ex->size[source];
525     *cone   = ex->cone[source];
526     *ornt   = ex->ornt[source];
527   }
528   PetscFunctionReturn(0);
529 }
530 
531 /* Computes new vertex along normal */
532 static PetscErrorCode DMPlexTransformMapCoordinates_Extrude(DMPlexTransform tr, DMPolytopeType pct, DMPolytopeType ct, PetscInt p, PetscInt r, PetscInt Nv, PetscInt dE, const PetscScalar in[], PetscScalar out[])
533 {
534   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
535   DM                       dm;
536   DMLabel                  active;
537   PetscReal                ones2[2]  = {0., 1.}, ones3[3] = { 0., 0., 1.};
538   PetscReal                normal[3] = {0., 0., 0.}, norm;
539   PetscBool                computeNormal;
540   PetscInt                 dim, dEx = ex->cdimEx, cStart, cEnd, d;
541 
542   PetscFunctionBeginHot;
543   PetscCheck(pct == DM_POLYTOPE_POINT,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parent point type %s",DMPolytopeTypes[pct]);
544   PetscCheck(ct  == DM_POLYTOPE_POINT,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for refined point type %s",DMPolytopeTypes[ct]);
545   PetscCheck(Nv == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertices should be produced from a single vertex, not %" PetscInt_FMT,Nv);
546   PetscCheck(dE == ex->cdim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Coordinate dim %" PetscInt_FMT " != %" PetscInt_FMT " original dimension", dE, ex->cdim);
547 
548   PetscCall(DMPlexTransformGetDM(tr, &dm));
549   PetscCall(DMGetDimension(dm, &dim));
550   PetscCall(DMPlexTransformGetActive(tr, &active));
551   computeNormal = dim != ex->cdim && !ex->useNormal ? PETSC_TRUE : PETSC_FALSE;
552   if (computeNormal) {
553     PetscInt *closure = NULL;
554     PetscInt  closureSize, cl;
555 
556     PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
557     PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_FALSE, &closureSize, &closure));
558     for (cl = 0; cl < closureSize*2; cl += 2) {
559       if ((closure[cl] >= cStart) && (closure[cl] < cEnd)) {
560         PetscReal cnormal[3] = {0, 0, 0};
561 
562         PetscCall(DMPlexComputeCellGeometryFVM(dm, closure[cl], NULL, NULL, cnormal));
563         for (d = 0; d < dEx; ++d) normal[d] += cnormal[d];
564       }
565     }
566     PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_FALSE, &closureSize, &closure));
567   } else if (ex->useNormal) {
568     for (d = 0; d < dEx; ++d) normal[d] = ex->normal[d];
569   } else if (active) { // find an active point in the closure of p and use its coordinate normal as the extrusion direction
570     PetscInt *closure = NULL;
571     PetscInt  closureSize, cl, pStart, pEnd;
572 
573     PetscCall(DMPlexGetDepthStratum(dm, ex->cdimEx-1, &pStart, &pEnd));
574     PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_FALSE, &closureSize, &closure));
575     for (cl = 0; cl < closureSize*2; cl += 2) {
576       if ((closure[cl] >= pStart) && (closure[cl] < pEnd)) {
577         PetscReal       cnormal[3] = {0, 0, 0};
578         const PetscInt *supp;
579         PetscInt        suppSize;
580 
581         PetscCall(DMPlexComputeCellGeometryFVM(dm, closure[cl], NULL, NULL, cnormal));
582         PetscCall(DMPlexGetSupportSize(dm, closure[cl], &suppSize));
583         PetscCall(DMPlexGetSupport(dm, closure[cl], &supp));
584         // Only use external faces, so I can get the orientation from any cell
585         if (suppSize == 1) {
586           const PetscInt *cone, *ornt;
587           PetscInt        coneSize, c;
588 
589           PetscCall(DMPlexGetConeSize(dm, supp[0], &coneSize));
590           PetscCall(DMPlexGetCone(dm, supp[0], &cone));
591           PetscCall(DMPlexGetConeOrientation(dm, supp[0], &ornt));
592           for (c = 0; c < coneSize; ++c) if (cone[c] == closure[cl]) break;
593           PetscCheck(c < coneSize, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Asymmetry in cone/support");
594           if (ornt[c] < 0) for (d = 0; d < dEx; ++d) cnormal[d] *= -1.;
595           for (d = 0; d < dEx; ++d) normal[d] += cnormal[d];
596         }
597       }
598     }
599     PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_FALSE, &closureSize, &closure));
600   } else if (ex->cdimEx == 2) {
601     for (d = 0; d < dEx; ++d) normal[d] = ones2[d];
602   } else if (ex->cdimEx == 3) {
603     for (d = 0; d < dEx; ++d) normal[d] = ones3[d];
604   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to determine normal for extrusion");
605   if (ex->normalFunc) {
606     PetscScalar n[3];
607     PetscReal   x[3], dot;
608 
609     for (d = 0; d < ex->cdim; ++d) x[d] = PetscRealPart(in[d]);
610     PetscCall((*ex->normalFunc)(ex->cdim, 0., x, r, n, NULL));
611     for (dot=0,d=0; d < dEx; d++) dot += PetscRealPart(n[d]) * normal[d];
612     for (d = 0; d < dEx; ++d) normal[d] = PetscSign(dot) * PetscRealPart(n[d]);
613   }
614 
615   for (d = 0, norm = 0.0; d < dEx; ++d) norm += PetscSqr(normal[d]);
616   for (d = 0; d < dEx; ++d) normal[d] *= norm == 0.0 ? 1.0 : 1./PetscSqrtReal(norm);
617   for (d = 0; d < dEx;      ++d) out[d]  = normal[d]*ex->layerPos[r];
618   for (d = 0; d < ex->cdim; ++d) out[d] += in[d];
619   PetscFunctionReturn(0);
620 }
621 
622 static PetscErrorCode DMPlexTransformInitialize_Extrude(DMPlexTransform tr)
623 {
624   PetscFunctionBegin;
625   tr->ops->view           = DMPlexTransformView_Extrude;
626   tr->ops->setfromoptions = DMPlexTransformSetFromOptions_Extrude;
627   tr->ops->setup          = DMPlexTransformSetUp_Extrude;
628   tr->ops->destroy        = DMPlexTransformDestroy_Extrude;
629   tr->ops->setdimensions  = DMPlexTransformSetDimensions_Extrude;
630   tr->ops->celltransform  = DMPlexTransformCellTransform_Extrude;
631   tr->ops->getsubcellorientation = DMPlexTransformGetSubcellOrientation_Extrude;
632   tr->ops->mapcoordinates = DMPlexTransformMapCoordinates_Extrude;
633   PetscFunctionReturn(0);
634 }
635 
636 PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Extrude(DMPlexTransform tr)
637 {
638   DMPlexTransform_Extrude *ex;
639   DM                       dm;
640   PetscInt                 dim;
641 
642   PetscFunctionBegin;
643   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
644   PetscCall(PetscNewLog(tr, &ex));
645   tr->data        = ex;
646   ex->thickness   = 1.;
647   ex->useTensor   = PETSC_TRUE;
648   ex->symmetric   = PETSC_FALSE;
649   ex->useNormal   = PETSC_FALSE;
650   ex->layerPos    = NULL;
651   PetscCall(DMPlexTransformGetDM(tr, &dm));
652   PetscCall(DMGetDimension(dm, &dim));
653   PetscCall(DMGetCoordinateDim(dm, &ex->cdim));
654   PetscCall(DMPlexTransformExtrudeSetLayers(tr, 1));
655   PetscCall(DMPlexTransformInitialize_Extrude(tr));
656   PetscFunctionReturn(0);
657 }
658 
659 /*@
660   DMPlexTransformExtrudeGetLayers - Get the number of extruded layers.
661 
662   Not collective
663 
664   Input Parameter:
665 . tr  - The DMPlexTransform
666 
667   Output Parameter:
668 . layers - The number of layers
669 
670   Level: intermediate
671 
672 .seealso: `DMPlexTransformExtrudeSetLayers()`
673 @*/
674 PetscErrorCode DMPlexTransformExtrudeGetLayers(DMPlexTransform tr, PetscInt *layers)
675 {
676   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
677 
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
680   PetscValidIntPointer(layers, 2);
681   *layers = ex->layers;
682   PetscFunctionReturn(0);
683 }
684 
685 /*@
686   DMPlexTransformExtrudeSetLayers - Set the number of extruded layers.
687 
688   Not collective
689 
690   Input Parameters:
691 + tr  - The DMPlexTransform
692 - layers - The number of layers
693 
694   Level: intermediate
695 
696 .seealso: `DMPlexTransformExtrudeGetLayers()`
697 @*/
698 PetscErrorCode DMPlexTransformExtrudeSetLayers(DMPlexTransform tr, PetscInt layers)
699 {
700   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
701 
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
704   ex->layers = layers;
705   PetscCall(PetscFree(ex->layerPos));
706   PetscCall(PetscCalloc1(ex->layers+1, &ex->layerPos));
707   PetscFunctionReturn(0);
708 }
709 
710 /*@
711   DMPlexTransformExtrudeGetThickness - Get the total thickness of the layers
712 
713   Not collective
714 
715   Input Parameter:
716 . tr  - The DMPlexTransform
717 
718   Output Parameter:
719 . thickness - The total thickness of the layers
720 
721   Level: intermediate
722 
723 .seealso: `DMPlexTransformExtrudeSetThickness()`
724 @*/
725 PetscErrorCode DMPlexTransformExtrudeGetThickness(DMPlexTransform tr, PetscReal *thickness)
726 {
727   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
728 
729   PetscFunctionBegin;
730   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
731   PetscValidRealPointer(thickness, 2);
732   *thickness = ex->thickness;
733   PetscFunctionReturn(0);
734 }
735 
736 /*@
737   DMPlexTransformExtrudeSetThickness - Set the total thickness of the layers
738 
739   Not collective
740 
741   Input Parameters:
742 + tr  - The DMPlexTransform
743 - thickness - The total thickness of the layers
744 
745   Level: intermediate
746 
747 .seealso: `DMPlexTransformExtrudeGetThickness()`
748 @*/
749 PetscErrorCode DMPlexTransformExtrudeSetThickness(DMPlexTransform tr, PetscReal thickness)
750 {
751   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
752 
753   PetscFunctionBegin;
754   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
755   PetscCheck(thickness > 0.,PetscObjectComm((PetscObject) tr), PETSC_ERR_ARG_OUTOFRANGE, "Height of layers %g must be positive", (double) thickness);
756   ex->thickness = thickness;
757   PetscFunctionReturn(0);
758 }
759 
760 /*@
761   DMPlexTransformExtrudeGetTensor - Get the flag to use tensor cells
762 
763   Not collective
764 
765   Input Parameter:
766 . tr  - The DMPlexTransform
767 
768   Output Parameter:
769 . useTensor - The flag to use tensor cells
770 
771   Note: This flag determines the orientation behavior of the created points. For example, if tensor is PETSC_TRUE, then DM_POLYTOPE_POINT_PRISM_TENSOR is made instead of DM_POLYTOPE_SEGMENT, DM_POLYTOPE_SEG_PRISM_TENSOR instead of DM_POLYTOPE_QUADRILATERAL, DM_POLYTOPE_TRI_PRISM_TENSOR instead of DM_POLYTOPE_TRI_PRISM, and DM_POLYTOPE_QUAD_PRISM_TENSOR instead of DM_POLYTOPE_HEXAHEDRON.
772 
773   Level: intermediate
774 
775 .seealso: `DMPlexTransformExtrudeSetTensor()`
776 @*/
777 PetscErrorCode DMPlexTransformExtrudeGetTensor(DMPlexTransform tr, PetscBool *useTensor)
778 {
779   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
780 
781   PetscFunctionBegin;
782   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
783   PetscValidBoolPointer(useTensor, 2);
784   *useTensor = ex->useTensor;
785   PetscFunctionReturn(0);
786 }
787 
788 /*@
789   DMPlexTransformExtrudeSetTensor - Set the flag to use tensor cells
790 
791   Not collective
792 
793   Input Parameters:
794 + tr  - The DMPlexTransform
795 - useTensor - The flag for tensor cells
796 
797   Note: This flag determines the orientation behavior of the created points. For example, if tensor is PETSC_TRUE, then DM_POLYTOPE_POINT_PRISM_TENSOR is made instead of DM_POLYTOPE_SEGMENT, DM_POLYTOPE_SEG_PRISM_TENSOR instead of DM_POLYTOPE_QUADRILATERAL, DM_POLYTOPE_TRI_PRISM_TENSOR instead of DM_POLYTOPE_TRI_PRISM, and DM_POLYTOPE_QUAD_PRISM_TENSOR instead of DM_POLYTOPE_HEXAHEDRON.
798 
799   Level: intermediate
800 
801 .seealso: `DMPlexTransformExtrudeGetTensor()`
802 @*/
803 PetscErrorCode DMPlexTransformExtrudeSetTensor(DMPlexTransform tr, PetscBool useTensor)
804 {
805   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
806 
807   PetscFunctionBegin;
808   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
809   ex->useTensor = useTensor;
810   PetscFunctionReturn(0);
811 }
812 
813 /*@
814   DMPlexTransformExtrudeGetSymmetric - Get the flag to extrude symmetrically from the initial surface
815 
816   Not collective
817 
818   Input Parameter:
819 . tr  - The DMPlexTransform
820 
821   Output Parameter:
822 . symmetric - The flag to extrude symmetrically
823 
824   Level: intermediate
825 
826 .seealso: `DMPlexTransformExtrudeSetSymmetric()`
827 @*/
828 PetscErrorCode DMPlexTransformExtrudeGetSymmetric(DMPlexTransform tr, PetscBool *symmetric)
829 {
830   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
831 
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
834   PetscValidBoolPointer(symmetric, 2);
835   *symmetric = ex->symmetric;
836   PetscFunctionReturn(0);
837 }
838 
839 /*@
840   DMPlexTransformExtrudeSetSymmetric - Set the flag to extrude symmetrically from the initial surface
841 
842   Not collective
843 
844   Input Parameters:
845 + tr  - The DMPlexTransform
846 - symmetric - The flag to extrude symmetrically
847 
848   Level: intermediate
849 
850 .seealso: `DMPlexTransformExtrudeGetSymmetric()`
851 @*/
852 PetscErrorCode DMPlexTransformExtrudeSetSymmetric(DMPlexTransform tr, PetscBool symmetric)
853 {
854   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
855 
856   PetscFunctionBegin;
857   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
858   ex->symmetric = symmetric;
859   PetscFunctionReturn(0);
860 }
861 
862 /*@
863   DMPlexTransformExtrudeGetNormal - Get the extrusion normal vector
864 
865   Not collective
866 
867   Input Parameter:
868 . tr  - The DMPlexTransform
869 
870   Output Parameter:
871 . normal - The extrusion direction
872 
873   Note: The user passes in an array, which is filled by the library.
874 
875   Level: intermediate
876 
877 .seealso: `DMPlexTransformExtrudeSetNormal()`
878 @*/
879 PetscErrorCode DMPlexTransformExtrudeGetNormal(DMPlexTransform tr, PetscReal normal[])
880 {
881   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
882   PetscInt                 d;
883 
884   PetscFunctionBegin;
885   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
886   if (ex->useNormal) {for (d = 0; d < ex->cdimEx; ++d) normal[d] = ex->normal[d];}
887   else               {for (d = 0; d < ex->cdimEx; ++d) normal[d] = 0.;}
888   PetscFunctionReturn(0);
889 }
890 
891 /*@
892   DMPlexTransformExtrudeSetNormal - Set the extrusion normal
893 
894   Not collective
895 
896   Input Parameters:
897 + tr     - The DMPlexTransform
898 - normal - The extrusion direction
899 
900   Level: intermediate
901 
902 .seealso: `DMPlexTransformExtrudeGetNormal()`
903 @*/
904 PetscErrorCode DMPlexTransformExtrudeSetNormal(DMPlexTransform tr, const PetscReal normal[])
905 {
906   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
907   PetscInt                 d;
908 
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
911   ex->useNormal = PETSC_TRUE;
912   for (d = 0; d < ex->cdimEx; ++d) ex->normal[d] = normal[d];
913   PetscFunctionReturn(0);
914 }
915 
916 /*@C
917   DMPlexTransformExtrudeSetNormalFunction - Set a function to determine the extrusion normal
918 
919   Not collective
920 
921   Input Parameters:
922 + tr     - The DMPlexTransform
923 - normalFunc - A function determining the extrusion direction
924 
925   Note: The calling sequence for the function is normalFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx)
926 $ dim  - The coordinate dimension of the original mesh (usually a surface)
927 $ time - The current time, or 0.
928 $ x    - The location of the current normal, in the coordinate space of the original mesh
929 $ r    - The extrusion replica number (layer number) of this point
930 $ u    - The user provides the computed normal on output; the sign and magnitude is not significant
931 $ ctx  - An optional user context
932 
933   Level: intermediate
934 
935 .seealso: `DMPlexTransformExtrudeGetNormal()`
936 @*/
937 PetscErrorCode DMPlexTransformExtrudeSetNormalFunction(DMPlexTransform tr, PetscSimplePointFunc normalFunc)
938 {
939   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
940 
941   PetscFunctionBegin;
942   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
943   ex->normalFunc = normalFunc;
944   PetscFunctionReturn(0);
945 }
946 
947 /*@
948   DMPlexTransformExtrudeSetThicknesses - Set the thickness of each layer
949 
950   Not collective
951 
952   Input Parameters:
953 + tr  - The DMPlexTransform
954 . Nth - The number of thicknesses
955 - thickness - The array of thicknesses
956 
957   Level: intermediate
958 
959 .seealso: `DMPlexTransformExtrudeSetThickness()`, `DMPlexTransformExtrudeGetThickness()`
960 @*/
961 PetscErrorCode DMPlexTransformExtrudeSetThicknesses(DMPlexTransform tr, PetscInt Nth, const PetscReal thicknesses[])
962 {
963   DMPlexTransform_Extrude *ex = (DMPlexTransform_Extrude *) tr->data;
964   PetscInt                 t;
965 
966   PetscFunctionBegin;
967   PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 1);
968   PetscCheck(Nth > 0,PetscObjectComm((PetscObject) tr), PETSC_ERR_ARG_OUTOFRANGE, "Number of thicknesses %" PetscInt_FMT " must be positive", Nth);
969   ex->Nth = PetscMin(Nth, ex->layers);
970   PetscCall(PetscFree(ex->thicknesses));
971   PetscCall(PetscMalloc1(ex->Nth, &ex->thicknesses));
972   for (t = 0; t < ex->Nth; ++t) {
973     PetscCheck(thicknesses[t] > 0.,PetscObjectComm((PetscObject) tr), PETSC_ERR_ARG_OUTOFRANGE, "Thickness %g of layer %" PetscInt_FMT " must be positive", (double) thicknesses[t], t);
974     ex->thicknesses[t] = thicknesses[t];
975   }
976   PetscFunctionReturn(0);
977 }
978