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 if (dim == meshDim - 1) { 132 /* TODO */ 133 } 134 ierr = PetscFEGetDimension(fe,&feDim);CHKERRQ(ierr); 135 ierr = PetscFEGetTabulation(fe,nq,qpoints,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 136 for (i = 0; i < numCells; i++) { 137 PetscInt c = isStride ? (sfirst + i * stride) : points[i]; 138 PetscInt closureSize; 139 PetscScalar *elem = NULL; 140 141 ierr = DMPlexVecGetClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 142 if (B) { 143 if (type == PETSC_SCALAR) { 144 PetscScalar *cB = &((PetscScalar *) B)[nc * nq * i]; 145 146 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,(PetscScalar)); 147 } else { 148 PetscReal *cB = &((PetscReal *) B)[nc * nq * i]; 149 150 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,PetscRealPart); 151 } 152 } 153 if (D) { 154 if (type == PETSC_SCALAR) { 155 PetscScalar *cD = &((PetscScalar *) D)[nc * nq * dim * i]; 156 157 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),(PetscScalar)); 158 } else { 159 PetscReal *cD = &((PetscReal *) D)[nc * nq * dim * i]; 160 161 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),PetscRealPart); 162 } 163 } 164 if (H) { 165 if (type == PETSC_SCALAR) { 166 PetscScalar *cH = &((PetscScalar *) H)[nc * nq * dim * dim * i]; 167 168 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),(PetscScalar)); 169 } else { 170 PetscReal *cH = &((PetscReal *) H)[nc * nq * dim * dim * i]; 171 172 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),PetscRealPart); 173 } 174 } 175 ierr = DMPlexVecRestoreClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 176 } 177 ierr = PetscFERestoreTabulation(fe,nq,qpoints,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 178 } else {SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented");} 179 if (!isStride) { 180 ierr = ISRestoreIndices(pointIS,&points);CHKERRQ(ierr); 181 } 182 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 183 PetscFunctionReturn(0); 184 } 185 186 static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 187 { 188 DMField_DS *dsfield = (DMField_DS *) field->data; 189 PetscSF cellSF = NULL; 190 const PetscSFNode *cells; 191 PetscInt c, nFound, numCells, feDim, nc; 192 const PetscInt *cellDegrees; 193 const PetscScalar *pointsArray; 194 PetscScalar *cellPoints; 195 PetscInt gatherSize, gatherMax; 196 PetscInt dim, dimR, offset; 197 MPI_Datatype pointType; 198 PetscObject cellDisc; 199 PetscFE cellFE; 200 PetscClassId discID; 201 PetscReal *coordsReal, *coordsRef; 202 PetscSection section; 203 PetscScalar *cellBs = NULL, *cellDs = NULL, *cellHs = NULL; 204 PetscReal *cellBr = NULL, *cellDr = NULL, *cellHr = NULL; 205 PetscReal *v, *J, *invJ, *detJ; 206 PetscErrorCode ierr; 207 208 PetscFunctionBegin; 209 nc = field->numComponents; 210 ierr = DMGetDefaultSection(field->dm,§ion);CHKERRQ(ierr); 211 ierr = DMFieldDSGetHeightDisc(field,0,&cellDisc);CHKERRQ(ierr); 212 ierr = PetscObjectGetClassId(cellDisc, &discID);CHKERRQ(ierr); 213 if (discID != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Discretization type not supported\n"); 214 cellFE = (PetscFE) cellDisc; 215 ierr = PetscFEGetDimension(cellFE,&feDim);CHKERRQ(ierr); 216 ierr = DMGetCoordinateDim(field->dm, &dim);CHKERRQ(ierr); 217 ierr = DMGetDimension(field->dm, &dimR);CHKERRQ(ierr); 218 ierr = DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF);CHKERRQ(ierr); 219 ierr = PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells);CHKERRQ(ierr); 220 for (c = 0; c < nFound; c++) { 221 if (cells[c].index < 0) SETERRQ1(PetscObjectComm((PetscObject)points),PETSC_ERR_ARG_WRONG, "Point %D could not be located\n", c); 222 } 223 ierr = PetscSFComputeDegreeBegin(cellSF,&cellDegrees);CHKERRQ(ierr); 224 ierr = PetscSFComputeDegreeEnd(cellSF,&cellDegrees);CHKERRQ(ierr); 225 for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) { 226 gatherMax = PetscMax(gatherMax,cellDegrees[c]); 227 gatherSize += cellDegrees[c]; 228 } 229 ierr = PetscMalloc3(gatherSize*dim,&cellPoints,gatherMax*dim,&coordsReal,gatherMax*dimR,&coordsRef);CHKERRQ(ierr); 230 ierr = PetscMalloc4(gatherMax*dimR,&v,gatherMax*dimR*dimR,&J,gatherMax*dimR*dimR,&invJ,gatherMax,&detJ);CHKERRQ(ierr); 231 if (datatype == PETSC_SCALAR) { 232 ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBs, D ? nc * dim * gatherSize : 0, &cellDs, H ? nc * dim * dim * gatherSize : 0, &cellHs);CHKERRQ(ierr); 233 } else { 234 ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr);CHKERRQ(ierr); 235 } 236 237 ierr = MPI_Type_contiguous(dim,MPIU_SCALAR,&pointType);CHKERRQ(ierr); 238 ierr = MPI_Type_commit(&pointType);CHKERRQ(ierr); 239 ierr = VecGetArrayRead(points,&pointsArray);CHKERRQ(ierr); 240 ierr = PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr); 241 ierr = PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr); 242 ierr = VecRestoreArrayRead(points,&pointsArray);CHKERRQ(ierr); 243 for (c = 0, offset = 0; c < numCells; c++) { 244 PetscInt nq = cellDegrees[c], p; 245 246 if (nq) { 247 PetscReal *fB, *fD, *fH; 248 PetscInt closureSize; 249 PetscScalar *elem = NULL; 250 PetscReal *quadPoints; 251 PetscQuadrature quad; 252 PetscInt d, e, f, g; 253 254 for (p = 0; p < dim * nq; p++) coordsReal[p] = cellPoints[dim * offset + p]; 255 ierr = DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef);CHKERRQ(ierr); 256 ierr = PetscFEGetTabulation(cellFE,nq,coordsRef,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 257 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &quad);CHKERRQ(ierr); 258 ierr = PetscMalloc1(dimR * nq, &quadPoints);CHKERRQ(ierr); 259 for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p]; 260 ierr = PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL);CHKERRQ(ierr); 261 ierr = DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ);CHKERRQ(ierr); 262 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 263 ierr = DMPlexVecGetClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 264 if (B) { 265 if (datatype == PETSC_SCALAR) { 266 PetscScalar *cB = &cellBs[nc * offset]; 267 268 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,(PetscScalar)); 269 } else { 270 PetscReal *cB = &cellBr[nc * offset]; 271 272 DMFieldDSdot(cB,fB,elem,nq,feDim,nc,PetscRealPart); 273 } 274 } 275 if (D) { 276 if (datatype == PETSC_SCALAR) { 277 PetscScalar *cD = &cellDs[nc * dim * offset]; 278 279 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),(PetscScalar)); 280 for (p = 0; p < nq; p++) { 281 for (g = 0; g < nc; g++) { 282 for (d = 0; d < dimR; d++) { 283 v[d] = 0.; 284 for (e = 0; e < dimR; e++) { 285 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 286 } 287 } 288 for (d = 0; d < dimR; d++) { 289 cD[(nc * p + g) * dimR + d] = v[d]; 290 } 291 } 292 } 293 } else { 294 PetscReal *cD = &cellDr[nc * dim * offset]; 295 296 DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),PetscRealPart); 297 for (p = 0; p < nq; p++) { 298 for (g = 0; g < nc; g++) { 299 for (d = 0; d < dimR; d++) { 300 v[d] = 0.; 301 for (e = 0; e < dimR; e++) { 302 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 303 } 304 } 305 for (d = 0; d < dimR; d++) { 306 cD[(nc * p + g) * dimR + d] = v[d]; 307 } 308 } 309 } 310 } 311 } 312 if (H) { 313 if (datatype == PETSC_SCALAR) { 314 PetscScalar *cH = &cellHs[nc * dim * dim * offset]; 315 316 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),(PetscScalar)); 317 for (p = 0; p < nq; p++) { 318 for (g = 0; g < nc * dimR; g++) { 319 for (d = 0; d < dimR; d++) { 320 v[d] = 0.; 321 for (e = 0; e < dimR; e++) { 322 v[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] = v[d]; 327 } 328 } 329 for (g = 0; g < nc; g++) { 330 for (f = 0; f < dimR; f++) { 331 for (d = 0; d < dimR; d++) { 332 v[d] = 0.; 333 for (e = 0; e < dimR; e++) { 334 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 335 } 336 } 337 for (d = 0; d < dimR; d++) { 338 cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 339 } 340 } 341 } 342 } 343 } else { 344 PetscReal *cH = &cellHr[nc * dim * dim * offset]; 345 346 DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),PetscRealPart); 347 for (p = 0; p < nq; p++) { 348 for (g = 0; g < nc * dimR; g++) { 349 for (d = 0; d < dimR; d++) { 350 v[d] = 0.; 351 for (e = 0; e < dimR; e++) { 352 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e]; 353 } 354 } 355 for (d = 0; d < dimR; d++) { 356 cH[(nc * dimR * p + g) * dimR + d] = v[d]; 357 } 358 } 359 for (g = 0; g < nc; g++) { 360 for (f = 0; f < dimR; f++) { 361 for (d = 0; d < dimR; d++) { 362 v[d] = 0.; 363 for (e = 0; e < dimR; e++) { 364 v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 365 } 366 } 367 for (d = 0; d < dimR; d++) { 368 cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 369 } 370 } 371 } 372 } 373 } 374 } 375 ierr = DMPlexVecRestoreClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 376 ierr = PetscFERestoreTabulation(cellFE,nq,coordsRef,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 377 } 378 offset += nq; 379 } 380 { 381 MPI_Datatype origtype; 382 383 if (datatype == PETSC_SCALAR) { 384 origtype = MPIU_SCALAR; 385 } else { 386 origtype = MPIU_REAL; 387 } 388 if (B) { 389 MPI_Datatype Btype; 390 391 ierr = MPI_Type_contiguous(nc, origtype, &Btype);CHKERRQ(ierr); 392 ierr = MPI_Type_commit(&Btype);CHKERRQ(ierr); 393 ierr = PetscSFScatterBegin(cellSF,Btype,(datatype == PETSC_SCALAR) ? cellBs : cellBr, B);CHKERRQ(ierr); 394 ierr = PetscSFScatterEnd(cellSF,Btype,(datatype == PETSC_SCALAR) ? cellBs : cellBr, B);CHKERRQ(ierr); 395 ierr = MPI_Type_free(&Btype);CHKERRQ(ierr); 396 } 397 if (D) { 398 MPI_Datatype Dtype; 399 400 ierr = MPI_Type_contiguous(nc * dim, origtype, &Dtype);CHKERRQ(ierr); 401 ierr = MPI_Type_commit(&Dtype);CHKERRQ(ierr); 402 ierr = PetscSFScatterBegin(cellSF,Dtype,(datatype == PETSC_SCALAR) ? cellDs : cellDr, D);CHKERRQ(ierr); 403 ierr = PetscSFScatterEnd(cellSF,Dtype,(datatype == PETSC_SCALAR) ? cellDs : cellDr, D);CHKERRQ(ierr); 404 ierr = MPI_Type_free(&Dtype);CHKERRQ(ierr); 405 } 406 if (H) { 407 MPI_Datatype Htype; 408 409 ierr = MPI_Type_contiguous(nc * dim * dim, origtype, &Htype);CHKERRQ(ierr); 410 ierr = MPI_Type_commit(&Htype);CHKERRQ(ierr); 411 ierr = PetscSFScatterBegin(cellSF,Htype,(datatype == PETSC_SCALAR) ? cellHs : cellHr, H);CHKERRQ(ierr); 412 ierr = PetscSFScatterEnd(cellSF,Htype,(datatype == PETSC_SCALAR) ? cellHs : cellHr, H);CHKERRQ(ierr); 413 ierr = MPI_Type_free(&Htype);CHKERRQ(ierr); 414 } 415 } 416 ierr = PetscFree4(v,J,invJ,detJ);CHKERRQ(ierr); 417 ierr = PetscFree3(cellBr, cellDr, cellHr);CHKERRQ(ierr); 418 ierr = PetscFree3(cellBs, cellDs, cellHs);CHKERRQ(ierr); 419 ierr = PetscFree3(cellPoints,coordsReal,coordsRef);CHKERRQ(ierr); 420 ierr = MPI_Type_free(&pointType);CHKERRQ(ierr); 421 ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424 425 static PetscErrorCode DMFieldGetFEInvariance_DS(DMField field, IS pointIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic) 426 { 427 DMField_DS *dsfield; 428 PetscObject disc; 429 PetscInt h, imin; 430 PetscClassId id; 431 PetscErrorCode ierr; 432 433 PetscFunctionBegin; 434 dsfield = (DMField_DS *) field->data; 435 ierr = ISGetMinMax(pointIS,&imin,NULL);CHKERRQ(ierr); 436 for (h = 0; h < dsfield->height; h++) { 437 PetscInt hEnd; 438 439 ierr = DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd);CHKERRQ(ierr); 440 if (imin < hEnd) break; 441 } 442 ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 443 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 444 if (id == PETSCFE_CLASSID) { 445 PetscFE fe = (PetscFE) disc; 446 PetscInt order, maxOrder; 447 PetscBool tensor = PETSC_FALSE; 448 PetscSpace sp; 449 450 ierr = PetscFEGetBasisSpace(fe, &sp);CHKERRQ(ierr); 451 ierr = PetscSpaceGetOrder(sp,&order);CHKERRQ(ierr); 452 ierr = PetscSpacePolynomialGetTensor(sp,&tensor);CHKERRQ(ierr); 453 if (tensor) { 454 PetscInt dim; 455 456 ierr = DMGetDimension(field->dm,&dim);CHKERRQ(ierr); 457 maxOrder = order * dim; 458 } else { 459 maxOrder = order; 460 } 461 if (isConstant) *isConstant = (maxOrder < 1) ? PETSC_TRUE : PETSC_FALSE; 462 if (isAffine) *isAffine = (maxOrder < 2) ? PETSC_TRUE : PETSC_FALSE; 463 if (isQuadratic) *isQuadratic = (maxOrder < 3) ? PETSC_TRUE : PETSC_FALSE; 464 } 465 PetscFunctionReturn(0); 466 } 467 468 static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad) 469 { 470 PetscInt h, dim, imax, imin; 471 DM dm; 472 DMField_DS *dsfield; 473 PetscObject disc; 474 PetscFE fe; 475 PetscClassId id; 476 PetscErrorCode ierr; 477 478 479 PetscFunctionBegin; 480 dm = field->dm; 481 dsfield = (DMField_DS *) field->data; 482 ierr = ISGetMinMax(pointIS,&imax,&imin);CHKERRQ(ierr); 483 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 484 for (h = 0; h <= dim; h++) { 485 PetscInt hStart, hEnd; 486 487 ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 488 if (imin >= hStart && imax < hEnd) break; 489 } 490 *quad = NULL; 491 if (h < dsfield->height) { 492 ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 493 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 494 if (id != PETSCFE_CLASSID) PetscFunctionReturn(0); 495 fe = (PetscFE) disc; 496 ierr = PetscFEGetQuadrature(fe,quad);CHKERRQ(ierr); 497 ierr = PetscObjectReference((PetscObject)*quad);CHKERRQ(ierr); 498 } 499 PetscFunctionReturn(0); 500 } 501 502 static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom) 503 { 504 const PetscInt *points; 505 PetscInt p, dim, dE, numFaces, Nq; 506 PetscBool affineCells; 507 DMLabel depthLabel; 508 IS cellIS; 509 DM dm = field->dm; 510 PetscErrorCode ierr; 511 512 PetscFunctionBegin; 513 dim = geom->dim; 514 dE = geom->dimEmbed; 515 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 516 ierr = DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS);CHKERRQ(ierr); 517 ierr = DMFieldGetFEInvariance(field,cellIS,NULL,&affineCells,NULL);CHKERRQ(ierr); 518 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 519 numFaces = geom->numCells; 520 Nq = geom->numPoints; 521 if (affineCells) { 522 PetscInt numCells, offset, *cells; 523 PetscFEGeom *cellGeom; 524 IS suppIS; 525 PetscQuadrature cellQuad = NULL; 526 527 ierr = DMFieldCreateDefaultQuadrature(field,cellIS,&cellQuad);CHKERRQ(ierr); 528 for (p = 0, numCells = 0; p < numFaces; p++) { 529 PetscInt point = points[p]; 530 PetscInt numSupp, numChildren; 531 532 ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL); CHKERRQ(ierr); 533 if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 534 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 535 if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp); 536 numCells += numSupp; 537 } 538 ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr); 539 for (p = 0, offset = 0; p < numFaces; p++) { 540 PetscInt point = points[p]; 541 PetscInt numSupp, s; 542 const PetscInt *supp; 543 544 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 545 ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 546 for (s = 0; s < numSupp; s++, offset++) { 547 cells[offset] = supp[s]; 548 } 549 } 550 ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr); 551 ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr); 552 for (p = 0, offset = 0; p < numFaces; p++) { 553 PetscInt point = points[p]; 554 PetscInt numSupp, s, q; 555 const PetscInt *supp; 556 557 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 558 ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 559 for (s = 0; s < numSupp; s++, offset++) { 560 for (q = 0; q < Nq * dE * dE; q++) { 561 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 562 } 563 } 564 } 565 ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr); 566 ierr = ISDestroy(&suppIS);CHKERRQ(ierr); 567 ierr = PetscFree(cells);CHKERRQ(ierr); 568 ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr); 569 } else { 570 PetscObject faceDisc, cellDisc; 571 PetscClassId faceId, cellId; 572 PetscDualSpace dsp; 573 DM K; 574 PetscInt (*co)[2][3]; 575 PetscInt coneSize; 576 PetscInt **counts; 577 PetscInt f, i, o, q, s; 578 const PetscInt *coneK; 579 PetscInt minOrient, maxOrient, numOrient; 580 PetscInt *orients; 581 PetscReal **orientPoints; 582 PetscReal *cellPoints; 583 PetscReal *dummyWeights; 584 PetscQuadrature cellQuad = NULL; 585 586 ierr = DMFieldDSGetHeightDisc(field, 1, &faceDisc);CHKERRQ(ierr); 587 ierr = DMFieldDSGetHeightDisc(field, 0, &cellDisc);CHKERRQ(ierr); 588 ierr = PetscObjectGetClassId(faceDisc,&faceId);CHKERRQ(ierr); 589 ierr = PetscObjectGetClassId(cellDisc,&cellId);CHKERRQ(ierr); 590 if (faceId != PETSCFE_CLASSID || cellId != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported\n"); 591 ierr = PetscFEGetDualSpace((PetscFE)cellDisc, &dsp);CHKERRQ(ierr); 592 ierr = PetscDualSpaceGetDM(dsp, &K); CHKERRQ(ierr); 593 ierr = DMPlexGetConeSize(K,0,&coneSize);CHKERRQ(ierr); 594 ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr); 595 ierr = PetscMalloc4(numFaces,&co,dE*Nq,&cellPoints,coneSize,&counts,Nq,&dummyWeights);CHKERRQ(ierr); 596 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)field), &cellQuad);CHKERRQ(ierr); 597 ierr = PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights);CHKERRQ(ierr); 598 minOrient = PETSC_MAX_INT; 599 maxOrient = PETSC_MIN_INT; 600 for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */ 601 PetscInt point = points[p]; 602 PetscInt numSupp, numChildren; 603 const PetscInt *supp; 604 605 ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL); CHKERRQ(ierr); 606 if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 607 ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 608 if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp); 609 ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 610 for (s = 0; s < numSupp; s++) { 611 PetscInt cell = supp[s]; 612 PetscInt numCone; 613 const PetscInt *cone, *orient; 614 615 ierr = DMPlexGetConeSize(dm, cell, &numCone);CHKERRQ(ierr); 616 if (numCone != coneSize) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element"); 617 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 618 ierr = DMPlexGetConeOrientation(dm, cell, &orient);CHKERRQ(ierr); 619 for (f = 0; f < coneSize; f++) { 620 if (cone[f] == point) break; 621 } 622 co[p][s][0] = f; 623 co[p][s][1] = orient[f]; 624 co[p][s][2] = cell; 625 minOrient = PetscMin(minOrient, orient[f]); 626 maxOrient = PetscMin(maxOrient, orient[f]); 627 } 628 for (; s < 2; s++) { 629 co[p][s][0] = -1; 630 co[p][s][1] = -1; 631 co[p][s][2] = -1; 632 } 633 } 634 numOrient = maxOrient + 1 - minOrient; 635 ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr); 636 /* count all (face,orientation) doubles that appear */ 637 ierr = PetscCalloc2(numOrient,&orients,numOrient,&orientPoints);CHKERRQ(ierr); 638 for (f = 0; f < coneSize; f++) {ierr = PetscCalloc1(numOrient, &counts[f]);CHKERRQ(ierr);} 639 for (p = 0; p < numFaces; p++) { 640 for (s = 0; s < 2; s++) { 641 if (co[p][s][0] >= 0) { 642 counts[co[p][s][0]][co[p][s][1] - minOrient]++; 643 orients[co[p][s][1] - minOrient]++; 644 } 645 } 646 } 647 for (o = 0; o < numOrient; o++) { 648 if (orients[o]) { 649 PetscInt orient = o + minOrient; 650 PetscInt q; 651 652 ierr = PetscMalloc1(Nq * dim, &orientPoints[o]);CHKERRQ(ierr); 653 /* rotate the quadrature points appropriately */ 654 switch (dim) { 655 case 0: 656 break; 657 case 1: 658 if (orient == -2 || orient == 1) { 659 for (q = 0; q < Nq; q++) { 660 orientPoints[o][q] = -geom->xi[q]; 661 } 662 } else { 663 for (q = 0; q < Nq; q++) { 664 orientPoints[o][q] = geom->xi[q]; 665 } 666 } 667 break; 668 case 2: 669 switch (coneSize) { 670 case 3: 671 for (q = 0; q < Nq; q++) { 672 PetscReal lambda[3]; 673 PetscReal lambdao[3]; 674 675 /* convert to barycentric */ 676 lambda[0] = - (geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.; 677 lambda[1] = (geom->xi[2 * q] + 1.) / 2.; 678 lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.; 679 if (orient >= 0) { 680 for (i = 0; i < 3; i++) { 681 lambdao[i] = lambda[(orient + i) % 3]; 682 } 683 } else { 684 for (i = 0; i < 3; i++) { 685 lambdao[i] = lambda[(-(orient + i) + 3) % 3]; 686 } 687 } 688 /* convert to coordinates */ 689 orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1]; 690 orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2]; 691 } 692 break; 693 case 4: 694 for (q = 0; q < Nq; q++) { 695 PetscReal xi[2], xio[2]; 696 PetscInt oabs = (orient >= 0) ? orient : -(orient + 1); 697 698 xi[0] = geom->xi[2 * q]; 699 xi[1] = geom->xi[2 * q + 1]; 700 switch (oabs) { 701 case 0: 702 xio[0] = xi[0]; 703 xio[1] = xi[1]; 704 break; 705 case 1: 706 xio[0] = xi[1]; 707 xio[1] = -xi[0]; 708 break; 709 case 2: 710 xio[0] = -xi[0]; 711 xio[1] = -xi[1]; 712 case 3: 713 xio[0] = -xi[1]; 714 xio[1] = xi[0]; 715 } 716 if (orient < 0) { 717 xio[0] = -xio[0]; 718 } 719 orientPoints[o][2 * q + 0] = xio[0]; 720 orientPoints[o][2 * q + 1] = xio[1]; 721 } 722 break; 723 default: 724 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not implemented yet\n"); 725 } 726 default: 727 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not implemented yet\n"); 728 } 729 } 730 } 731 for (f = 0; f < coneSize; f++) { 732 PetscInt face = coneK[f]; 733 PetscScalar v0[3]; 734 PetscScalar J[9]; 735 PetscInt numCells, offset; 736 PetscInt *cells; 737 IS suppIS; 738 739 ierr = DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, NULL);CHKERRQ(ierr); 740 for (o = 0; o <= numOrient; o++) { 741 PetscFEGeom *cellGeom; 742 743 if (!counts[f][o]) continue; 744 /* If this (face,orientation) double appears, 745 * convert the face quadrature points into volume quadrature points */ 746 for (q = 0; q < Nq; q++) { 747 PetscReal xi0[3] = {-1., -1., -1.}; 748 749 CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]); 750 } 751 for (p = 0, numCells = 0; p < numFaces; p++) { 752 for (s = 0; s < 2; s++) { 753 if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++; 754 } 755 } 756 ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr); 757 for (p = 0, offset = 0; p < numFaces; p++) { 758 for (s = 0; s < 2; s++) { 759 if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 760 cells[offset++] = co[p][s][2]; 761 } 762 } 763 } 764 ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr); 765 ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr); 766 for (p = 0, offset = 0; p < numFaces; p++) { 767 for (s = 0; s < 2; s++) { 768 if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 769 for (q = 0; q < Nq * dE * dE; q++) { 770 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 771 } 772 offset++; 773 } 774 } 775 } 776 ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr); 777 ierr = ISDestroy(&suppIS);CHKERRQ(ierr); 778 ierr = PetscFree(cells);CHKERRQ(ierr); 779 } 780 } 781 for (o = 0; o < numOrient; o++) { 782 if (orients[o]) { 783 ierr = PetscFree(orientPoints[o]);CHKERRQ(ierr); 784 } 785 } 786 ierr = PetscFree2(orients,orientPoints);CHKERRQ(ierr); 787 ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr); 788 ierr = PetscFree4(co,cellPoints,counts,dummyWeights);CHKERRQ(ierr); 789 } 790 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 791 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 792 PetscFunctionReturn(0); 793 } 794 795 static PetscErrorCode DMFieldInitialize_DS(DMField field) 796 { 797 PetscFunctionBegin; 798 field->ops->destroy = DMFieldDestroy_DS; 799 field->ops->evaluate = DMFieldEvaluate_DS; 800 field->ops->evaluateFE = DMFieldEvaluateFE_DS; 801 field->ops->getFEInvariance = DMFieldGetFEInvariance_DS; 802 field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS; 803 field->ops->view = DMFieldView_DS; 804 field->ops->computeFaceData = DMFieldComputeFaceData_DS; 805 PetscFunctionReturn(0); 806 } 807 808 PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field) 809 { 810 DMField_DS *dsfield; 811 PetscErrorCode ierr; 812 813 PetscFunctionBegin; 814 ierr = PetscNewLog(field,&dsfield);CHKERRQ(ierr); 815 field->data = dsfield; 816 ierr = DMFieldInitialize_DS(field);CHKERRQ(ierr); 817 PetscFunctionReturn(0); 818 } 819 820 PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec,DMField *field) 821 { 822 DMField b; 823 DMField_DS *dsfield; 824 PetscObject disc = NULL; 825 PetscBool isContainer = PETSC_FALSE; 826 PetscClassId id = -1; 827 PetscInt numComponents = -1, dsNumFields; 828 PetscSection section; 829 PetscErrorCode ierr; 830 831 PetscFunctionBegin; 832 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 833 ierr = PetscSectionGetFieldComponents(section,fieldNum,&numComponents);CHKERRQ(ierr); 834 ierr = DMGetNumFields(dm,&dsNumFields);CHKERRQ(ierr); 835 if (dsNumFields) {ierr = DMGetField(dm,fieldNum,&disc);CHKERRQ(ierr);} 836 if (disc) { 837 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 838 isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE; 839 } 840 if (!disc || isContainer) { 841 MPI_Comm comm = PetscObjectComm((PetscObject) dm); 842 PetscInt cStart, cEnd, dim; 843 PetscInt localConeSize = 0, coneSize; 844 PetscFE fe; 845 PetscDualSpace Q; 846 PetscSpace P; 847 DM K; 848 PetscQuadrature quad, fquad; 849 PetscBool isSimplex; 850 851 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 852 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 853 if (cEnd > cStart) { 854 ierr = DMPlexGetConeSize(dm, cStart, &localConeSize);CHKERRQ(ierr); 855 } 856 ierr = MPI_Allreduce(&localConeSize,&coneSize,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 857 isSimplex = (coneSize == (dim + 1)) ? PETSC_TRUE : PETSC_FALSE; 858 ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 859 ierr = PetscSpaceSetOrder(P, 1);CHKERRQ(ierr); 860 ierr = PetscSpaceSetNumComponents(P, numComponents);CHKERRQ(ierr); 861 ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 862 ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 863 ierr = PetscSpacePolynomialSetTensor(P, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr); 864 ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 865 ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 866 ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 867 ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 868 ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 869 ierr = DMDestroy(&K);CHKERRQ(ierr); 870 ierr = PetscDualSpaceSetNumComponents(Q, numComponents);CHKERRQ(ierr); 871 ierr = PetscDualSpaceSetOrder(Q, 1);CHKERRQ(ierr); 872 ierr = PetscDualSpaceLagrangeSetTensor(Q, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr); 873 ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 874 ierr = PetscFECreate(comm, &fe);CHKERRQ(ierr); 875 ierr = PetscFESetType(fe,PETSCFEBASIC);CHKERRQ(ierr); 876 ierr = PetscFESetBasisSpace(fe, P);CHKERRQ(ierr); 877 ierr = PetscFESetDualSpace(fe, Q);CHKERRQ(ierr); 878 ierr = PetscFESetNumComponents(fe, numComponents);CHKERRQ(ierr); 879 ierr = PetscFESetUp(fe);CHKERRQ(ierr); 880 ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 881 ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 882 if (isSimplex) { 883 ierr = PetscDTGaussJacobiQuadrature(dim, 1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr); 884 ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr); 885 } 886 else { 887 ierr = PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr); 888 ierr = PetscDTGaussTensorQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr); 889 } 890 ierr = PetscFESetQuadrature(fe, quad);CHKERRQ(ierr); 891 ierr = PetscFESetFaceQuadrature(fe, fquad);CHKERRQ(ierr); 892 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 893 ierr = PetscQuadratureDestroy(&fquad);CHKERRQ(ierr); 894 disc = (PetscObject) fe; 895 } else { 896 ierr = PetscObjectReference(disc);CHKERRQ(ierr); 897 } 898 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 899 if (id == PETSCFE_CLASSID) { 900 PetscFE fe = (PetscFE) disc; 901 902 ierr = PetscFEGetNumComponents(fe,&numComponents);CHKERRQ(ierr); 903 } else {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented");} 904 ierr = DMFieldCreate(dm,numComponents,DMFIELD_VERTEX,&b);CHKERRQ(ierr); 905 ierr = DMFieldSetType(b,DMFIELDDS);CHKERRQ(ierr); 906 dsfield = (DMField_DS *) b->data; 907 dsfield->fieldNum = fieldNum; 908 ierr = DMGetDimension(dm,&dsfield->height);CHKERRQ(ierr); 909 dsfield->height++; 910 ierr = PetscCalloc1(dsfield->height,&dsfield->disc);CHKERRQ(ierr); 911 dsfield->disc[0] = disc; 912 ierr = PetscObjectReference((PetscObject)vec);CHKERRQ(ierr); 913 dsfield->vec = vec; 914 *field = b; 915 916 PetscFunctionReturn(0); 917 } 918