1 #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/ 2 #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/ 3 #include <petscfe.h> 4 #include <petscdmplex.h> 5 #include <petscds.h> 6 7 typedef struct _n_DMField_DS 8 { 9 PetscInt fieldNum; 10 Vec vec; 11 PetscInt height; 12 PetscObject *disc; 13 PetscBool multifieldVec; 14 } 15 DMField_DS; 16 17 static PetscErrorCode DMFieldDestroy_DS(DMField field) 18 { 19 DMField_DS *dsfield; 20 PetscInt i; 21 PetscErrorCode ierr; 22 23 PetscFunctionBegin; 24 dsfield = (DMField_DS *) field->data; 25 ierr = VecDestroy(&dsfield->vec);CHKERRQ(ierr); 26 for (i = 0; i < dsfield->height; i++) { 27 ierr = PetscObjectDereference(dsfield->disc[i]);CHKERRQ(ierr); 28 } 29 ierr = PetscFree(dsfield->disc);CHKERRQ(ierr); 30 ierr = PetscFree(dsfield);CHKERRQ(ierr); 31 PetscFunctionReturn(0); 32 } 33 34 static PetscErrorCode DMFieldView_DS(DMField field,PetscViewer viewer) 35 { 36 DMField_DS *dsfield = (DMField_DS *) field->data; 37 PetscBool iascii; 38 PetscObject disc; 39 PetscErrorCode ierr; 40 41 PetscFunctionBegin; 42 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 43 disc = dsfield->disc[0]; 44 if (iascii) { 45 PetscViewerASCIIPrintf(viewer, "PetscDS field %D\n",dsfield->fieldNum);CHKERRQ(ierr); 46 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 47 ierr = PetscObjectView(disc,viewer);CHKERRQ(ierr); 48 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 49 } 50 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 51 if (dsfield->multifieldVec) { 52 SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"View of subfield not implemented yet"); 53 } else { 54 ierr = VecView(dsfield->vec,viewer);CHKERRQ(ierr); 55 } 56 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 57 PetscFunctionReturn(0); 58 } 59 60 static PetscErrorCode DMFieldDSGetHeightDisc(DMField field, PetscInt height, PetscObject *disc) 61 { 62 DMField_DS *dsfield = (DMField_DS *) field->data; 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 if (!dsfield->disc[height]) { 67 PetscClassId id; 68 69 ierr = PetscObjectGetClassId(dsfield->disc[0],&id);CHKERRQ(ierr); 70 if (id == PETSCFE_CLASSID) { 71 PetscFE fe = (PetscFE) dsfield->disc[0]; 72 73 ierr = PetscFECreateHeightTrace(fe,height,(PetscFE *)&dsfield->disc[height]);CHKERRQ(ierr); 74 } 75 } 76 *disc = dsfield->disc[height]; 77 PetscFunctionReturn(0); 78 } 79 80 #define DMFieldDSdot(y,A,b,m,n,c,cast) \ 81 do { \ 82 PetscInt _i, _j, _k; \ 83 for (_i = 0; _i < (m); _i++) { \ 84 for (_k = 0; _k < (c); _k++) { \ 85 (y)[_i * (c) + _k] = 0.; \ 86 } \ 87 for (_j = 0; _j < (n); _j++) { \ 88 for (_k = 0; _k < (c); _k++) { \ 89 (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \ 90 } \ 91 } \ 92 } \ 93 } while (0) 94 95 static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H) 96 { 97 DMField_DS *dsfield = (DMField_DS *) field->data; 98 DM dm; 99 PetscObject disc; 100 PetscClassId classid; 101 PetscInt nq, nc, dim, meshDim, numCells; 102 PetscSection section; 103 const PetscReal *qpoints; 104 PetscBool isStride; 105 const PetscInt *points = NULL; 106 PetscInt sfirst = -1, stride = -1; 107 PetscErrorCode ierr; 108 109 PetscFunctionBeginHot; 110 dm = field->dm; 111 nc = field->numComponents; 112 ierr = PetscQuadratureGetData(quad,&dim,NULL,&nq,&qpoints,NULL);CHKERRQ(ierr); 113 ierr = DMFieldDSGetHeightDisc(field,dsfield->height - 1 - dim,&disc);CHKERRQ(ierr); 114 ierr = DMGetDimension(dm,&meshDim);CHKERRQ(ierr); 115 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 116 ierr = PetscSectionGetField(section,dsfield->fieldNum,§ion);CHKERRQ(ierr); 117 ierr = PetscObjectGetClassId(disc,&classid);CHKERRQ(ierr); 118 /* TODO: batch */ 119 ierr = PetscObjectTypeCompare((PetscObject)pointIS,ISSTRIDE,&isStride);CHKERRQ(ierr); 120 ierr = ISGetLocalSize(pointIS,&numCells);CHKERRQ(ierr); 121 if (isStride) { 122 ierr = ISStrideGetInfo(pointIS,&sfirst,&stride);CHKERRQ(ierr); 123 } else { 124 ierr = ISGetIndices(pointIS,&points);CHKERRQ(ierr); 125 } 126 if (classid == PETSCFE_CLASSID) { 127 PetscFE fe = (PetscFE) disc; 128 PetscInt feDim, i; 129 PetscReal *fB = NULL, *fD = NULL, *fH = NULL; 130 131 ierr = PetscFEGetDimension(fe,&feDim);CHKERRQ(ierr); 132 ierr = PetscFEGetTabulation(fe,nq,qpoints,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 133 for (i = 0; i < numCells; i++) { 134 PetscInt c = isStride ? (sfirst + i * stride) : points[i]; 135 PetscInt closureSize; 136 PetscScalar *elem = NULL; 137 138 ierr = DMPlexVecGetClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 139 if (B) { 140 if (type == PETSC_SCALAR) { 141 PetscScalar *cB = &((PetscScalar *) B)[nc * nq * i]; 142 143 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,(PetscScalar)); 144 } else { 145 PetscReal *cB = &((PetscReal *) B)[nc * nq * i]; 146 147 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,PetscRealPart); 148 } 149 } 150 if (D) { 151 if (type == PETSC_SCALAR) { 152 PetscScalar *cD = &((PetscScalar *) D)[nc * nq * dim * i]; 153 154 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),(PetscScalar)); 155 } else { 156 PetscReal *cD = &((PetscReal *) D)[nc * nq * dim * i]; 157 158 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),PetscRealPart); 159 } 160 } 161 if (H) { 162 if (type == PETSC_SCALAR) { 163 PetscScalar *cH = &((PetscScalar *) H)[nc * nq * dim * dim * i]; 164 165 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),(PetscScalar)); 166 } else { 167 PetscReal *cH = &((PetscReal *) H)[nc * nq * dim * dim * i]; 168 169 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),PetscRealPart); 170 } 171 } 172 ierr = DMPlexVecRestoreClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 173 } 174 ierr = PetscFERestoreTabulation(fe,nq,qpoints,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 175 } else {SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented");} 176 if (!isStride) { 177 ierr = ISRestoreIndices(pointIS,&points);CHKERRQ(ierr); 178 } 179 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 180 PetscFunctionReturn(0); 181 } 182 183 static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 184 { 185 DMField_DS *dsfield = (DMField_DS *) field->data; 186 PetscSF cellSF = NULL; 187 const PetscSFNode *cells; 188 PetscInt c, nFound, numCells, feDim, nc; 189 const PetscInt *cellDegrees; 190 const PetscScalar *pointsArray; 191 PetscScalar *cellPoints; 192 PetscInt gatherSize, gatherMax; 193 PetscInt dim, dimR, offset; 194 MPI_Datatype pointType; 195 PetscObject cellDisc; 196 PetscFE cellFE; 197 PetscClassId discID; 198 PetscReal *coordsReal, *coordsRef; 199 PetscSection section; 200 PetscScalar *cellBs = NULL, *cellDs = NULL, *cellHs = NULL; 201 PetscReal *cellBr = NULL, *cellDr = NULL, *cellHr = NULL; 202 PetscReal *v, *J, *invJ, *detJ; 203 PetscErrorCode ierr; 204 205 PetscFunctionBegin; 206 nc = field->numComponents; 207 ierr = DMGetDefaultSection(field->dm,§ion);CHKERRQ(ierr); 208 ierr = DMFieldDSGetHeightDisc(field,0,&cellDisc);CHKERRQ(ierr); 209 ierr = PetscObjectGetClassId(cellDisc, &discID);CHKERRQ(ierr); 210 if (discID != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Discretization type not supported\n"); 211 cellFE = (PetscFE) cellDisc; 212 ierr = PetscFEGetDimension(cellFE,&feDim);CHKERRQ(ierr); 213 ierr = DMGetCoordinateDim(field->dm, &dim);CHKERRQ(ierr); 214 ierr = DMGetDimension(field->dm, &dimR);CHKERRQ(ierr); 215 ierr = DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF);CHKERRQ(ierr); 216 ierr = PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells);CHKERRQ(ierr); 217 for (c = 0; c < nFound; c++) { 218 if (cells[c].index < 0) SETERRQ1(PetscObjectComm((PetscObject)points),PETSC_ERR_ARG_WRONG, "Point %D could not be located\n", c); 219 } 220 ierr = PetscSFComputeDegreeBegin(cellSF,&cellDegrees);CHKERRQ(ierr); 221 ierr = PetscSFComputeDegreeEnd(cellSF,&cellDegrees);CHKERRQ(ierr); 222 for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) { 223 gatherMax = PetscMax(gatherMax,cellDegrees[c]); 224 gatherSize += cellDegrees[c]; 225 } 226 ierr = PetscMalloc3(gatherSize*dim,&cellPoints,gatherMax*dim,&coordsReal,gatherMax*dimR,&coordsRef);CHKERRQ(ierr); 227 ierr = PetscMalloc4(gatherMax*dimR,&v,gatherMax*dimR*dimR,&J,gatherMax*dimR*dimR,&invJ,gatherMax,&detJ);CHKERRQ(ierr); 228 if (datatype == PETSC_SCALAR) { 229 ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBs, D ? nc * dim * gatherSize : 0, &cellDs, H ? nc * dim * dim * gatherSize : 0, &cellHs);CHKERRQ(ierr); 230 } else { 231 ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr);CHKERRQ(ierr); 232 } 233 234 ierr = MPI_Type_contiguous(dim,MPIU_SCALAR,&pointType);CHKERRQ(ierr); 235 ierr = MPI_Type_commit(&pointType);CHKERRQ(ierr); 236 ierr = VecGetArrayRead(points,&pointsArray);CHKERRQ(ierr); 237 ierr = PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr); 238 ierr = PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr); 239 ierr = VecRestoreArrayRead(points,&pointsArray);CHKERRQ(ierr); 240 for (c = 0, offset = 0; c < numCells; c++) { 241 PetscInt nq = cellDegrees[c], p; 242 243 if (nq) { 244 PetscReal *fB, *fD, *fH; 245 PetscInt closureSize; 246 PetscScalar *elem = NULL; 247 PetscReal *quadPoints; 248 PetscQuadrature quad; 249 PetscInt d, e, f, g; 250 251 for (p = 0; p < dim * nq; p++) coordsReal[p] = cellPoints[dim * offset + p]; 252 ierr = DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef);CHKERRQ(ierr); 253 ierr = PetscFEGetTabulation(cellFE,nq,coordsRef,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 254 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &quad);CHKERRQ(ierr); 255 ierr = PetscMalloc1(dimR * nq, &quadPoints);CHKERRQ(ierr); 256 for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p]; 257 ierr = PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL);CHKERRQ(ierr); 258 ierr = DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ);CHKERRQ(ierr); 259 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 260 ierr = DMPlexVecGetClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 261 if (B) { 262 if (datatype == PETSC_SCALAR) { 263 PetscScalar *cB = &cellBs[nc * offset]; 264 265 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,(PetscScalar)); 266 } else { 267 PetscReal *cB = &cellBr[nc * offset]; 268 269 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,PetscRealPart); 270 } 271 } 272 if (D) { 273 if (datatype == PETSC_SCALAR) { 274 PetscScalar *cD = &cellDs[nc * dim * offset]; 275 276 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),(PetscScalar)); 277 for (p = 0; p < nq; p++) { 278 for (g = 0; g < nc; g++) { 279 for (d = 0; d < dimR; d++) { 280 v[d] = 0.; 281 for (e = 0; e < dimR; e++) { 282 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 283 } 284 } 285 for (d = 0; d < dimR; d++) { 286 cD[(nc * p + g) * dimR + d] = v[d]; 287 } 288 } 289 } 290 } else { 291 PetscReal *cD = &cellDr[nc * dim * offset]; 292 293 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),PetscRealPart); 294 for (p = 0; p < nq; p++) { 295 for (g = 0; g < nc; g++) { 296 for (d = 0; d < dimR; d++) { 297 v[d] = 0.; 298 for (e = 0; e < dimR; e++) { 299 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 300 } 301 } 302 for (d = 0; d < dimR; d++) { 303 cD[(nc * p + g) * dimR + d] = v[d]; 304 } 305 } 306 } 307 } 308 } 309 if (H) { 310 if (datatype == PETSC_SCALAR) { 311 PetscScalar *cH = &cellHs[nc * dim * dim * offset]; 312 313 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),(PetscScalar)); 314 for (p = 0; p < nq; p++) { 315 for (g = 0; g < nc * dimR; g++) { 316 for (d = 0; d < dimR; d++) { 317 v[d] = 0.; 318 for (e = 0; e < dimR; e++) { 319 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e]; 320 } 321 } 322 for (d = 0; d < dimR; d++) { 323 cH[(nc * dimR * p + g) * dimR + d] = v[d]; 324 } 325 } 326 for (g = 0; g < nc; g++) { 327 for (f = 0; f < dimR; f++) { 328 for (d = 0; d < dimR; d++) { 329 v[d] = 0.; 330 for (e = 0; e < dimR; e++) { 331 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 332 } 333 } 334 for (d = 0; d < dimR; d++) { 335 cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 336 } 337 } 338 } 339 } 340 } else { 341 PetscReal *cH = &cellHr[nc * dim * dim * offset]; 342 343 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),PetscRealPart); 344 for (p = 0; p < nq; p++) { 345 for (g = 0; g < nc * dimR; g++) { 346 for (d = 0; d < dimR; d++) { 347 v[d] = 0.; 348 for (e = 0; e < dimR; e++) { 349 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e]; 350 } 351 } 352 for (d = 0; d < dimR; d++) { 353 cH[(nc * dimR * p + g) * dimR + d] = v[d]; 354 } 355 } 356 for (g = 0; g < nc; g++) { 357 for (f = 0; f < dimR; f++) { 358 for (d = 0; d < dimR; d++) { 359 v[d] = 0.; 360 for (e = 0; e < dimR; e++) { 361 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 362 } 363 } 364 for (d = 0; d < dimR; d++) { 365 cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 366 } 367 } 368 } 369 } 370 } 371 } 372 ierr = DMPlexVecRestoreClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 373 ierr = PetscFERestoreTabulation(cellFE,nq,coordsRef,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 374 } 375 offset += nq; 376 } 377 { 378 MPI_Datatype origtype; 379 380 if (datatype == PETSC_SCALAR) { 381 origtype = MPIU_SCALAR; 382 } else { 383 origtype = MPIU_REAL; 384 } 385 if (B) { 386 MPI_Datatype Btype; 387 388 ierr = MPI_Type_contiguous(nc, origtype, &Btype);CHKERRQ(ierr); 389 ierr = MPI_Type_commit(&Btype);CHKERRQ(ierr); 390 ierr = PetscSFScatterBegin(cellSF,Btype,(datatype == PETSC_SCALAR) ? cellBs : cellBr, B);CHKERRQ(ierr); 391 ierr = PetscSFScatterEnd(cellSF,Btype,(datatype == PETSC_SCALAR) ? cellBs : cellBr, B);CHKERRQ(ierr); 392 ierr = MPI_Type_free(&Btype);CHKERRQ(ierr); 393 } 394 if (D) { 395 MPI_Datatype Dtype; 396 397 ierr = MPI_Type_contiguous(nc * dim, origtype, &Dtype);CHKERRQ(ierr); 398 ierr = MPI_Type_commit(&Dtype);CHKERRQ(ierr); 399 ierr = PetscSFScatterBegin(cellSF,Dtype,(datatype == PETSC_SCALAR) ? cellDs : cellDr, D);CHKERRQ(ierr); 400 ierr = PetscSFScatterEnd(cellSF,Dtype,(datatype == PETSC_SCALAR) ? cellDs : cellDr, D);CHKERRQ(ierr); 401 ierr = MPI_Type_free(&Dtype);CHKERRQ(ierr); 402 } 403 if (H) { 404 MPI_Datatype Htype; 405 406 ierr = MPI_Type_contiguous(nc * dim * dim, origtype, &Htype);CHKERRQ(ierr); 407 ierr = MPI_Type_commit(&Htype);CHKERRQ(ierr); 408 ierr = PetscSFScatterBegin(cellSF,Htype,(datatype == PETSC_SCALAR) ? cellHs : cellHr, H);CHKERRQ(ierr); 409 ierr = PetscSFScatterEnd(cellSF,Htype,(datatype == PETSC_SCALAR) ? cellHs : cellHr, H);CHKERRQ(ierr); 410 ierr = MPI_Type_free(&Htype);CHKERRQ(ierr); 411 } 412 } 413 ierr = PetscFree4(v,J,invJ,detJ);CHKERRQ(ierr); 414 ierr = PetscFree3(cellBr, cellDr, cellHr);CHKERRQ(ierr); 415 ierr = PetscFree3(cellBs, cellDs, cellHs);CHKERRQ(ierr); 416 ierr = PetscFree3(cellPoints,coordsReal,coordsRef);CHKERRQ(ierr); 417 ierr = MPI_Type_free(&pointType);CHKERRQ(ierr); 418 ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421 422 static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H) 423 { 424 DMField_DS *dsfield = (DMField_DS *) field->data; 425 PetscInt h, imin; 426 PetscInt dim; 427 PetscClassId id; 428 PetscQuadrature quad = NULL; 429 PetscBool isAffine; 430 PetscFEGeom *geom; 431 PetscInt Nq, Nc, dimC, qNc, N; 432 PetscInt numPoints; 433 void *qB = NULL, *qD = NULL, *qH = NULL; 434 const PetscReal *weights; 435 MPI_Datatype mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL; 436 PetscObject disc; 437 DMField coordField; 438 PetscErrorCode ierr; 439 440 PetscFunctionBegin; 441 Nc = field->numComponents; 442 dsfield = (DMField_DS *) field->data; 443 ierr = DMGetCoordinateDim(field->dm, &dimC);CHKERRQ(ierr); 444 ierr = DMGetDimension(field->dm, &dim);CHKERRQ(ierr); 445 ierr = ISGetLocalSize(pointIS, &numPoints); 446 ierr = ISGetMinMax(pointIS,&imin,NULL);CHKERRQ(ierr); 447 for (h = 0; h < dsfield->height; h++) { 448 PetscInt hEnd; 449 450 ierr = DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd);CHKERRQ(ierr); 451 if (imin < hEnd) break; 452 } 453 dim -= h; 454 ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 455 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 456 if (id != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported\n"); 457 ierr = DMGetCoordinateField(field->dm, &coordField);CHKERRQ(ierr); 458 ierr = DMFieldGetFEInvariance(coordField, pointIS, NULL, &isAffine, NULL);CHKERRQ(ierr); 459 if (isAffine) { 460 ierr = DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad);CHKERRQ(ierr); 461 } 462 if (!quad) {ierr = DMFieldCreateDefaultQuadrature(field, pointIS, &quad);CHKERRQ(ierr);} 463 if (!quad) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages\n"); 464 ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&geom);CHKERRQ(ierr); 465 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights);CHKERRQ(ierr); 466 if (qNc != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components\n"); 467 N = numPoints * Nq * Nc; 468 if (B) ierr = DMGetWorkArray(field->dm, N, mpitype, &qB);CHKERRQ(ierr); 469 if (D) ierr = DMGetWorkArray(field->dm, N * dimC, mpitype, &qD);CHKERRQ(ierr); 470 if (H) ierr = DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH);CHKERRQ(ierr); 471 ierr = DMFieldEvaluateFE(field,pointIS,quad,type,qB,qD,qH);CHKERRQ(ierr); 472 if (B) { 473 PetscInt i, j, k; 474 475 if (type == PETSC_SCALAR) { 476 PetscScalar * sB = (PetscScalar *) B; 477 PetscScalar * sqB = (PetscScalar *) qB; 478 479 for (i = 0; i < numPoints; i++) { 480 PetscReal vol = 0.; 481 482 for (j = 0; j < Nc; j++) {sB[i * Nc + j] = 0.;} 483 for (k = 0; k < Nq; k++) { 484 vol += geom->detJ[i * Nq + k] * weights[k]; 485 for (j = 0; j < Nc; j++) { 486 sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[ (i * Nq + k) * Nc + j]; 487 } 488 } 489 for (k = 0; k < Nq * Nc; k++) sB[i * Nq * Nc + k] /= vol; 490 } 491 } else { 492 PetscReal * rB = (PetscReal *) B; 493 PetscReal * rqB = (PetscReal *) qB; 494 495 for (i = 0; i < numPoints; i++) { 496 PetscReal vol = 0.; 497 498 for (j = 0; j < Nc; j++) {rB[i * Nc + j] = 0.;} 499 for (k = 0; k < Nq; k++) { 500 vol += geom->detJ[i * Nq + k] * weights[k]; 501 for (j = 0; j < Nc; j++) { 502 rB[i * Nc + j] += weights[k] * rqB[ (i * Nq + k) * Nc + j]; 503 } 504 } 505 for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol; 506 } 507 } 508 } 509 if (D) { 510 PetscInt i, j, k, l, m; 511 512 if (type == PETSC_SCALAR) { 513 PetscScalar * sD = (PetscScalar *) D; 514 PetscScalar * sqD = (PetscScalar *) qD; 515 516 for (i = 0; i < numPoints; i++) { 517 PetscReal vol = 0.; 518 519 for (j = 0; j < Nc * dimC; j++) {sD[i * Nc * dimC + j] = 0.;} 520 for (k = 0; k < Nq; k++) { 521 vol += geom->detJ[i * Nq + k] * weights[k]; 522 for (j = 0; j < Nc; j++) { 523 PetscScalar pD[3] = {0.,0.,0.}; 524 525 for (l = 0; l < dimC; l++) { 526 for (m = 0; m < dim; m++) { 527 pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * sqD[((i * Nq + k) * Nc + j) * dim + m]; 528 } 529 } 530 for (l = 0; l < dimC; l++) { 531 sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l]; 532 } 533 } 534 } 535 for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol; 536 } 537 } else { 538 PetscReal * rD = (PetscReal *) D; 539 PetscReal * rqD = (PetscReal *) qD; 540 541 for (i = 0; i < numPoints; i++) { 542 PetscReal vol = 0.; 543 544 for (j = 0; j < Nc * dimC; j++) {rD[i * Nc * dimC + j] = 0.;} 545 for (k = 0; k < Nq; k++) { 546 vol += geom->detJ[i * Nq + k] * weights[k]; 547 for (j = 0; j < Nc; j++) { 548 PetscReal pD[3] = {0.,0.,0.}; 549 550 for (l = 0; l < dimC; l++) { 551 for (m = 0; m < dim; m++) { 552 pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * rqD[((i * Nq + k) * Nc + j) * dim + m]; 553 } 554 } 555 for (l = 0; l < dimC; l++) { 556 rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l]; 557 } 558 } 559 } 560 for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol; 561 } 562 } 563 } 564 if (H) { 565 PetscInt i, j, k, l, m, q, r; 566 567 if (type == PETSC_SCALAR) { 568 PetscScalar * sH = (PetscScalar *) H; 569 PetscScalar * sqH = (PetscScalar *) qH; 570 571 for (i = 0; i < numPoints; i++) { 572 PetscReal vol = 0.; 573 574 for (j = 0; j < Nc * dimC * dimC; j++) {sH[i * Nc * dimC * dimC + j] = 0.;} 575 for (k = 0; k < Nq; k++) { 576 const PetscScalar *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC]; 577 578 vol += geom->detJ[i * Nq + k] * weights[k]; 579 for (j = 0; j < Nc; j++) { 580 PetscScalar pH[3][3] = {{0.,0.,0.},{0.,0.,0.},{0.,0.,0.}}; 581 const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC]; 582 583 for (l = 0; l < dimC; l++) { 584 for (m = 0; m < dimC; m++) { 585 for (q = 0; q < dim; q++) { 586 for (r = 0; r < dim; r++) { 587 pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r]; 588 } 589 } 590 } 591 } 592 for (l = 0; l < dimC; l++) { 593 for (m = 0; m < dimC; m++) { 594 sH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m]; 595 } 596 } 597 } 598 } 599 for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol; 600 } 601 } else { 602 PetscReal * rH = (PetscReal *) H; 603 PetscReal * rqH = (PetscReal *) qH; 604 605 for (i = 0; i < numPoints; i++) { 606 PetscReal vol = 0.; 607 608 for (j = 0; j < Nc * dimC * dimC; j++) {rH[i * Nc * dimC * dimC + j] = 0.;} 609 for (k = 0; k < Nq; k++) { 610 const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC]; 611 612 vol += geom->detJ[i * Nq + k] * weights[k]; 613 for (j = 0; j < Nc; j++) { 614 PetscReal pH[3][3] = {{0.,0.,0.},{0.,0.,0.},{0.,0.,0.}}; 615 const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC]; 616 617 for (l = 0; l < dimC; l++) { 618 for (m = 0; m < dimC; m++) { 619 for (q = 0; q < dim; q++) { 620 for (r = 0; r < dim; r++) { 621 pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r]; 622 } 623 } 624 } 625 } 626 for (l = 0; l < dimC; l++) { 627 for (m = 0; m < dimC; m++) { 628 rH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m]; 629 } 630 } 631 } 632 } 633 for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol; 634 } 635 } 636 } 637 if (B) ierr = DMRestoreWorkArray(field->dm, N, mpitype, &qB);CHKERRQ(ierr); 638 if (D) ierr = DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD);CHKERRQ(ierr); 639 if (H) ierr = DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH);CHKERRQ(ierr); 640 ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr); 641 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 642 PetscFunctionReturn(0); 643 } 644 645 static PetscErrorCode DMFieldGetFEInvariance_DS(DMField field, IS pointIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic) 646 { 647 DMField_DS *dsfield; 648 PetscObject disc; 649 PetscInt h, imin, imax; 650 PetscClassId id; 651 PetscErrorCode ierr; 652 653 PetscFunctionBegin; 654 dsfield = (DMField_DS *) field->data; 655 ierr = ISGetMinMax(pointIS,&imin,&imax);CHKERRQ(ierr); 656 if (imin >= imax) { 657 if (isConstant) {*isConstant = PETSC_TRUE;} 658 if (isAffine) {*isAffine = PETSC_TRUE;} 659 if (isQuadratic) {*isQuadratic = PETSC_TRUE;} 660 PetscFunctionReturn(0); 661 } 662 for (h = 0; h < dsfield->height; h++) { 663 PetscInt hEnd; 664 665 ierr = DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd);CHKERRQ(ierr); 666 if (imin < hEnd) break; 667 } 668 ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 669 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 670 if (id == PETSCFE_CLASSID) { 671 PetscFE fe = (PetscFE) disc; 672 PetscInt order, maxOrder; 673 PetscBool tensor = PETSC_FALSE; 674 PetscSpace sp; 675 676 ierr = PetscFEGetBasisSpace(fe, &sp);CHKERRQ(ierr); 677 ierr = PetscSpaceGetOrder(sp,&order);CHKERRQ(ierr); 678 ierr = PetscSpacePolynomialGetTensor(sp,&tensor);CHKERRQ(ierr); 679 if (tensor) { 680 PetscInt dim; 681 682 ierr = DMGetDimension(field->dm,&dim);CHKERRQ(ierr); 683 maxOrder = order * dim; 684 } else { 685 maxOrder = order; 686 } 687 if (isConstant) *isConstant = (maxOrder < 1) ? PETSC_TRUE : PETSC_FALSE; 688 if (isAffine) *isAffine = (maxOrder < 2) ? PETSC_TRUE : PETSC_FALSE; 689 if (isQuadratic) *isQuadratic = (maxOrder < 3) ? PETSC_TRUE : PETSC_FALSE; 690 } 691 PetscFunctionReturn(0); 692 } 693 694 static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad) 695 { 696 PetscInt h, dim, imax, imin; 697 DM dm; 698 DMField_DS *dsfield; 699 PetscObject disc; 700 PetscFE fe; 701 PetscClassId id; 702 PetscErrorCode ierr; 703 704 705 PetscFunctionBegin; 706 dm = field->dm; 707 dsfield = (DMField_DS *) field->data; 708 ierr = ISGetMinMax(pointIS,&imax,&imin);CHKERRQ(ierr); 709 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 710 for (h = 0; h <= dim; h++) { 711 PetscInt hStart, hEnd; 712 713 ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 714 if (imin >= hStart && imax < hEnd) break; 715 } 716 *quad = NULL; 717 if (h < dsfield->height) { 718 ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 719 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 720 if (id != PETSCFE_CLASSID) PetscFunctionReturn(0); 721 fe = (PetscFE) disc; 722 ierr = PetscFEGetQuadrature(fe,quad);CHKERRQ(ierr); 723 ierr = PetscObjectReference((PetscObject)*quad);CHKERRQ(ierr); 724 } 725 PetscFunctionReturn(0); 726 } 727 728 static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom) 729 { 730 const PetscInt *points; 731 PetscInt p, dim, dE, numFaces, Nq; 732 PetscBool affineCells; 733 DMLabel depthLabel; 734 IS cellIS; 735 DM dm = field->dm; 736 PetscErrorCode ierr; 737 738 PetscFunctionBegin; 739 dim = geom->dim; 740 dE = geom->dimEmbed; 741 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 742 ierr = DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS);CHKERRQ(ierr); 743 ierr = DMFieldGetFEInvariance(field,cellIS,NULL,&affineCells,NULL);CHKERRQ(ierr); 744 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 745 numFaces = geom->numCells; 746 Nq = geom->numPoints; 747 if (affineCells) { 748 PetscInt numCells, offset, *cells; 749 PetscFEGeom *cellGeom; 750 IS suppIS; 751 PetscQuadrature cellQuad = NULL; 752 753 ierr = DMFieldCreateDefaultQuadrature(field,cellIS,&cellQuad);CHKERRQ(ierr); 754 for (p = 0, numCells = 0; p < numFaces; p++) { 755 PetscInt point = points[p]; 756 PetscInt numSupp, numChildren; 757 758 ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL); CHKERRQ(ierr); 759 if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 760 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 761 if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp); 762 numCells += numSupp; 763 } 764 ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr); 765 for (p = 0, offset = 0; p < numFaces; p++) { 766 PetscInt point = points[p]; 767 PetscInt numSupp, s; 768 const PetscInt *supp; 769 770 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 771 ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 772 for (s = 0; s < numSupp; s++, offset++) { 773 cells[offset] = supp[s]; 774 } 775 } 776 ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr); 777 ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr); 778 for (p = 0, offset = 0; p < numFaces; p++) { 779 PetscInt point = points[p]; 780 PetscInt numSupp, s, q; 781 const PetscInt *supp; 782 783 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 784 ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 785 for (s = 0; s < numSupp; s++, offset++) { 786 for (q = 0; q < Nq * dE * dE; q++) { 787 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 788 } 789 } 790 } 791 ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr); 792 ierr = ISDestroy(&suppIS);CHKERRQ(ierr); 793 ierr = PetscFree(cells);CHKERRQ(ierr); 794 ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr); 795 } else { 796 PetscObject faceDisc, cellDisc; 797 PetscClassId faceId, cellId; 798 PetscDualSpace dsp; 799 DM K; 800 PetscInt (*co)[2][3]; 801 PetscInt coneSize; 802 PetscInt **counts; 803 PetscInt f, i, o, q, s; 804 const PetscInt *coneK; 805 PetscInt minOrient, maxOrient, numOrient; 806 PetscInt *orients; 807 PetscReal **orientPoints; 808 PetscReal *cellPoints; 809 PetscReal *dummyWeights; 810 PetscQuadrature cellQuad = NULL; 811 812 ierr = DMFieldDSGetHeightDisc(field, 1, &faceDisc);CHKERRQ(ierr); 813 ierr = DMFieldDSGetHeightDisc(field, 0, &cellDisc);CHKERRQ(ierr); 814 ierr = PetscObjectGetClassId(faceDisc,&faceId);CHKERRQ(ierr); 815 ierr = PetscObjectGetClassId(cellDisc,&cellId);CHKERRQ(ierr); 816 if (faceId != PETSCFE_CLASSID || cellId != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported\n"); 817 ierr = PetscFEGetDualSpace((PetscFE)cellDisc, &dsp);CHKERRQ(ierr); 818 ierr = PetscDualSpaceGetDM(dsp, &K); CHKERRQ(ierr); 819 ierr = DMPlexGetConeSize(K,0,&coneSize);CHKERRQ(ierr); 820 ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr); 821 ierr = PetscMalloc4(numFaces,&co,dE*Nq,&cellPoints,coneSize,&counts,Nq,&dummyWeights);CHKERRQ(ierr); 822 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)field), &cellQuad);CHKERRQ(ierr); 823 ierr = PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights);CHKERRQ(ierr); 824 minOrient = PETSC_MAX_INT; 825 maxOrient = PETSC_MIN_INT; 826 for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */ 827 PetscInt point = points[p]; 828 PetscInt numSupp, numChildren; 829 const PetscInt *supp; 830 831 ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL); CHKERRQ(ierr); 832 if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 833 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 834 if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp); 835 ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 836 for (s = 0; s < numSupp; s++) { 837 PetscInt cell = supp[s]; 838 PetscInt numCone; 839 const PetscInt *cone, *orient; 840 841 ierr = DMPlexGetConeSize(dm, cell, &numCone);CHKERRQ(ierr); 842 if (numCone != coneSize) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element"); 843 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 844 ierr = DMPlexGetConeOrientation(dm, cell, &orient);CHKERRQ(ierr); 845 for (f = 0; f < coneSize; f++) { 846 if (cone[f] == point) break; 847 } 848 co[p][s][0] = f; 849 co[p][s][1] = orient[f]; 850 co[p][s][2] = cell; 851 minOrient = PetscMin(minOrient, orient[f]); 852 maxOrient = PetscMin(maxOrient, orient[f]); 853 } 854 for (; s < 2; s++) { 855 co[p][s][0] = -1; 856 co[p][s][1] = -1; 857 co[p][s][2] = -1; 858 } 859 } 860 numOrient = maxOrient + 1 - minOrient; 861 ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr); 862 /* count all (face,orientation) doubles that appear */ 863 ierr = PetscCalloc2(numOrient,&orients,numOrient,&orientPoints);CHKERRQ(ierr); 864 for (f = 0; f < coneSize; f++) {ierr = PetscCalloc1(numOrient, &counts[f]);CHKERRQ(ierr);} 865 for (p = 0; p < numFaces; p++) { 866 for (s = 0; s < 2; s++) { 867 if (co[p][s][0] >= 0) { 868 counts[co[p][s][0]][co[p][s][1] - minOrient]++; 869 orients[co[p][s][1] - minOrient]++; 870 } 871 } 872 } 873 for (o = 0; o < numOrient; o++) { 874 if (orients[o]) { 875 PetscInt orient = o + minOrient; 876 PetscInt q; 877 878 ierr = PetscMalloc1(Nq * dim, &orientPoints[o]);CHKERRQ(ierr); 879 /* rotate the quadrature points appropriately */ 880 switch (dim) { 881 case 0: 882 break; 883 case 1: 884 if (orient == -2 || orient == 1) { 885 for (q = 0; q < Nq; q++) { 886 orientPoints[o][q] = -geom->xi[q]; 887 } 888 } else { 889 for (q = 0; q < Nq; q++) { 890 orientPoints[o][q] = geom->xi[q]; 891 } 892 } 893 break; 894 case 2: 895 switch (coneSize) { 896 case 3: 897 for (q = 0; q < Nq; q++) { 898 PetscReal lambda[3]; 899 PetscReal lambdao[3]; 900 901 /* convert to barycentric */ 902 lambda[0] = - (geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.; 903 lambda[1] = (geom->xi[2 * q] + 1.) / 2.; 904 lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.; 905 if (orient >= 0) { 906 for (i = 0; i < 3; i++) { 907 lambdao[i] = lambda[(orient + i) % 3]; 908 } 909 } else { 910 for (i = 0; i < 3; i++) { 911 lambdao[i] = lambda[(-(orient + i) + 3) % 3]; 912 } 913 } 914 /* convert to coordinates */ 915 orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1]; 916 orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2]; 917 } 918 break; 919 case 4: 920 for (q = 0; q < Nq; q++) { 921 PetscReal xi[2], xio[2]; 922 PetscInt oabs = (orient >= 0) ? orient : -(orient + 1); 923 924 xi[0] = geom->xi[2 * q]; 925 xi[1] = geom->xi[2 * q + 1]; 926 switch (oabs) { 927 case 0: 928 xio[0] = xi[0]; 929 xio[1] = xi[1]; 930 break; 931 case 1: 932 xio[0] = xi[1]; 933 xio[1] = -xi[0]; 934 break; 935 case 2: 936 xio[0] = -xi[0]; 937 xio[1] = -xi[1]; 938 case 3: 939 xio[0] = -xi[1]; 940 xio[1] = xi[0]; 941 } 942 if (orient < 0) { 943 xio[0] = -xio[0]; 944 } 945 orientPoints[o][2 * q + 0] = xio[0]; 946 orientPoints[o][2 * q + 1] = xio[1]; 947 } 948 break; 949 default: 950 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not implemented yet\n"); 951 } 952 default: 953 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not implemented yet\n"); 954 } 955 } 956 } 957 for (f = 0; f < coneSize; f++) { 958 PetscInt face = coneK[f]; 959 PetscScalar v0[3]; 960 PetscScalar J[9]; 961 PetscInt numCells, offset; 962 PetscInt *cells; 963 IS suppIS; 964 965 ierr = DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, NULL);CHKERRQ(ierr); 966 for (o = 0; o <= numOrient; o++) { 967 PetscFEGeom *cellGeom; 968 969 if (!counts[f][o]) continue; 970 /* If this (face,orientation) double appears, 971 * convert the face quadrature points into volume quadrature points */ 972 for (q = 0; q < Nq; q++) { 973 PetscReal xi0[3] = {-1., -1., -1.}; 974 975 CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]); 976 } 977 for (p = 0, numCells = 0; p < numFaces; p++) { 978 for (s = 0; s < 2; s++) { 979 if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++; 980 } 981 } 982 ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr); 983 for (p = 0, offset = 0; p < numFaces; p++) { 984 for (s = 0; s < 2; s++) { 985 if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 986 cells[offset++] = co[p][s][2]; 987 } 988 } 989 } 990 ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr); 991 ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr); 992 for (p = 0, offset = 0; p < numFaces; p++) { 993 for (s = 0; s < 2; s++) { 994 if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 995 for (q = 0; q < Nq * dE * dE; q++) { 996 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 997 } 998 offset++; 999 } 1000 } 1001 } 1002 ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr); 1003 ierr = ISDestroy(&suppIS);CHKERRQ(ierr); 1004 ierr = PetscFree(cells);CHKERRQ(ierr); 1005 } 1006 } 1007 for (o = 0; o < numOrient; o++) { 1008 if (orients[o]) { 1009 ierr = PetscFree(orientPoints[o]);CHKERRQ(ierr); 1010 } 1011 } 1012 ierr = PetscFree2(orients,orientPoints);CHKERRQ(ierr); 1013 ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr); 1014 ierr = PetscFree4(co,cellPoints,counts,dummyWeights);CHKERRQ(ierr); 1015 } 1016 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1017 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 1018 PetscFunctionReturn(0); 1019 } 1020 1021 static PetscErrorCode DMFieldInitialize_DS(DMField field) 1022 { 1023 PetscFunctionBegin; 1024 field->ops->destroy = DMFieldDestroy_DS; 1025 field->ops->evaluate = DMFieldEvaluate_DS; 1026 field->ops->evaluateFE = DMFieldEvaluateFE_DS; 1027 field->ops->evaluateFV = DMFieldEvaluateFV_DS; 1028 field->ops->getFEInvariance = DMFieldGetFEInvariance_DS; 1029 field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS; 1030 field->ops->view = DMFieldView_DS; 1031 field->ops->computeFaceData = DMFieldComputeFaceData_DS; 1032 PetscFunctionReturn(0); 1033 } 1034 1035 PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field) 1036 { 1037 DMField_DS *dsfield; 1038 PetscErrorCode ierr; 1039 1040 PetscFunctionBegin; 1041 ierr = PetscNewLog(field,&dsfield);CHKERRQ(ierr); 1042 field->data = dsfield; 1043 ierr = DMFieldInitialize_DS(field);CHKERRQ(ierr); 1044 PetscFunctionReturn(0); 1045 } 1046 1047 PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec,DMField *field) 1048 { 1049 DMField b; 1050 DMField_DS *dsfield; 1051 PetscObject disc = NULL; 1052 PetscBool isContainer = PETSC_FALSE; 1053 PetscClassId id = -1; 1054 PetscInt numComponents = -1, dsNumFields; 1055 PetscSection section; 1056 PetscErrorCode ierr; 1057 1058 PetscFunctionBegin; 1059 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 1060 ierr = PetscSectionGetFieldComponents(section,fieldNum,&numComponents);CHKERRQ(ierr); 1061 ierr = DMGetNumFields(dm,&dsNumFields);CHKERRQ(ierr); 1062 if (dsNumFields) {ierr = DMGetField(dm,fieldNum,&disc);CHKERRQ(ierr);} 1063 if (disc) { 1064 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 1065 isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE; 1066 } 1067 if (!disc || isContainer) { 1068 MPI_Comm comm = PetscObjectComm((PetscObject) dm); 1069 PetscInt cStart, cEnd, dim; 1070 PetscInt localConeSize = 0, coneSize; 1071 PetscFE fe; 1072 PetscDualSpace Q; 1073 PetscSpace P; 1074 DM K; 1075 PetscQuadrature quad, fquad; 1076 PetscBool isSimplex; 1077 1078 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1079 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1080 if (cEnd > cStart) { 1081 ierr = DMPlexGetConeSize(dm, cStart, &localConeSize);CHKERRQ(ierr); 1082 } 1083 ierr = MPI_Allreduce(&localConeSize,&coneSize,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 1084 isSimplex = (coneSize == (dim + 1)) ? PETSC_TRUE : PETSC_FALSE; 1085 ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 1086 ierr = PetscSpaceSetOrder(P, 1);CHKERRQ(ierr); 1087 ierr = PetscSpaceSetNumComponents(P, numComponents);CHKERRQ(ierr); 1088 ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 1089 ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 1090 ierr = PetscSpacePolynomialSetTensor(P, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr); 1091 ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 1092 ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 1093 ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 1094 ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 1095 ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 1096 ierr = DMDestroy(&K);CHKERRQ(ierr); 1097 ierr = PetscDualSpaceSetNumComponents(Q, numComponents);CHKERRQ(ierr); 1098 ierr = PetscDualSpaceSetOrder(Q, 1);CHKERRQ(ierr); 1099 ierr = PetscDualSpaceLagrangeSetTensor(Q, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr); 1100 ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 1101 ierr = PetscFECreate(comm, &fe);CHKERRQ(ierr); 1102 ierr = PetscFESetType(fe,PETSCFEBASIC);CHKERRQ(ierr); 1103 ierr = PetscFESetBasisSpace(fe, P);CHKERRQ(ierr); 1104 ierr = PetscFESetDualSpace(fe, Q);CHKERRQ(ierr); 1105 ierr = PetscFESetNumComponents(fe, numComponents);CHKERRQ(ierr); 1106 ierr = PetscFESetUp(fe);CHKERRQ(ierr); 1107 ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 1108 ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 1109 if (isSimplex) { 1110 ierr = PetscDTGaussJacobiQuadrature(dim, 1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr); 1111 ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr); 1112 } 1113 else { 1114 ierr = PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr); 1115 ierr = PetscDTGaussTensorQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr); 1116 } 1117 ierr = PetscFESetQuadrature(fe, quad);CHKERRQ(ierr); 1118 ierr = PetscFESetFaceQuadrature(fe, fquad);CHKERRQ(ierr); 1119 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 1120 ierr = PetscQuadratureDestroy(&fquad);CHKERRQ(ierr); 1121 disc = (PetscObject) fe; 1122 } else { 1123 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 1124 } 1125 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 1126 if (id == PETSCFE_CLASSID) { 1127 PetscFE fe = (PetscFE) disc; 1128 1129 ierr = PetscFEGetNumComponents(fe,&numComponents);CHKERRQ(ierr); 1130 } else {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented");} 1131 ierr = DMFieldCreate(dm,numComponents,DMFIELD_VERTEX,&b);CHKERRQ(ierr); 1132 ierr = DMFieldSetType(b,DMFIELDDS);CHKERRQ(ierr); 1133 dsfield = (DMField_DS *) b->data; 1134 dsfield->fieldNum = fieldNum; 1135 ierr = DMGetDimension(dm,&dsfield->height);CHKERRQ(ierr); 1136 dsfield->height++; 1137 ierr = PetscCalloc1(dsfield->height,&dsfield->disc);CHKERRQ(ierr); 1138 dsfield->disc[0] = disc; 1139 ierr = PetscObjectReference((PetscObject)vec);CHKERRQ(ierr); 1140 dsfield->vec = vec; 1141 *field = b; 1142 1143 PetscFunctionReturn(0); 1144 } 1145