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