1 #include <petscdm.h> 2 #include <petscdmplex.h> 3 #include <petscdmswarm.h> 4 #include "../src/dm/impls/swarm/data_bucket.h" 5 6 PetscBool SwarmProjcite = PETSC_FALSE; 7 const char SwarmProjCitation[] = "@article{PusztayKnepleyAdams2022,\n" 8 "title = {Conservative Projection Between FEM and Particle Bases},\n" 9 "author = {Joseph V. Pusztay and Matthew G. Knepley and Mark F. Adams},\n" 10 "journal = {SIAM Journal on Scientific Computing},\n" 11 "volume = {44},\n" 12 "number = {4},\n" 13 "pages = {C310--C319},\n" 14 "doi = {10.1137/21M145407},\n" 15 "year = {2022}\n}\n"; 16 17 PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *xi); 18 19 static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem) 20 { 21 const PetscInt Nc = 1; 22 PetscQuadrature q, fq; 23 DM K; 24 PetscSpace P; 25 PetscDualSpace Q; 26 PetscInt order, quadPointsPerEdge; 27 PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 28 29 PetscFunctionBegin; 30 /* Create space */ 31 PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)dm), &P)); 32 /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); */ 33 PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); 34 /* PetscCall(PetscSpaceSetFromOptions(P)); */ 35 PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 36 PetscCall(PetscSpaceSetDegree(P, 1, PETSC_DETERMINE)); 37 PetscCall(PetscSpaceSetNumComponents(P, Nc)); 38 PetscCall(PetscSpaceSetNumVariables(P, dim)); 39 PetscCall(PetscSpaceSetUp(P)); 40 PetscCall(PetscSpaceGetDegree(P, &order, NULL)); 41 PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); 42 /* Create dual space */ 43 PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)dm), &Q)); 44 PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 45 /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); */ 46 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); 47 PetscCall(PetscDualSpaceSetDM(Q, K)); 48 PetscCall(DMDestroy(&K)); 49 PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 50 PetscCall(PetscDualSpaceSetOrder(Q, order)); 51 PetscCall(PetscDualSpaceLagrangeSetTensor(Q, tensor)); 52 /* PetscCall(PetscDualSpaceSetFromOptions(Q)); */ 53 PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 54 PetscCall(PetscDualSpaceSetUp(Q)); 55 /* Create element */ 56 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)dm), fem)); 57 /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); */ 58 /* PetscCall(PetscFESetFromOptions(*fem)); */ 59 PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 60 PetscCall(PetscFESetBasisSpace(*fem, P)); 61 PetscCall(PetscFESetDualSpace(*fem, Q)); 62 PetscCall(PetscFESetNumComponents(*fem, Nc)); 63 PetscCall(PetscFESetUp(*fem)); 64 PetscCall(PetscSpaceDestroy(&P)); 65 PetscCall(PetscDualSpaceDestroy(&Q)); 66 /* Create quadrature (with specified order if given) */ 67 qorder = qorder >= 0 ? qorder : order; 68 quadPointsPerEdge = PetscMax(qorder + 1, 1); 69 if (isSimplex) { 70 PetscCall(PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); 71 PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); 72 } else { 73 PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); 74 PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); 75 } 76 PetscCall(PetscFESetQuadrature(*fem, q)); 77 PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 78 PetscCall(PetscQuadratureDestroy(&q)); 79 PetscCall(PetscQuadratureDestroy(&fq)); 80 PetscFunctionReturn(PETSC_SUCCESS); 81 } 82 83 static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm, DM dmc, PetscInt nsub) 84 { 85 PetscInt dim, nfaces, nbasis; 86 PetscInt q, npoints_q, e, nel, pcnt, ps, pe, d, k, r; 87 PetscTabulation T; 88 Vec coorlocal; 89 PetscSection coordSection; 90 PetscScalar *elcoor = NULL; 91 PetscReal *swarm_coor; 92 PetscInt *swarm_cellid; 93 const PetscReal *xiq; 94 PetscQuadrature quadrature; 95 PetscFE fe, feRef; 96 PetscBool is_simplex; 97 98 PetscFunctionBegin; 99 PetscCall(DMGetDimension(dmc, &dim)); 100 is_simplex = PETSC_FALSE; 101 PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 102 PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 103 if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 104 105 PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 106 107 for (r = 0; r < nsub; r++) { 108 PetscCall(PetscFERefine(fe, &feRef)); 109 PetscCall(PetscFECopyQuadrature(feRef, fe)); 110 PetscCall(PetscFEDestroy(&feRef)); 111 } 112 113 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 114 PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL)); 115 PetscCall(PetscFEGetDimension(fe, &nbasis)); 116 PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 117 118 /* 0->cell, 1->edge, 2->vert */ 119 PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 120 nel = pe - ps; 121 122 PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 123 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 124 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 125 126 PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 127 PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 128 129 pcnt = 0; 130 for (e = 0; e < nel; e++) { 131 PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 132 133 for (q = 0; q < npoints_q; q++) { 134 for (d = 0; d < dim; d++) { 135 swarm_coor[dim * pcnt + d] = 0.0; 136 for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][q * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); 137 } 138 swarm_cellid[pcnt] = e; 139 pcnt++; 140 } 141 PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 142 } 143 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 144 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 145 146 PetscCall(PetscFEDestroy(&fe)); 147 PetscFunctionReturn(PETSC_SUCCESS); 148 } 149 150 static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm, DM dmc, PetscInt npoints) 151 { 152 PetscInt dim; 153 PetscInt ii, jj, q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, nfaces; 154 PetscReal *xi, ds, ds2; 155 PetscReal **basis; 156 Vec coorlocal; 157 PetscSection coordSection; 158 PetscScalar *elcoor = NULL; 159 PetscReal *swarm_coor; 160 PetscInt *swarm_cellid; 161 PetscBool is_simplex; 162 163 PetscFunctionBegin; 164 PetscCall(DMGetDimension(dmc, &dim)); 165 PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only 2D is supported"); 166 is_simplex = PETSC_FALSE; 167 PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 168 PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 169 if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 170 PetscCheck(is_simplex, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only the simplex is supported"); 171 172 PetscCall(PetscMalloc1(dim * npoints * npoints, &xi)); 173 pcnt = 0; 174 ds = 1.0 / ((PetscReal)(npoints - 1)); 175 ds2 = 1.0 / ((PetscReal)(npoints)); 176 for (jj = 0; jj < npoints; jj++) { 177 for (ii = 0; ii < npoints - jj; ii++) { 178 xi[dim * pcnt + 0] = ii * ds; 179 xi[dim * pcnt + 1] = jj * ds; 180 181 xi[dim * pcnt + 0] *= (1.0 - 1.2 * ds2); 182 xi[dim * pcnt + 1] *= (1.0 - 1.2 * ds2); 183 184 xi[dim * pcnt + 0] += 0.35 * ds2; 185 xi[dim * pcnt + 1] += 0.35 * ds2; 186 pcnt++; 187 } 188 } 189 npoints_q = pcnt; 190 191 npe = 3; /* nodes per element (triangle) */ 192 PetscCall(PetscMalloc1(npoints_q, &basis)); 193 for (q = 0; q < npoints_q; q++) { 194 PetscCall(PetscMalloc1(npe, &basis[q])); 195 196 basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1]; 197 basis[q][1] = xi[dim * q + 0]; 198 basis[q][2] = xi[dim * q + 1]; 199 } 200 201 /* 0->cell, 1->edge, 2->vert */ 202 PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 203 nel = pe - ps; 204 205 PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 206 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 207 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 208 209 PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 210 PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 211 212 pcnt = 0; 213 for (e = 0; e < nel; e++) { 214 PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 215 216 for (q = 0; q < npoints_q; q++) { 217 for (d = 0; d < dim; d++) { 218 swarm_coor[dim * pcnt + d] = 0.0; 219 for (k = 0; k < npe; k++) swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]); 220 } 221 swarm_cellid[pcnt] = e; 222 pcnt++; 223 } 224 PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 225 } 226 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 227 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 228 229 PetscCall(PetscFree(xi)); 230 for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q])); 231 PetscCall(PetscFree(basis)); 232 PetscFunctionReturn(PETSC_SUCCESS); 233 } 234 235 PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param) 236 { 237 PetscInt dim; 238 239 PetscFunctionBegin; 240 PetscCall(DMGetDimension(celldm, &dim)); 241 switch (layout) { 242 case DMSWARMPIC_LAYOUT_REGULAR: 243 PetscCheck(dim != 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No 3D support for REGULAR+PLEX"); 244 PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm, celldm, layout_param)); 245 break; 246 case DMSWARMPIC_LAYOUT_GAUSS: { 247 PetscQuadrature quad, facequad; 248 const PetscReal *xi; 249 DMPolytopeType ct; 250 PetscInt cStart, Nq; 251 252 PetscCall(DMPlexGetHeightStratum(celldm, 0, &cStart, NULL)); 253 PetscCall(DMPlexGetCellType(celldm, cStart, &ct)); 254 PetscCall(PetscDTCreateDefaultQuadrature(ct, layout_param, &quad, &facequad)); 255 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xi, NULL)); 256 PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, Nq, (PetscReal *)xi)); 257 PetscCall(PetscQuadratureDestroy(&quad)); 258 PetscCall(PetscQuadratureDestroy(&facequad)); 259 } break; 260 case DMSWARMPIC_LAYOUT_SUBDIVISION: 261 PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm, celldm, layout_param)); 262 break; 263 } 264 PetscFunctionReturn(PETSC_SUCCESS); 265 } 266 267 /* 268 typedef struct { 269 PetscReal x,y; 270 } Point2d; 271 272 static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 273 { 274 PetscFunctionBegin; 275 *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 276 PetscFunctionReturn(PETSC_SUCCESS); 277 } 278 */ 279 /* 280 static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 281 { 282 PetscReal s1,s2,s3; 283 PetscBool b1, b2, b3; 284 285 PetscFunctionBegin; 286 signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 287 signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 288 signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 289 290 *v = ((b1 == b2) && (b2 == b3)); 291 PetscFunctionReturn(PETSC_SUCCESS); 292 } 293 */ 294 /* 295 static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 296 { 297 PetscReal x1,y1,x2,y2,x3,y3; 298 PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 299 300 PetscFunctionBegin; 301 x1 = coords[2*0+0]; 302 x2 = coords[2*1+0]; 303 x3 = coords[2*2+0]; 304 305 y1 = coords[2*0+1]; 306 y2 = coords[2*1+1]; 307 y3 = coords[2*2+1]; 308 309 c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 310 b[0] = xp[0] - c; 311 c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 312 b[1] = xp[1] - c; 313 314 A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 315 A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 316 317 detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 318 *dJ = PetscAbsReal(detJ); 319 od = 1.0/detJ; 320 321 inv[0][0] = A[1][1] * od; 322 inv[0][1] = -A[0][1] * od; 323 inv[1][0] = -A[1][0] * od; 324 inv[1][1] = A[0][0] * od; 325 326 xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 327 xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 328 PetscFunctionReturn(PETSC_SUCCESS); 329 } 330 */ 331 332 static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[], PetscScalar coords[], PetscReal xip[], PetscReal *dJ) 333 { 334 PetscReal x1, y1, x2, y2, x3, y3; 335 PetscReal b[2], A[2][2], inv[2][2], detJ, od; 336 337 PetscFunctionBegin; 338 x1 = PetscRealPart(coords[2 * 0 + 0]); 339 x2 = PetscRealPart(coords[2 * 1 + 0]); 340 x3 = PetscRealPart(coords[2 * 2 + 0]); 341 342 y1 = PetscRealPart(coords[2 * 0 + 1]); 343 y2 = PetscRealPart(coords[2 * 1 + 1]); 344 y3 = PetscRealPart(coords[2 * 2 + 1]); 345 346 b[0] = xp[0] - x1; 347 b[1] = xp[1] - y1; 348 349 A[0][0] = x2 - x1; 350 A[0][1] = x3 - x1; 351 A[1][0] = y2 - y1; 352 A[1][1] = y3 - y1; 353 354 detJ = A[0][0] * A[1][1] - A[0][1] * A[1][0]; 355 *dJ = PetscAbsReal(detJ); 356 od = 1.0 / detJ; 357 358 inv[0][0] = A[1][1] * od; 359 inv[0][1] = -A[0][1] * od; 360 inv[1][0] = -A[1][0] * od; 361 inv[1][1] = A[0][0] * od; 362 363 xip[0] = inv[0][0] * b[0] + inv[0][1] * b[1]; 364 xip[1] = inv[1][0] * b[0] + inv[1][1] * b[1]; 365 PetscFunctionReturn(PETSC_SUCCESS); 366 } 367 368 static PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm, PetscReal *swarm_field, DM dm, Vec v_field) 369 { 370 const PetscReal PLEX_C_EPS = 1.0e-8; 371 Vec v_field_l, denom_l, coor_l, denom; 372 PetscInt k, p, e, npoints; 373 PetscInt *mpfield_cell; 374 PetscReal *mpfield_coor; 375 PetscReal xi_p[2]; 376 PetscScalar Ni[3]; 377 PetscSection coordSection; 378 PetscScalar *elcoor = NULL; 379 380 PetscFunctionBegin; 381 PetscCall(VecZeroEntries(v_field)); 382 383 PetscCall(DMGetLocalVector(dm, &v_field_l)); 384 PetscCall(DMGetGlobalVector(dm, &denom)); 385 PetscCall(DMGetLocalVector(dm, &denom_l)); 386 PetscCall(VecZeroEntries(v_field_l)); 387 PetscCall(VecZeroEntries(denom)); 388 PetscCall(VecZeroEntries(denom_l)); 389 390 PetscCall(DMGetCoordinatesLocal(dm, &coor_l)); 391 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 392 393 PetscCall(DMSwarmGetLocalSize(swarm, &npoints)); 394 PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 395 PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 396 397 for (p = 0; p < npoints; p++) { 398 PetscReal *coor_p, dJ; 399 PetscScalar elfield[3]; 400 PetscBool point_located; 401 402 e = mpfield_cell[p]; 403 coor_p = &mpfield_coor[2 * p]; 404 405 PetscCall(DMPlexVecGetClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 406 407 /* 408 while (!point_located && (failed_counter < 25)) { 409 PetscCall(PointInTriangle(point, coords[0], coords[1], coords[2], &point_located)); 410 point.x = coor_p[0]; 411 point.y = coor_p[1]; 412 point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 413 point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 414 failed_counter++; 415 } 416 417 if (!point_located) { 418 PetscPrintf(PETSC_COMM_SELF,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ") with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e) in %" PetscInt_FMT " iterations\n",point.x,point.y,e,coords[0].x,coords[0].y,coords[1].x,coords[1].y,coords[2].x,coords[2].y,failed_counter); 419 } 420 421 PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ")",point.x,point.y,e); 422 else { 423 PetscCall(_ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ)); 424 xi_p[0] = 0.5*(xi_p[0] + 1.0); 425 xi_p[1] = 0.5*(xi_p[1] + 1.0); 426 427 PetscPrintf(PETSC_COMM_SELF,"[p=%" PetscInt_FMT "] x(%+1.4e,%+1.4e) -> mapped to element %" PetscInt_FMT " xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]); 428 429 } 430 */ 431 432 PetscCall(ComputeLocalCoordinateAffine2d(coor_p, elcoor, xi_p, &dJ)); 433 /* 434 PetscPrintf(PETSC_COMM_SELF,"[p=%" PetscInt_FMT "] x(%+1.4e,%+1.4e) -> mapped to element %" PetscInt_FMT " xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]); 435 */ 436 /* 437 point_located = PETSC_TRUE; 438 if (xi_p[0] < 0.0) { 439 if (xi_p[0] > -PLEX_C_EPS) { 440 xi_p[0] = 0.0; 441 } else { 442 point_located = PETSC_FALSE; 443 } 444 } 445 if (xi_p[1] < 0.0) { 446 if (xi_p[1] > -PLEX_C_EPS) { 447 xi_p[1] = 0.0; 448 } else { 449 point_located = PETSC_FALSE; 450 } 451 } 452 if (xi_p[1] > (1.0-xi_p[0])) { 453 if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 454 xi_p[1] = 1.0 - xi_p[0]; 455 } else { 456 point_located = PETSC_FALSE; 457 } 458 } 459 if (!point_located) { 460 PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 461 PetscPrintf(PETSC_COMM_SELF,"[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ") with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e)\n",coor_p[0],coor_p[1],e,elcoor[0],elcoor[1],elcoor[2],elcoor[3],elcoor[4],elcoor[5]); 462 } 463 PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ")",coor_p[0],coor_p[1],e); 464 */ 465 466 Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 467 Ni[1] = xi_p[0]; 468 Ni[2] = xi_p[1]; 469 470 point_located = PETSC_TRUE; 471 for (k = 0; k < 3; k++) { 472 if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 473 if (PetscRealPart(Ni[k]) > (1.0 + PLEX_C_EPS)) point_located = PETSC_FALSE; 474 } 475 if (!point_located) { 476 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[Error] xi,eta = %+1.8e, %+1.8e\n", (double)xi_p[0], (double)xi_p[1])); 477 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ") with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e)\n", (double)coor_p[0], (double)coor_p[1], e, (double)PetscRealPart(elcoor[0]), (double)PetscRealPart(elcoor[1]), (double)PetscRealPart(elcoor[2]), (double)PetscRealPart(elcoor[3]), (double)PetscRealPart(elcoor[4]), (double)PetscRealPart(elcoor[5]))); 478 } 479 PetscCheck(point_located, PETSC_COMM_SELF, PETSC_ERR_SUP, "Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ")", (double)coor_p[0], (double)coor_p[1], e); 480 481 for (k = 0; k < 3; k++) { 482 Ni[k] = Ni[k] * dJ; 483 elfield[k] = Ni[k] * swarm_field[p]; 484 } 485 PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 486 487 PetscCall(DMPlexVecSetClosure(dm, NULL, v_field_l, e, elfield, ADD_VALUES)); 488 PetscCall(DMPlexVecSetClosure(dm, NULL, denom_l, e, Ni, ADD_VALUES)); 489 } 490 491 PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 492 PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 493 494 PetscCall(DMLocalToGlobalBegin(dm, v_field_l, ADD_VALUES, v_field)); 495 PetscCall(DMLocalToGlobalEnd(dm, v_field_l, ADD_VALUES, v_field)); 496 PetscCall(DMLocalToGlobalBegin(dm, denom_l, ADD_VALUES, denom)); 497 PetscCall(DMLocalToGlobalEnd(dm, denom_l, ADD_VALUES, denom)); 498 499 PetscCall(VecPointwiseDivide(v_field, v_field, denom)); 500 501 PetscCall(DMRestoreLocalVector(dm, &v_field_l)); 502 PetscCall(DMRestoreLocalVector(dm, &denom_l)); 503 PetscCall(DMRestoreGlobalVector(dm, &denom)); 504 PetscFunctionReturn(PETSC_SUCCESS); 505 } 506 507 PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]) 508 { 509 PetscInt f, dim; 510 511 PetscFunctionBegin; 512 PetscCall(DMGetDimension(swarm, &dim)); 513 switch (dim) { 514 case 2: 515 for (f = 0; f < nfields; f++) { 516 PetscReal *swarm_field; 517 518 PetscCall(DMSwarmDataFieldGetEntries(dfield[f], (void **)&swarm_field)); 519 PetscCall(DMSwarmProjectField_ApproxP1_PLEX_2D(swarm, swarm_field, celldm, vecs[f])); 520 } 521 break; 522 case 3: 523 SETERRQ(PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "No support for 3D"); 524 default: 525 break; 526 } 527 PetscFunctionReturn(PETSC_SUCCESS); 528 } 529 530 PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[]) 531 { 532 PetscBool is_simplex, is_tensorcell; 533 PetscInt dim, nfaces, ps, pe, p, d, nbasis, pcnt, e, k, nel; 534 PetscFE fe; 535 PetscQuadrature quadrature; 536 PetscTabulation T; 537 PetscReal *xiq; 538 Vec coorlocal; 539 PetscSection coordSection; 540 PetscScalar *elcoor = NULL; 541 PetscReal *swarm_coor; 542 PetscInt *swarm_cellid; 543 544 PetscFunctionBegin; 545 PetscCall(DMGetDimension(dmc, &dim)); 546 547 is_simplex = PETSC_FALSE; 548 is_tensorcell = PETSC_FALSE; 549 PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 550 PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 551 552 if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 553 554 switch (dim) { 555 case 2: 556 if (nfaces == 4) is_tensorcell = PETSC_TRUE; 557 break; 558 case 3: 559 if (nfaces == 6) is_tensorcell = PETSC_TRUE; 560 break; 561 default: 562 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D"); 563 } 564 565 /* check points provided fail inside the reference cell */ 566 if (is_simplex) { 567 for (p = 0; p < npoints; p++) { 568 PetscReal sum; 569 for (d = 0; d < dim; d++) PetscCheck(xi[dim * p + d] >= -1.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain"); 570 sum = 0.0; 571 for (d = 0; d < dim; d++) sum += xi[dim * p + d]; 572 PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain"); 573 } 574 } else if (is_tensorcell) { 575 for (p = 0; p < npoints; p++) { 576 for (d = 0; d < dim; d++) PetscCheck(PetscAbsReal(xi[dim * p + d]) <= 1.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the tensor domain [-1,1]^d"); 577 } 578 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell"); 579 580 PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm), &quadrature)); 581 PetscCall(PetscMalloc1(npoints * dim, &xiq)); 582 PetscCall(PetscArraycpy(xiq, xi, npoints * dim)); 583 PetscCall(PetscQuadratureSetData(quadrature, dim, 1, npoints, (const PetscReal *)xiq, NULL)); 584 PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 585 PetscCall(PetscFESetQuadrature(fe, quadrature)); 586 PetscCall(PetscFEGetDimension(fe, &nbasis)); 587 PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 588 589 /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 590 /* 0->cell, 1->edge, 2->vert */ 591 PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 592 nel = pe - ps; 593 594 PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, -1)); 595 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 596 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 597 598 PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 599 PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 600 601 pcnt = 0; 602 for (e = 0; e < nel; e++) { 603 PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 604 605 for (p = 0; p < npoints; p++) { 606 for (d = 0; d < dim; d++) { 607 swarm_coor[dim * pcnt + d] = 0.0; 608 for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][p * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); 609 } 610 swarm_cellid[pcnt] = e; 611 pcnt++; 612 } 613 PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 614 } 615 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 616 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 617 618 PetscCall(PetscQuadratureDestroy(&quadrature)); 619 PetscCall(PetscFEDestroy(&fe)); 620 PetscFunctionReturn(PETSC_SUCCESS); 621 } 622