-
Notifications
You must be signed in to change notification settings - Fork 45
/
petscutils.c
376 lines (344 loc) · 15.2 KB
/
petscutils.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
#include "../include/petscutils.h"
// -----------------------------------------------------------------------------
// Convert PETSc MemType to libCEED MemType
// -----------------------------------------------------------------------------
CeedMemType MemTypeP2C(PetscMemType mem_type) {
return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
}
// -----------------------------------------------------------------------------
// Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
// -----------------------------------------------------------------------------
PetscErrorCode ProjectToUnitSphere(DM dm) {
Vec coordinates;
PetscScalar *coords;
PetscInt Nv, v, dim, d;
PetscErrorCode ierr;
PetscFunctionBeginUser;
ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
Nv /= dim;
ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
for (v = 0; v < Nv; ++v) {
PetscReal r = 0.0;
for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
r = PetscSqrtReal(r);
for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
}
ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
PetscFunctionReturn(0);
};
// -----------------------------------------------------------------------------
// Apply 3D Kershaw mesh transformation
// -----------------------------------------------------------------------------
// Transition from a value of "a" for x=0, to a value of "b" for x=1. Optionally
// smooth -- see the commented versions at the end.
static double step(const double a, const double b, double x) {
if (x <= 0) return a;
if (x >= 1) return b;
return a + (b-a) * (x);
}
// 1D transformation at the right boundary
static double right(const double eps, const double x) {
return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
}
// 1D transformation at the left boundary
static double left(const double eps, const double x) {
return 1-right(eps,1-x);
}
// Apply 3D Kershaw mesh transformation
// The eps parameters are in (0, 1]
// Uniform mesh is recovered for eps=1
PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
PetscErrorCode ierr;
Vec coord;
PetscInt ncoord;
PetscScalar *c;
PetscFunctionBeginUser;
ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
for (PetscInt i = 0; i < ncoord; i += 3) {
PetscScalar x = c[i], y = c[i+1], z = c[i+2];
PetscInt layer = x*6;
PetscScalar lambda = (x-layer/6.0)*6;
c[i] = x;
switch (layer) {
case 0:
c[i+1] = left(eps, y);
c[i+2] = left(eps, z);
break;
case 1:
case 4:
c[i+1] = step(left(eps, y), right(eps, y), lambda);
c[i+2] = step(left(eps, z), right(eps, z), lambda);
break;
case 2:
c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
break;
case 3:
c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
break;
default:
c[i+1] = right(eps, y);
c[i+2] = right(eps, z);
}
}
ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
// -----------------------------------------------------------------------------
// PETSc FE Boilerplate
// -----------------------------------------------------------------------------
PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc,
PetscBool isSimplex, const char prefix[],
PetscInt order, PetscFE *fem) {
PetscQuadrature q, fq;
DM K;
PetscSpace P;
PetscDualSpace Q;
PetscInt quadPointsPerEdge;
PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
PetscErrorCode ierr;
PetscFunctionBeginUser;
/* Create space */
ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P); CHKERRQ(ierr);
ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix); CHKERRQ(ierr);
ierr = PetscSpacePolynomialSetTensor(P, tensor); CHKERRQ(ierr);
ierr = PetscSpaceSetFromOptions(P); CHKERRQ(ierr);
ierr = PetscSpaceSetNumComponents(P, Nc); CHKERRQ(ierr);
ierr = PetscSpaceSetNumVariables(P, dim); CHKERRQ(ierr);
ierr = PetscSpaceSetDegree(P, order, order); CHKERRQ(ierr);
ierr = PetscSpaceSetUp(P); CHKERRQ(ierr);
ierr = PetscSpacePolynomialGetTensor(P, &tensor); CHKERRQ(ierr);
/* Create dual space */
ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);
CHKERRQ(ierr);
ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE); CHKERRQ(ierr);
ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix); CHKERRQ(ierr);
ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K); CHKERRQ(ierr);
ierr = PetscDualSpaceSetDM(Q, K); CHKERRQ(ierr);
ierr = DMDestroy(&K); CHKERRQ(ierr);
ierr = PetscDualSpaceSetNumComponents(Q, Nc); CHKERRQ(ierr);
ierr = PetscDualSpaceSetOrder(Q, order); CHKERRQ(ierr);
ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor); CHKERRQ(ierr);
ierr = PetscDualSpaceSetFromOptions(Q); CHKERRQ(ierr);
ierr = PetscDualSpaceSetUp(Q); CHKERRQ(ierr);
/* Create element */
ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem); CHKERRQ(ierr);
ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix); CHKERRQ(ierr);
ierr = PetscFESetFromOptions(*fem); CHKERRQ(ierr);
ierr = PetscFESetBasisSpace(*fem, P); CHKERRQ(ierr);
ierr = PetscFESetDualSpace(*fem, Q); CHKERRQ(ierr);
ierr = PetscFESetNumComponents(*fem, Nc); CHKERRQ(ierr);
ierr = PetscFESetUp(*fem); CHKERRQ(ierr);
ierr = PetscSpaceDestroy(&P); CHKERRQ(ierr);
ierr = PetscDualSpaceDestroy(&Q); CHKERRQ(ierr);
/* Create quadrature */
quadPointsPerEdge = PetscMax(order + 1,1);
if (isSimplex) {
ierr = PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0,
&q); CHKERRQ(ierr);
ierr = PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0,
&fq); CHKERRQ(ierr);
} else {
ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0,
&q); CHKERRQ(ierr);
ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0,
&fq); CHKERRQ(ierr);
}
ierr = PetscFESetQuadrature(*fem, q); CHKERRQ(ierr);
ierr = PetscFESetFaceQuadrature(*fem, fq); CHKERRQ(ierr);
ierr = PetscQuadratureDestroy(&q); CHKERRQ(ierr);
ierr = PetscQuadratureDestroy(&fq); CHKERRQ(ierr);
PetscFunctionReturn(0);
};
// -----------------------------------------------------------------------------
// Create BC label
// -----------------------------------------------------------------------------
static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
int ierr;
DMLabel label;
PetscFunctionBeginUser;
ierr = DMCreateLabel(dm, name); CHKERRQ(ierr);
ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr);
ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr);
ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr);
PetscFunctionReturn(0);
};
// -----------------------------------------------------------------------------
// This function sets up a DM for a given degree
// -----------------------------------------------------------------------------
PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u,
PetscInt dim, bool enforce_bc, BCFunction bc_func) {
PetscInt ierr, marker_ids[1] = {1};
PetscFE fe;
PetscFunctionBeginUser;
// Setup FE
ierr = PetscFECreateByDegree(dm, dim, num_comp_u, PETSC_FALSE, NULL, degree,
&fe);
CHKERRQ(ierr);
ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
// Setup DM
ierr = DMCreateDS(dm); CHKERRQ(ierr);
if (enforce_bc) {
PetscBool has_label;
DMHasLabel(dm, "marker", &has_label);
if (!has_label) {CreateBCLabel(dm, "marker");}
ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL,
(void(*)(void))bc_func, NULL, 1, marker_ids, NULL);
CHKERRQ(ierr);
}
ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
CHKERRQ(ierr);
ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
PetscFunctionReturn(0);
};
// -----------------------------------------------------------------------------
// This function sets up a BDDC vertex only DM from an existing fine DM
// -----------------------------------------------------------------------------
PetscErrorCode SetupVertexDMFromDM(DM dm, DM dm_vertex, PetscInt num_comp_u,
bool enforce_bc, BCFunction bc_func) {
PetscInt ierr, dim;
PetscFunctionBeginUser;
ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
ierr = SetupDMByDegree(dm_vertex, 1, num_comp_u, dim, enforce_bc, bc_func);
CHKERRQ(ierr);
PetscFunctionReturn(0);
};
// -----------------------------------------------------------------------------
// Utility function - essential BC dofs are encoded in closure indices as -(i+1)
// -----------------------------------------------------------------------------
PetscInt Involute(PetscInt i) {
return i >= 0 ? i : -(i + 1);
};
// -----------------------------------------------------------------------------
// Get CEED restriction data from DMPlex
// -----------------------------------------------------------------------------
PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
CeedInt topo_dim, CeedInt height, DMLabel domain_label, CeedInt value,
CeedElemRestriction *elem_restr) {
PetscSection section;
PetscInt p, num_elem, num_dof, *elem_restr_offsets, e_offset, num_fields, dim,
depth;
DMLabel depth_label;
IS depth_is, iter_is;
Vec U_loc;
const PetscInt *iter_indices;
PetscErrorCode ierr;
PetscFunctionBeginUser;
ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
dim -= height;
ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr);
ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr);
PetscInt num_comp[num_fields], field_off[num_fields+1];
field_off[0] = 0;
for (PetscInt f = 0; f < num_fields; f++) {
ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr);
field_off[f+1] = field_off[f] + num_comp[f];
}
ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is);
CHKERRQ(ierr);
if (domain_label) {
IS domain_is;
ierr = DMLabelGetStratumIS(domain_label, value, &domain_is); CHKERRQ(ierr);
if (domain_is) { // domain_is is non-empty
ierr = ISIntersect(depth_is, domain_is, &iter_is); CHKERRQ(ierr);
ierr = ISDestroy(&domain_is); CHKERRQ(ierr);
} else { // domain_is is NULL (empty)
iter_is = NULL;
}
ierr = ISDestroy(&depth_is); CHKERRQ(ierr);
} else {
iter_is = depth_is;
}
if (iter_is) {
ierr = ISGetLocalSize(iter_is, &num_elem); CHKERRQ(ierr);
ierr = ISGetIndices(iter_is, &iter_indices); CHKERRQ(ierr);
} else {
num_elem = 0;
iter_indices = NULL;
}
ierr = PetscMalloc1(num_elem*PetscPowInt(P, topo_dim), &elem_restr_offsets);
CHKERRQ(ierr);
for (p = 0, e_offset = 0; p < num_elem; p++) {
PetscInt c = iter_indices[p];
PetscInt num_indices, *indices, num_nodes;
ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
&num_indices, &indices, NULL, NULL);
CHKERRQ(ierr);
bool flip = false;
if (height > 0) {
PetscInt num_cells, num_faces, start = -1;
const PetscInt *orients, *faces, *cells;
ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr);
if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
"Expected one cell in support of exterior face, but got %D cells",
num_cells);
ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr);
for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
"Could not find face %D in cone of its support",
c);
ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
if (orients[start] < 0) flip = true;
}
if (num_indices % field_off[num_fields]) SETERRQ1(PETSC_COMM_SELF,
PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
c);
num_nodes = num_indices / field_off[num_fields];
for (PetscInt i = 0; i < num_nodes; i++) {
PetscInt ii = i;
if (flip) {
if (P == num_nodes) ii = num_nodes - 1 - i;
else if (P*P == num_nodes) {
PetscInt row = i / P, col = i % P;
ii = row + col * P;
} else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
"No support for flipping point with %D nodes != P (%D) or P^2",
num_nodes, P);
}
// Check that indices are blocked by node and thus can be coalesced as a single field with
// field_off[num_fields] = sum(num_comp) components.
for (PetscInt f = 0; f < num_fields; f++) {
for (PetscInt j = 0; j < num_comp[f]; j++) {
if (Involute(indices[field_off[f]*num_nodes + ii*num_comp[f] + j])
!= Involute(indices[ii*num_comp[0]]) + field_off[f] + j)
SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
"Cell %D closure indices not interlaced for node %D field %D component %D",
c, ii, f, j);
}
}
// Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
PetscInt loc = Involute(indices[ii*num_comp[0]]);
elem_restr_offsets[e_offset++] = loc;
}
ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
&num_indices, &indices, NULL, NULL);
CHKERRQ(ierr);
}
if (e_offset != num_elem*PetscPowInt(P, topo_dim))
SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
"ElemRestriction of size (%D,%D) initialized %D nodes", num_elem,
PetscPowInt(P, topo_dim),e_offset);
if (iter_is) {
ierr = ISRestoreIndices(iter_is, &iter_indices); CHKERRQ(ierr);
}
ierr = ISDestroy(&iter_is); CHKERRQ(ierr);
ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr);
ierr = VecGetLocalSize(U_loc, &num_dof); CHKERRQ(ierr);
ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr);
CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, topo_dim),
field_off[num_fields], 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
elem_restr_offsets, elem_restr);
ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
PetscFunctionReturn(0);
};
// -----------------------------------------------------------------------------