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