151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/ 20145028aSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/ 351a454edSToby Isaac #include <petscfe.h> 451a454edSToby Isaac #include <petscdmplex.h> 551a454edSToby Isaac #include <petscds.h> 651a454edSToby Isaac 751a454edSToby Isaac typedef struct _n_DMField_DS 851a454edSToby Isaac { 951a454edSToby Isaac PetscInt fieldNum; 1051a454edSToby Isaac Vec vec; 1151a454edSToby Isaac PetscInt height; 1251a454edSToby Isaac PetscObject *disc; 1351a454edSToby Isaac PetscBool multifieldVec; 1451a454edSToby Isaac } 1551a454edSToby Isaac DMField_DS; 1651a454edSToby Isaac 1751a454edSToby Isaac static PetscErrorCode DMFieldDestroy_DS(DMField field) 1851a454edSToby Isaac { 1951a454edSToby Isaac DMField_DS *dsfield; 2051a454edSToby Isaac PetscInt i; 2151a454edSToby Isaac PetscErrorCode ierr; 2251a454edSToby Isaac 2351a454edSToby Isaac PetscFunctionBegin; 2451a454edSToby Isaac dsfield = (DMField_DS *) field->data; 2551a454edSToby Isaac ierr = VecDestroy(&dsfield->vec);CHKERRQ(ierr); 2651a454edSToby Isaac for (i = 0; i < dsfield->height; i++) { 2751a454edSToby Isaac ierr = PetscObjectDereference(dsfield->disc[i]);CHKERRQ(ierr); 2851a454edSToby Isaac } 2951a454edSToby Isaac ierr = PetscFree(dsfield->disc);CHKERRQ(ierr); 3051a454edSToby Isaac ierr = PetscFree(dsfield);CHKERRQ(ierr); 3151a454edSToby Isaac PetscFunctionReturn(0); 3251a454edSToby Isaac } 3351a454edSToby Isaac 3451a454edSToby Isaac static PetscErrorCode DMFieldView_DS(DMField field,PetscViewer viewer) 3551a454edSToby Isaac { 3651a454edSToby Isaac DMField_DS *dsfield = (DMField_DS *) field->data; 3751a454edSToby Isaac PetscBool iascii; 3851a454edSToby Isaac PetscObject disc; 3951a454edSToby Isaac PetscErrorCode ierr; 4051a454edSToby Isaac 4151a454edSToby Isaac PetscFunctionBegin; 4251a454edSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4351a454edSToby Isaac disc = dsfield->disc[0]; 4451a454edSToby Isaac if (iascii) { 4551a454edSToby Isaac PetscViewerASCIIPrintf(viewer, "PetscDS field %D\n",dsfield->fieldNum);CHKERRQ(ierr); 4651a454edSToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4751a454edSToby Isaac ierr = PetscObjectView(disc,viewer);CHKERRQ(ierr); 4851a454edSToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4951a454edSToby Isaac } 5051a454edSToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 5151a454edSToby Isaac if (dsfield->multifieldVec) { 5251a454edSToby Isaac SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"View of subfield not implemented yet"); 5351a454edSToby Isaac } else { 5451a454edSToby Isaac ierr = VecView(dsfield->vec,viewer);CHKERRQ(ierr); 5551a454edSToby Isaac } 5651a454edSToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 5751a454edSToby Isaac PetscFunctionReturn(0); 5851a454edSToby Isaac } 5951a454edSToby Isaac 6051a454edSToby Isaac static PetscErrorCode DMFieldDSGetHeightDisc(DMField field, PetscInt height, PetscObject *disc) 6151a454edSToby Isaac { 6251a454edSToby Isaac DMField_DS *dsfield = (DMField_DS *) field->data; 6351a454edSToby Isaac PetscErrorCode ierr; 6451a454edSToby Isaac 6551a454edSToby Isaac PetscFunctionBegin; 6651a454edSToby Isaac if (!dsfield->disc[height]) { 6751a454edSToby Isaac PetscClassId id; 6851a454edSToby Isaac 6951a454edSToby Isaac ierr = PetscObjectGetClassId(dsfield->disc[0],&id);CHKERRQ(ierr); 7051a454edSToby Isaac if (id == PETSCFE_CLASSID) { 7151a454edSToby Isaac PetscFE fe = (PetscFE) dsfield->disc[0]; 7251a454edSToby Isaac 7351a454edSToby Isaac ierr = PetscFECreateHeightTrace(fe,height,(PetscFE *)&dsfield->disc[height]);CHKERRQ(ierr); 7451a454edSToby Isaac } 7551a454edSToby Isaac } 7651a454edSToby Isaac *disc = dsfield->disc[height]; 7751a454edSToby Isaac PetscFunctionReturn(0); 7851a454edSToby Isaac } 7951a454edSToby Isaac 8051a454edSToby Isaac #define DMFieldDSdot(y,A,b,m,n,c,cast) \ 8151a454edSToby Isaac do { \ 8251a454edSToby Isaac PetscInt _i, _j, _k; \ 8351a454edSToby Isaac for (_i = 0; _i < (m); _i++) { \ 8451a454edSToby Isaac for (_k = 0; _k < (c); _k++) { \ 8551a454edSToby Isaac (y)[_i * (c) + _k] = 0.; \ 8651a454edSToby Isaac } \ 8751a454edSToby Isaac for (_j = 0; _j < (n); _j++) { \ 8851a454edSToby Isaac for (_k = 0; _k < (c); _k++) { \ 8951a454edSToby Isaac (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \ 9051a454edSToby Isaac } \ 9151a454edSToby Isaac } \ 9251a454edSToby Isaac } \ 9351a454edSToby Isaac } while (0) 9451a454edSToby Isaac 9551a454edSToby Isaac static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H) 9651a454edSToby Isaac { 9751a454edSToby Isaac DMField_DS *dsfield = (DMField_DS *) field->data; 9851a454edSToby Isaac DM dm; 9951a454edSToby Isaac PetscObject disc; 10051a454edSToby Isaac PetscClassId classid; 10151a454edSToby Isaac PetscInt nq, nc, dim, meshDim, numCells; 10251a454edSToby Isaac PetscSection section; 10351a454edSToby Isaac const PetscReal *qpoints; 10451a454edSToby Isaac PetscBool isStride; 10551a454edSToby Isaac const PetscInt *points = NULL; 10651a454edSToby Isaac PetscInt sfirst = -1, stride = -1; 10751a454edSToby Isaac PetscErrorCode ierr; 10851a454edSToby Isaac 10951a454edSToby Isaac PetscFunctionBeginHot; 11051a454edSToby Isaac dm = field->dm; 11151a454edSToby Isaac nc = field->numComponents; 11251a454edSToby Isaac ierr = PetscQuadratureGetData(quad,&dim,NULL,&nq,&qpoints,NULL);CHKERRQ(ierr); 113f99c8401SToby Isaac ierr = DMFieldDSGetHeightDisc(field,dsfield->height - 1 - dim,&disc);CHKERRQ(ierr); 11451a454edSToby Isaac ierr = DMGetDimension(dm,&meshDim);CHKERRQ(ierr); 11551a454edSToby Isaac ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 11651a454edSToby Isaac ierr = PetscSectionGetField(section,dsfield->fieldNum,§ion);CHKERRQ(ierr); 11751a454edSToby Isaac ierr = PetscObjectGetClassId(disc,&classid);CHKERRQ(ierr); 11851a454edSToby Isaac /* TODO: batch */ 11951a454edSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)pointIS,ISSTRIDE,&isStride);CHKERRQ(ierr); 12051a454edSToby Isaac ierr = ISGetLocalSize(pointIS,&numCells);CHKERRQ(ierr); 12151a454edSToby Isaac if (isStride) { 12251a454edSToby Isaac ierr = ISStrideGetInfo(pointIS,&sfirst,&stride);CHKERRQ(ierr); 12351a454edSToby Isaac } else { 12451a454edSToby Isaac ierr = ISGetIndices(pointIS,&points);CHKERRQ(ierr); 12551a454edSToby Isaac } 12651a454edSToby Isaac if (classid == PETSCFE_CLASSID) { 12751a454edSToby Isaac PetscFE fe = (PetscFE) disc; 12851a454edSToby Isaac PetscInt feDim, i; 12951a454edSToby Isaac PetscReal *fB = NULL, *fD = NULL, *fH = NULL; 13051a454edSToby Isaac 13151a454edSToby Isaac if (dim == meshDim - 1) { 13251a454edSToby Isaac /* TODO */ 13351a454edSToby Isaac } 13451a454edSToby Isaac ierr = PetscFEGetDimension(fe,&feDim);CHKERRQ(ierr); 13551a454edSToby Isaac ierr = PetscFEGetTabulation(fe,nq,qpoints,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 13651a454edSToby Isaac for (i = 0; i < numCells; i++) { 13751a454edSToby Isaac PetscInt c = isStride ? (sfirst + i * stride) : points[i]; 1382710589cSToby Isaac PetscInt closureSize; 1392710589cSToby Isaac PetscScalar *elem = NULL; 14051a454edSToby Isaac 14151a454edSToby Isaac ierr = DMPlexVecGetClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 14251a454edSToby Isaac if (B) { 14351a454edSToby Isaac if (type == PETSC_SCALAR) { 14451a454edSToby Isaac PetscScalar *cB = &((PetscScalar *) B)[nc * nq * i]; 14551a454edSToby Isaac 14651a454edSToby Isaac DMFieldDSdot(cB,fB,elem,nq,feDim,nc,(PetscScalar)); 14751a454edSToby Isaac } else { 14851a454edSToby Isaac PetscReal *cB = &((PetscReal *) B)[nc * nq * i]; 14951a454edSToby Isaac 15051a454edSToby Isaac DMFieldDSdot(cB,fB,elem,nq,feDim,nc,PetscRealPart); 15151a454edSToby Isaac } 15251a454edSToby Isaac } 15351a454edSToby Isaac if (D) { 15451a454edSToby Isaac if (type == PETSC_SCALAR) { 15551a454edSToby Isaac PetscScalar *cD = &((PetscScalar *) D)[nc * nq * dim * i]; 15651a454edSToby Isaac 15751a454edSToby Isaac DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),(PetscScalar)); 15851a454edSToby Isaac } else { 15951a454edSToby Isaac PetscReal *cD = &((PetscReal *) D)[nc * nq * dim * i]; 16051a454edSToby Isaac 16151a454edSToby Isaac DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),PetscRealPart); 16251a454edSToby Isaac } 16351a454edSToby Isaac } 16451a454edSToby Isaac if (H) { 16551a454edSToby Isaac if (type == PETSC_SCALAR) { 16651a454edSToby Isaac PetscScalar *cH = &((PetscScalar *) H)[nc * nq * dim * dim * i]; 16751a454edSToby Isaac 16851a454edSToby Isaac DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),(PetscScalar)); 16951a454edSToby Isaac } else { 17051a454edSToby Isaac PetscReal *cH = &((PetscReal *) H)[nc * nq * dim * dim * i]; 17151a454edSToby Isaac 17251a454edSToby Isaac DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),PetscRealPart); 17351a454edSToby Isaac } 17451a454edSToby Isaac } 1752710589cSToby Isaac ierr = DMPlexVecRestoreClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 17651a454edSToby Isaac } 17751a454edSToby Isaac ierr = PetscFERestoreTabulation(fe,nq,qpoints,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 17851a454edSToby Isaac } else {SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented");} 17951a454edSToby Isaac if (!isStride) { 18051a454edSToby Isaac ierr = ISRestoreIndices(pointIS,&points);CHKERRQ(ierr); 18151a454edSToby Isaac } 18251a454edSToby Isaac ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 18351a454edSToby Isaac PetscFunctionReturn(0); 18451a454edSToby Isaac } 18551a454edSToby Isaac 186a6216909SToby Isaac static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 187a6216909SToby Isaac { 188a6216909SToby Isaac DMField_DS *dsfield = (DMField_DS *) field->data; 189a6216909SToby Isaac PetscSF cellSF = NULL; 190a6216909SToby Isaac const PetscSFNode *cells; 191a6216909SToby Isaac PetscInt c, nFound, numCells, feDim, nc; 192a6216909SToby Isaac const PetscInt *cellDegrees; 193a6216909SToby Isaac const PetscScalar *pointsArray; 194a6216909SToby Isaac PetscScalar *cellPoints; 195a6216909SToby Isaac PetscInt gatherSize, gatherMax; 196a6216909SToby Isaac PetscInt dim, dimR, offset; 197a6216909SToby Isaac MPI_Datatype pointType; 198a6216909SToby Isaac PetscObject cellDisc; 199a6216909SToby Isaac PetscFE cellFE; 200a6216909SToby Isaac PetscClassId discID; 201a6216909SToby Isaac PetscReal *coordsReal, *coordsRef; 202a6216909SToby Isaac PetscSection section; 203a6216909SToby Isaac PetscScalar *cellBs = NULL, *cellDs = NULL, *cellHs = NULL; 204a6216909SToby Isaac PetscReal *cellBr = NULL, *cellDr = NULL, *cellHr = NULL; 205*4cbe5319SToby Isaac PetscReal *v, *J, *invJ, *detJ; 206a6216909SToby Isaac PetscErrorCode ierr; 207a6216909SToby Isaac 208a6216909SToby Isaac PetscFunctionBegin; 209a6216909SToby Isaac nc = field->numComponents; 210a6216909SToby Isaac ierr = DMGetDefaultSection(field->dm,§ion);CHKERRQ(ierr); 211a6216909SToby Isaac ierr = DMFieldDSGetHeightDisc(field,0,&cellDisc);CHKERRQ(ierr); 212a6216909SToby Isaac ierr = PetscObjectGetClassId(cellDisc, &discID);CHKERRQ(ierr); 213a6216909SToby Isaac if (discID != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Discretization type not supported\n"); 214a6216909SToby Isaac cellFE = (PetscFE) cellDisc; 215a6216909SToby Isaac ierr = PetscFEGetDimension(cellFE,&feDim);CHKERRQ(ierr); 216a6216909SToby Isaac ierr = DMGetCoordinateDim(field->dm, &dim);CHKERRQ(ierr); 217a6216909SToby Isaac ierr = DMGetDimension(field->dm, &dimR);CHKERRQ(ierr); 218a6216909SToby Isaac ierr = DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF);CHKERRQ(ierr); 219a6216909SToby Isaac ierr = PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells);CHKERRQ(ierr); 220a6216909SToby Isaac for (c = 0; c < nFound; c++) { 221a6216909SToby Isaac if (cells[c].index < 0) SETERRQ1(PetscObjectComm((PetscObject)points),PETSC_ERR_ARG_WRONG, "Point %D could not be located\n", c); 222a6216909SToby Isaac } 223a6216909SToby Isaac ierr = PetscSFComputeDegreeBegin(cellSF,&cellDegrees);CHKERRQ(ierr); 224a6216909SToby Isaac ierr = PetscSFComputeDegreeEnd(cellSF,&cellDegrees);CHKERRQ(ierr); 225a6216909SToby Isaac for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) { 226a6216909SToby Isaac gatherMax = PetscMax(gatherMax,cellDegrees[c]); 227a6216909SToby Isaac gatherSize += cellDegrees[c]; 228a6216909SToby Isaac } 229a6216909SToby Isaac ierr = PetscMalloc3(gatherSize*dim,&cellPoints,gatherMax*dim,&coordsReal,gatherMax*dimR,&coordsRef);CHKERRQ(ierr); 230*4cbe5319SToby Isaac ierr = PetscMalloc4(gatherMax*dimR,&v,gatherMax*dimR*dimR,&J,gatherMax*dimR*dimR,&invJ,gatherMax,&detJ);CHKERRQ(ierr); 231a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 232a6216909SToby Isaac ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBs, D ? nc * dim * gatherSize : 0, &cellDs, H ? nc * dim * dim * gatherSize : 0, &cellHs);CHKERRQ(ierr); 233a6216909SToby Isaac } else { 234a6216909SToby Isaac ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr);CHKERRQ(ierr); 235a6216909SToby Isaac } 236a6216909SToby Isaac 237a6216909SToby Isaac ierr = MPI_Type_contiguous(dim,MPIU_SCALAR,&pointType);CHKERRQ(ierr); 238a6216909SToby Isaac ierr = MPI_Type_commit(&pointType);CHKERRQ(ierr); 239a6216909SToby Isaac ierr = VecGetArrayRead(points,&pointsArray);CHKERRQ(ierr); 240a6216909SToby Isaac ierr = PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr); 241a6216909SToby Isaac ierr = PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr); 242a6216909SToby Isaac ierr = VecRestoreArrayRead(points,&pointsArray);CHKERRQ(ierr); 243a6216909SToby Isaac for (c = 0, offset = 0; c < numCells; c++) { 244a6216909SToby Isaac PetscInt nq = cellDegrees[c], p; 245a6216909SToby Isaac 246a6216909SToby Isaac if (nq) { 247a6216909SToby Isaac PetscReal *fB, *fD, *fH; 248a6216909SToby Isaac PetscInt closureSize; 249a6216909SToby Isaac PetscScalar *elem = NULL; 250*4cbe5319SToby Isaac PetscReal *quadPoints; 251*4cbe5319SToby Isaac PetscQuadrature quad; 252*4cbe5319SToby Isaac PetscInt d, e, f, g; 253a6216909SToby Isaac 254a6216909SToby Isaac for (p = 0; p < dim * nq; p++) coordsReal[p] = cellPoints[dim * offset + p]; 255a6216909SToby Isaac ierr = DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef);CHKERRQ(ierr); 256a6216909SToby Isaac ierr = PetscFEGetTabulation(cellFE,nq,coordsRef,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 257*4cbe5319SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &quad);CHKERRQ(ierr); 258*4cbe5319SToby Isaac ierr = PetscMalloc1(dimR * nq, &quadPoints);CHKERRQ(ierr); 259*4cbe5319SToby Isaac for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p]; 260*4cbe5319SToby Isaac ierr = PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL);CHKERRQ(ierr); 261*4cbe5319SToby Isaac ierr = DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ);CHKERRQ(ierr); 262*4cbe5319SToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 263a6216909SToby Isaac ierr = DMPlexVecGetClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 264a6216909SToby Isaac if (B) { 265a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 266a6216909SToby Isaac PetscScalar *cB = &cellBs[nc * offset]; 267a6216909SToby Isaac 268a6216909SToby Isaac DMFieldDSdot(cB,fB,elem,nq,feDim,nc,(PetscScalar)); 269a6216909SToby Isaac } else { 270a6216909SToby Isaac PetscReal *cB = &cellBr[nc * offset]; 271a6216909SToby Isaac 272a6216909SToby Isaac DMFieldDSdot(cB,fB,elem,nq,feDim,nc,PetscRealPart); 273a6216909SToby Isaac } 274a6216909SToby Isaac } 275a6216909SToby Isaac if (D) { 276a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 277a6216909SToby Isaac PetscScalar *cD = &cellDs[nc * dim * offset]; 278a6216909SToby Isaac 279a6216909SToby Isaac DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),(PetscScalar)); 280*4cbe5319SToby Isaac for (p = 0; p < nq; p++) { 281*4cbe5319SToby Isaac for (g = 0; g < nc; g++) { 282*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 283*4cbe5319SToby Isaac v[d] = 0.; 284*4cbe5319SToby Isaac for (e = 0; e < dimR; e++) { 285*4cbe5319SToby Isaac v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 286*4cbe5319SToby Isaac } 287*4cbe5319SToby Isaac } 288*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 289*4cbe5319SToby Isaac cD[(nc * p + g) * dimR + d] = v[d]; 290*4cbe5319SToby Isaac } 291*4cbe5319SToby Isaac } 292*4cbe5319SToby Isaac } 293a6216909SToby Isaac } else { 294a6216909SToby Isaac PetscReal *cD = &cellDr[nc * dim * offset]; 295a6216909SToby Isaac 296a6216909SToby Isaac DMFieldDSdot(cD,fD,elem,nq,feDim,(nc * dim),PetscRealPart); 297*4cbe5319SToby Isaac for (p = 0; p < nq; p++) { 298*4cbe5319SToby Isaac for (g = 0; g < nc; g++) { 299*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 300*4cbe5319SToby Isaac v[d] = 0.; 301*4cbe5319SToby Isaac for (e = 0; e < dimR; e++) { 302*4cbe5319SToby Isaac v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 303*4cbe5319SToby Isaac } 304*4cbe5319SToby Isaac } 305*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 306*4cbe5319SToby Isaac cD[(nc * p + g) * dimR + d] = v[d]; 307*4cbe5319SToby Isaac } 308*4cbe5319SToby Isaac } 309*4cbe5319SToby Isaac } 310a6216909SToby Isaac } 311a6216909SToby Isaac } 312a6216909SToby Isaac if (H) { 313a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 314a6216909SToby Isaac PetscScalar *cH = &cellHs[nc * dim * dim * offset]; 315a6216909SToby Isaac 316a6216909SToby Isaac DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),(PetscScalar)); 317*4cbe5319SToby Isaac for (p = 0; p < nq; p++) { 318*4cbe5319SToby Isaac for (g = 0; g < nc * dimR; g++) { 319*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 320*4cbe5319SToby Isaac v[d] = 0.; 321*4cbe5319SToby Isaac for (e = 0; e < dimR; e++) { 322*4cbe5319SToby Isaac v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e]; 323*4cbe5319SToby Isaac } 324*4cbe5319SToby Isaac } 325*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 326*4cbe5319SToby Isaac cH[(nc * dimR * p + g) * dimR + d] = v[d]; 327*4cbe5319SToby Isaac } 328*4cbe5319SToby Isaac } 329*4cbe5319SToby Isaac for (g = 0; g < nc; g++) { 330*4cbe5319SToby Isaac for (f = 0; f < dimR; f++) { 331*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 332*4cbe5319SToby Isaac v[d] = 0.; 333*4cbe5319SToby Isaac for (e = 0; e < dimR; e++) { 334*4cbe5319SToby Isaac v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 335*4cbe5319SToby Isaac } 336*4cbe5319SToby Isaac } 337*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 338*4cbe5319SToby Isaac cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 339*4cbe5319SToby Isaac } 340*4cbe5319SToby Isaac } 341*4cbe5319SToby Isaac } 342*4cbe5319SToby Isaac } 343a6216909SToby Isaac } else { 344a6216909SToby Isaac PetscReal *cH = &cellHr[nc * dim * dim * offset]; 345a6216909SToby Isaac 346a6216909SToby Isaac DMFieldDSdot(cH,fH,elem,nq,feDim,(nc * dim * dim),PetscRealPart); 347*4cbe5319SToby Isaac for (p = 0; p < nq; p++) { 348*4cbe5319SToby Isaac for (g = 0; g < nc * dimR; g++) { 349*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 350*4cbe5319SToby Isaac v[d] = 0.; 351*4cbe5319SToby Isaac for (e = 0; e < dimR; e++) { 352*4cbe5319SToby Isaac v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e]; 353*4cbe5319SToby Isaac } 354*4cbe5319SToby Isaac } 355*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 356*4cbe5319SToby Isaac cH[(nc * dimR * p + g) * dimR + d] = v[d]; 357*4cbe5319SToby Isaac } 358*4cbe5319SToby Isaac } 359*4cbe5319SToby Isaac for (g = 0; g < nc; g++) { 360*4cbe5319SToby Isaac for (f = 0; f < dimR; f++) { 361*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 362*4cbe5319SToby Isaac v[d] = 0.; 363*4cbe5319SToby Isaac for (e = 0; e < dimR; e++) { 364*4cbe5319SToby Isaac v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 365*4cbe5319SToby Isaac } 366*4cbe5319SToby Isaac } 367*4cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 368*4cbe5319SToby Isaac cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 369*4cbe5319SToby Isaac } 370*4cbe5319SToby Isaac } 371*4cbe5319SToby Isaac } 372*4cbe5319SToby Isaac } 373a6216909SToby Isaac } 374a6216909SToby Isaac } 375a6216909SToby Isaac ierr = DMPlexVecRestoreClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr); 376*4cbe5319SToby Isaac ierr = PetscFERestoreTabulation(cellFE,nq,coordsRef,B ? &fB : NULL,D ? &fD : NULL,H ? &fH : NULL);CHKERRQ(ierr); 377a6216909SToby Isaac } 378a6216909SToby Isaac offset += nq; 379a6216909SToby Isaac } 380a6216909SToby Isaac { 381a6216909SToby Isaac MPI_Datatype origtype; 382*4cbe5319SToby Isaac 383a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 384a6216909SToby Isaac origtype = MPIU_SCALAR; 385a6216909SToby Isaac } else { 386a6216909SToby Isaac origtype = MPIU_REAL; 387a6216909SToby Isaac } 388a6216909SToby Isaac if (B) { 389a6216909SToby Isaac MPI_Datatype Btype; 390a6216909SToby Isaac 391*4cbe5319SToby Isaac ierr = MPI_Type_contiguous(nc, origtype, &Btype);CHKERRQ(ierr); 392a6216909SToby Isaac ierr = MPI_Type_commit(&Btype);CHKERRQ(ierr); 393a6216909SToby Isaac ierr = PetscSFScatterBegin(cellSF,Btype,(datatype == PETSC_SCALAR) ? cellBs : cellBr, B);CHKERRQ(ierr); 394a6216909SToby Isaac ierr = PetscSFScatterEnd(cellSF,Btype,(datatype == PETSC_SCALAR) ? cellBs : cellBr, B);CHKERRQ(ierr); 395a6216909SToby Isaac ierr = MPI_Type_free(&Btype);CHKERRQ(ierr); 396a6216909SToby Isaac } 397a6216909SToby Isaac if (D) { 398a6216909SToby Isaac MPI_Datatype Dtype; 399a6216909SToby Isaac 400*4cbe5319SToby Isaac ierr = MPI_Type_contiguous(nc * dim, origtype, &Dtype);CHKERRQ(ierr); 401a6216909SToby Isaac ierr = MPI_Type_commit(&Dtype);CHKERRQ(ierr); 402a6216909SToby Isaac ierr = PetscSFScatterBegin(cellSF,Dtype,(datatype == PETSC_SCALAR) ? cellDs : cellDr, D);CHKERRQ(ierr); 403a6216909SToby Isaac ierr = PetscSFScatterEnd(cellSF,Dtype,(datatype == PETSC_SCALAR) ? cellDs : cellDr, D);CHKERRQ(ierr); 404a6216909SToby Isaac ierr = MPI_Type_free(&Dtype);CHKERRQ(ierr); 405a6216909SToby Isaac } 406a6216909SToby Isaac if (H) { 407a6216909SToby Isaac MPI_Datatype Htype; 408a6216909SToby Isaac 409*4cbe5319SToby Isaac ierr = MPI_Type_contiguous(nc * dim * dim, origtype, &Htype);CHKERRQ(ierr); 410a6216909SToby Isaac ierr = MPI_Type_commit(&Htype);CHKERRQ(ierr); 411a6216909SToby Isaac ierr = PetscSFScatterBegin(cellSF,Htype,(datatype == PETSC_SCALAR) ? cellHs : cellHr, H);CHKERRQ(ierr); 412a6216909SToby Isaac ierr = PetscSFScatterEnd(cellSF,Htype,(datatype == PETSC_SCALAR) ? cellHs : cellHr, H);CHKERRQ(ierr); 413a6216909SToby Isaac ierr = MPI_Type_free(&Htype);CHKERRQ(ierr); 414a6216909SToby Isaac } 415a6216909SToby Isaac } 416*4cbe5319SToby Isaac ierr = PetscFree4(v,J,invJ,detJ);CHKERRQ(ierr); 417a6216909SToby Isaac ierr = PetscFree3(cellBr, cellDr, cellHr);CHKERRQ(ierr); 418a6216909SToby Isaac ierr = PetscFree3(cellBs, cellDs, cellHs);CHKERRQ(ierr); 419a6216909SToby Isaac ierr = PetscFree3(cellPoints,coordsReal,coordsRef);CHKERRQ(ierr); 420a6216909SToby Isaac ierr = MPI_Type_free(&pointType);CHKERRQ(ierr); 421a6216909SToby Isaac ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr); 422a6216909SToby Isaac PetscFunctionReturn(0); 423a6216909SToby Isaac } 424a6216909SToby Isaac 42551a454edSToby Isaac static PetscErrorCode DMFieldGetFEInvariance_DS(DMField field, IS pointIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic) 42651a454edSToby Isaac { 42751a454edSToby Isaac DMField_DS *dsfield; 42851a454edSToby Isaac PetscObject disc; 42951a454edSToby Isaac PetscInt h, imin; 43051a454edSToby Isaac PetscClassId id; 43151a454edSToby Isaac PetscErrorCode ierr; 43251a454edSToby Isaac 43351a454edSToby Isaac PetscFunctionBegin; 43451a454edSToby Isaac dsfield = (DMField_DS *) field->data; 43551a454edSToby Isaac ierr = ISGetMinMax(pointIS,&imin,NULL);CHKERRQ(ierr); 436f99c8401SToby Isaac for (h = 0; h < dsfield->height; h++) { 43751a454edSToby Isaac PetscInt hEnd; 43851a454edSToby Isaac 43951a454edSToby Isaac ierr = DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd);CHKERRQ(ierr); 44051a454edSToby Isaac if (imin < hEnd) break; 44151a454edSToby Isaac } 44251a454edSToby Isaac ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 44351a454edSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 44451a454edSToby Isaac if (id == PETSCFE_CLASSID) { 44551a454edSToby Isaac PetscFE fe = (PetscFE) disc; 44651a454edSToby Isaac PetscInt order, maxOrder; 44751a454edSToby Isaac PetscBool tensor = PETSC_FALSE; 44851a454edSToby Isaac PetscSpace sp; 44951a454edSToby Isaac 45051a454edSToby Isaac ierr = PetscFEGetBasisSpace(fe, &sp);CHKERRQ(ierr); 45151a454edSToby Isaac ierr = PetscSpaceGetOrder(sp,&order);CHKERRQ(ierr); 45251a454edSToby Isaac ierr = PetscSpacePolynomialGetTensor(sp,&tensor);CHKERRQ(ierr); 45351a454edSToby Isaac if (tensor) { 45451a454edSToby Isaac PetscInt dim; 45551a454edSToby Isaac 45651a454edSToby Isaac ierr = DMGetDimension(field->dm,&dim);CHKERRQ(ierr); 45751a454edSToby Isaac maxOrder = order * dim; 45851a454edSToby Isaac } else { 45951a454edSToby Isaac maxOrder = order; 46051a454edSToby Isaac } 46151a454edSToby Isaac if (isConstant) *isConstant = (maxOrder < 1) ? PETSC_TRUE : PETSC_FALSE; 46251a454edSToby Isaac if (isAffine) *isAffine = (maxOrder < 2) ? PETSC_TRUE : PETSC_FALSE; 46351a454edSToby Isaac if (isQuadratic) *isQuadratic = (maxOrder < 3) ? PETSC_TRUE : PETSC_FALSE; 46451a454edSToby Isaac } 46551a454edSToby Isaac PetscFunctionReturn(0); 46651a454edSToby Isaac } 46751a454edSToby Isaac 46851a454edSToby Isaac static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad) 46951a454edSToby Isaac { 47051a454edSToby Isaac PetscInt h, dim, imax, imin; 47151a454edSToby Isaac DM dm; 47251a454edSToby Isaac DMField_DS *dsfield; 47351a454edSToby Isaac PetscObject disc; 47451a454edSToby Isaac PetscFE fe; 47551a454edSToby Isaac PetscClassId id; 47651a454edSToby Isaac PetscErrorCode ierr; 47751a454edSToby Isaac 47851a454edSToby Isaac 47951a454edSToby Isaac PetscFunctionBegin; 48051a454edSToby Isaac dm = field->dm; 48151a454edSToby Isaac dsfield = (DMField_DS *) field->data; 48251a454edSToby Isaac ierr = ISGetMinMax(pointIS,&imax,&imin);CHKERRQ(ierr); 48351a454edSToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 48451a454edSToby Isaac for (h = 0; h <= dim; h++) { 48551a454edSToby Isaac PetscInt hStart, hEnd; 48651a454edSToby Isaac 48751a454edSToby Isaac ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 48851a454edSToby Isaac if (imin >= hStart && imax < hEnd) break; 48951a454edSToby Isaac } 49051a454edSToby Isaac *quad = NULL; 491f99c8401SToby Isaac if (h < dsfield->height) { 49251a454edSToby Isaac ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr); 49351a454edSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 49451a454edSToby Isaac if (id != PETSCFE_CLASSID) PetscFunctionReturn(0); 49551a454edSToby Isaac fe = (PetscFE) disc; 49651a454edSToby Isaac ierr = PetscFEGetQuadrature(fe,quad);CHKERRQ(ierr); 49751a454edSToby Isaac ierr = PetscObjectReference((PetscObject)*quad);CHKERRQ(ierr); 49851a454edSToby Isaac } 49951a454edSToby Isaac PetscFunctionReturn(0); 50051a454edSToby Isaac } 50151a454edSToby Isaac 5020145028aSToby Isaac static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom) 5030145028aSToby Isaac { 5040145028aSToby Isaac const PetscInt *points; 5050145028aSToby Isaac PetscInt p, dim, dE, numFaces, Nq; 5060145028aSToby Isaac PetscBool affineCells; 5070145028aSToby Isaac DMLabel depthLabel; 5080145028aSToby Isaac IS cellIS; 5090145028aSToby Isaac DM dm = field->dm; 5100145028aSToby Isaac PetscErrorCode ierr; 5110145028aSToby Isaac 5120145028aSToby Isaac PetscFunctionBegin; 5130145028aSToby Isaac dim = geom->dim; 5140145028aSToby Isaac dE = geom->dimEmbed; 5150145028aSToby Isaac ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 5160145028aSToby Isaac ierr = DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS);CHKERRQ(ierr); 5170145028aSToby Isaac ierr = DMFieldGetFEInvariance(field,cellIS,NULL,&affineCells,NULL);CHKERRQ(ierr); 5180145028aSToby Isaac ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 5190145028aSToby Isaac numFaces = geom->numCells; 5200145028aSToby Isaac Nq = geom->numPoints; 5210145028aSToby Isaac if (affineCells) { 5220145028aSToby Isaac PetscInt numCells, offset, *cells; 5230145028aSToby Isaac PetscFEGeom *cellGeom; 5240145028aSToby Isaac IS suppIS; 5250145028aSToby Isaac PetscQuadrature cellQuad = NULL; 5260145028aSToby Isaac 5270145028aSToby Isaac ierr = DMFieldCreateDefaultQuadrature(field,cellIS,&cellQuad);CHKERRQ(ierr); 5280145028aSToby Isaac for (p = 0, numCells = 0; p < numFaces; p++) { 5290145028aSToby Isaac PetscInt point = points[p]; 5300145028aSToby Isaac PetscInt numSupp, numChildren; 5310145028aSToby Isaac 5320145028aSToby Isaac ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL); CHKERRQ(ierr); 5330145028aSToby Isaac if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 5340145028aSToby Isaac ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 5350145028aSToby Isaac if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp); 5360145028aSToby Isaac numCells += numSupp; 5370145028aSToby Isaac } 5380145028aSToby Isaac ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr); 5390145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 5400145028aSToby Isaac PetscInt point = points[p]; 5410145028aSToby Isaac PetscInt numSupp, s; 5420145028aSToby Isaac const PetscInt *supp; 5430145028aSToby Isaac 5440145028aSToby Isaac ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 5450145028aSToby Isaac ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 5460145028aSToby Isaac for (s = 0; s < numSupp; s++, offset++) { 5470145028aSToby Isaac cells[offset] = supp[s]; 5480145028aSToby Isaac } 5490145028aSToby Isaac } 5500145028aSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr); 5510145028aSToby Isaac ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr); 5520145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 5530145028aSToby Isaac PetscInt point = points[p]; 5540145028aSToby Isaac PetscInt numSupp, s, q; 5550145028aSToby Isaac const PetscInt *supp; 5560145028aSToby Isaac 5570145028aSToby Isaac ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 5580145028aSToby Isaac ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 5590145028aSToby Isaac for (s = 0; s < numSupp; s++, offset++) { 5600145028aSToby Isaac for (q = 0; q < Nq * dE * dE; q++) { 5610145028aSToby Isaac geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 5620145028aSToby Isaac } 5630145028aSToby Isaac } 5640145028aSToby Isaac } 5650145028aSToby Isaac ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr); 5660145028aSToby Isaac ierr = ISDestroy(&suppIS);CHKERRQ(ierr); 5670145028aSToby Isaac ierr = PetscFree(cells);CHKERRQ(ierr); 5680145028aSToby Isaac ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr); 5690145028aSToby Isaac } else { 5700145028aSToby Isaac PetscObject faceDisc, cellDisc; 5710145028aSToby Isaac PetscClassId faceId, cellId; 5720145028aSToby Isaac PetscDualSpace dsp; 5730145028aSToby Isaac DM K; 5740145028aSToby Isaac PetscInt (*co)[2][3]; 5750145028aSToby Isaac PetscInt coneSize; 5760145028aSToby Isaac PetscInt **counts; 5770145028aSToby Isaac PetscInt f, i, o, q, s; 5780145028aSToby Isaac const PetscInt *coneK; 5790145028aSToby Isaac PetscInt minOrient, maxOrient, numOrient; 5800145028aSToby Isaac PetscInt *orients; 5810145028aSToby Isaac PetscReal **orientPoints; 5820145028aSToby Isaac PetscReal *cellPoints; 5830145028aSToby Isaac PetscReal *dummyWeights; 5840145028aSToby Isaac PetscQuadrature cellQuad = NULL; 5850145028aSToby Isaac 5860145028aSToby Isaac ierr = DMFieldDSGetHeightDisc(field, 1, &faceDisc);CHKERRQ(ierr); 5870145028aSToby Isaac ierr = DMFieldDSGetHeightDisc(field, 0, &cellDisc);CHKERRQ(ierr); 5880145028aSToby Isaac ierr = PetscObjectGetClassId(faceDisc,&faceId);CHKERRQ(ierr); 5890145028aSToby Isaac ierr = PetscObjectGetClassId(cellDisc,&cellId);CHKERRQ(ierr); 5900145028aSToby Isaac if (faceId != PETSCFE_CLASSID || cellId != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported\n"); 5910145028aSToby Isaac ierr = PetscFEGetDualSpace((PetscFE)cellDisc, &dsp);CHKERRQ(ierr); 5920145028aSToby Isaac ierr = PetscDualSpaceGetDM(dsp, &K); CHKERRQ(ierr); 5930145028aSToby Isaac ierr = DMPlexGetConeSize(K,0,&coneSize);CHKERRQ(ierr); 5940145028aSToby Isaac ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr); 5950145028aSToby Isaac ierr = PetscMalloc4(numFaces,&co,dE*Nq,&cellPoints,coneSize,&counts,Nq,&dummyWeights);CHKERRQ(ierr); 5960145028aSToby Isaac ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)field), &cellQuad);CHKERRQ(ierr); 5970145028aSToby Isaac ierr = PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights);CHKERRQ(ierr); 5980145028aSToby Isaac minOrient = PETSC_MAX_INT; 5990145028aSToby Isaac maxOrient = PETSC_MIN_INT; 6000145028aSToby Isaac for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */ 6010145028aSToby Isaac PetscInt point = points[p]; 6020145028aSToby Isaac PetscInt numSupp, numChildren; 6030145028aSToby Isaac const PetscInt *supp; 6040145028aSToby Isaac 6050145028aSToby Isaac ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL); CHKERRQ(ierr); 6060145028aSToby Isaac if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 6070145028aSToby Isaac ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr); 6080145028aSToby Isaac if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp); 6090145028aSToby Isaac ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr); 6100145028aSToby Isaac for (s = 0; s < numSupp; s++) { 6110145028aSToby Isaac PetscInt cell = supp[s]; 6120145028aSToby Isaac PetscInt numCone; 6130145028aSToby Isaac const PetscInt *cone, *orient; 6140145028aSToby Isaac 6150145028aSToby Isaac ierr = DMPlexGetConeSize(dm, cell, &numCone);CHKERRQ(ierr); 6160145028aSToby Isaac if (numCone != coneSize) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element"); 6170145028aSToby Isaac ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 6180145028aSToby Isaac ierr = DMPlexGetConeOrientation(dm, cell, &orient);CHKERRQ(ierr); 6190145028aSToby Isaac for (f = 0; f < coneSize; f++) { 6200145028aSToby Isaac if (cone[f] == point) break; 6210145028aSToby Isaac } 6220145028aSToby Isaac co[p][s][0] = f; 6230145028aSToby Isaac co[p][s][1] = orient[f]; 6240145028aSToby Isaac co[p][s][2] = cell; 6250145028aSToby Isaac minOrient = PetscMin(minOrient, orient[f]); 6260145028aSToby Isaac maxOrient = PetscMin(maxOrient, orient[f]); 6270145028aSToby Isaac } 6280145028aSToby Isaac for (; s < 2; s++) { 6290145028aSToby Isaac co[p][s][0] = -1; 6300145028aSToby Isaac co[p][s][1] = -1; 6310145028aSToby Isaac co[p][s][2] = -1; 6320145028aSToby Isaac } 6330145028aSToby Isaac } 6340145028aSToby Isaac numOrient = maxOrient + 1 - minOrient; 6350145028aSToby Isaac ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr); 6360145028aSToby Isaac /* count all (face,orientation) doubles that appear */ 6370145028aSToby Isaac ierr = PetscCalloc2(numOrient,&orients,numOrient,&orientPoints);CHKERRQ(ierr); 6380145028aSToby Isaac for (f = 0; f < coneSize; f++) {ierr = PetscCalloc1(numOrient, &counts[f]);CHKERRQ(ierr);} 6390145028aSToby Isaac for (p = 0; p < numFaces; p++) { 6400145028aSToby Isaac for (s = 0; s < 2; s++) { 6410145028aSToby Isaac if (co[p][s][0] >= 0) { 6420145028aSToby Isaac counts[co[p][s][0]][co[p][s][1] - minOrient]++; 6430145028aSToby Isaac orients[co[p][s][1] - minOrient]++; 6440145028aSToby Isaac } 6450145028aSToby Isaac } 6460145028aSToby Isaac } 6470145028aSToby Isaac for (o = 0; o < numOrient; o++) { 6480145028aSToby Isaac if (orients[o]) { 6490145028aSToby Isaac PetscInt orient = o + minOrient; 6500145028aSToby Isaac PetscInt q; 6510145028aSToby Isaac 6520145028aSToby Isaac ierr = PetscMalloc1(Nq * dim, &orientPoints[o]);CHKERRQ(ierr); 6530145028aSToby Isaac /* rotate the quadrature points appropriately */ 6540145028aSToby Isaac switch (dim) { 6550145028aSToby Isaac case 0: 6560145028aSToby Isaac break; 6570145028aSToby Isaac case 1: 6580145028aSToby Isaac if (orient == -2 || orient == 1) { 6590145028aSToby Isaac for (q = 0; q < Nq; q++) { 6600145028aSToby Isaac orientPoints[o][q] = -geom->xi[q]; 6610145028aSToby Isaac } 6620145028aSToby Isaac } else { 6630145028aSToby Isaac for (q = 0; q < Nq; q++) { 6640145028aSToby Isaac orientPoints[o][q] = geom->xi[q]; 6650145028aSToby Isaac } 6660145028aSToby Isaac } 6670145028aSToby Isaac break; 6680145028aSToby Isaac case 2: 6690145028aSToby Isaac switch (coneSize) { 6700145028aSToby Isaac case 3: 6710145028aSToby Isaac for (q = 0; q < Nq; q++) { 6720145028aSToby Isaac PetscReal lambda[3]; 6730145028aSToby Isaac PetscReal lambdao[3]; 6740145028aSToby Isaac 6750145028aSToby Isaac /* convert to barycentric */ 6760145028aSToby Isaac lambda[0] = - (geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.; 6770145028aSToby Isaac lambda[1] = (geom->xi[2 * q] + 1.) / 2.; 6780145028aSToby Isaac lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.; 6790145028aSToby Isaac if (orient >= 0) { 6800145028aSToby Isaac for (i = 0; i < 3; i++) { 6810145028aSToby Isaac lambdao[i] = lambda[(orient + i) % 3]; 6820145028aSToby Isaac } 6830145028aSToby Isaac } else { 6840145028aSToby Isaac for (i = 0; i < 3; i++) { 6850145028aSToby Isaac lambdao[i] = lambda[(-(orient + i) + 3) % 3]; 6860145028aSToby Isaac } 6870145028aSToby Isaac } 6880145028aSToby Isaac /* convert to coordinates */ 6890145028aSToby Isaac orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1]; 6900145028aSToby Isaac orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2]; 6910145028aSToby Isaac } 6920145028aSToby Isaac break; 6930145028aSToby Isaac case 4: 6940145028aSToby Isaac for (q = 0; q < Nq; q++) { 6950145028aSToby Isaac PetscReal xi[2], xio[2]; 6960145028aSToby Isaac PetscInt oabs = (orient >= 0) ? orient : -(orient + 1); 6970145028aSToby Isaac 6980145028aSToby Isaac xi[0] = geom->xi[2 * q]; 6990145028aSToby Isaac xi[1] = geom->xi[2 * q + 1]; 7000145028aSToby Isaac switch (oabs) { 7010145028aSToby Isaac case 0: 7020145028aSToby Isaac xio[0] = xi[0]; 7030145028aSToby Isaac xio[1] = xi[1]; 7040145028aSToby Isaac break; 7050145028aSToby Isaac case 1: 7060145028aSToby Isaac xio[0] = xi[1]; 7070145028aSToby Isaac xio[1] = -xi[0]; 7080145028aSToby Isaac break; 7090145028aSToby Isaac case 2: 7100145028aSToby Isaac xio[0] = -xi[0]; 7110145028aSToby Isaac xio[1] = -xi[1]; 7120145028aSToby Isaac case 3: 7130145028aSToby Isaac xio[0] = -xi[1]; 7140145028aSToby Isaac xio[1] = xi[0]; 7150145028aSToby Isaac } 7160145028aSToby Isaac if (orient < 0) { 7170145028aSToby Isaac xio[0] = -xio[0]; 7180145028aSToby Isaac } 7190145028aSToby Isaac orientPoints[o][2 * q + 0] = xio[0]; 7200145028aSToby Isaac orientPoints[o][2 * q + 1] = xio[1]; 7210145028aSToby Isaac } 7220145028aSToby Isaac break; 7230145028aSToby Isaac default: 7240145028aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not implemented yet\n"); 7250145028aSToby Isaac } 7260145028aSToby Isaac default: 7270145028aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not implemented yet\n"); 7280145028aSToby Isaac } 7290145028aSToby Isaac } 7300145028aSToby Isaac } 7310145028aSToby Isaac for (f = 0; f < coneSize; f++) { 7320145028aSToby Isaac PetscInt face = coneK[f]; 7330145028aSToby Isaac PetscScalar v0[3]; 7340145028aSToby Isaac PetscScalar J[9]; 7350145028aSToby Isaac PetscInt numCells, offset; 7360145028aSToby Isaac PetscInt *cells; 7370145028aSToby Isaac IS suppIS; 7380145028aSToby Isaac 7390145028aSToby Isaac ierr = DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, NULL);CHKERRQ(ierr); 7400145028aSToby Isaac for (o = 0; o <= numOrient; o++) { 7410145028aSToby Isaac PetscFEGeom *cellGeom; 7420145028aSToby Isaac 7430145028aSToby Isaac if (!counts[f][o]) continue; 7440145028aSToby Isaac /* If this (face,orientation) double appears, 7450145028aSToby Isaac * convert the face quadrature points into volume quadrature points */ 7460145028aSToby Isaac for (q = 0; q < Nq; q++) { 7470145028aSToby Isaac PetscReal xi0[3] = {-1., -1., -1.}; 7480145028aSToby Isaac 7490145028aSToby Isaac CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]); 7500145028aSToby Isaac } 7510145028aSToby Isaac for (p = 0, numCells = 0; p < numFaces; p++) { 7520145028aSToby Isaac for (s = 0; s < 2; s++) { 7530145028aSToby Isaac if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++; 7540145028aSToby Isaac } 7550145028aSToby Isaac } 7560145028aSToby Isaac ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr); 7570145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 7580145028aSToby Isaac for (s = 0; s < 2; s++) { 7590145028aSToby Isaac if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 7600145028aSToby Isaac cells[offset++] = co[p][s][2]; 7610145028aSToby Isaac } 7620145028aSToby Isaac } 7630145028aSToby Isaac } 7640145028aSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr); 7650145028aSToby Isaac ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr); 7660145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 7670145028aSToby Isaac for (s = 0; s < 2; s++) { 7680145028aSToby Isaac if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 7690145028aSToby Isaac for (q = 0; q < Nq * dE * dE; q++) { 7700145028aSToby Isaac geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 7710145028aSToby Isaac } 7720145028aSToby Isaac offset++; 7730145028aSToby Isaac } 7740145028aSToby Isaac } 7750145028aSToby Isaac } 7760145028aSToby Isaac ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr); 7770145028aSToby Isaac ierr = ISDestroy(&suppIS);CHKERRQ(ierr); 7780145028aSToby Isaac ierr = PetscFree(cells);CHKERRQ(ierr); 7790145028aSToby Isaac } 7800145028aSToby Isaac } 7810145028aSToby Isaac for (o = 0; o < numOrient; o++) { 7820145028aSToby Isaac if (orients[o]) { 7830145028aSToby Isaac ierr = PetscFree(orientPoints[o]);CHKERRQ(ierr); 7840145028aSToby Isaac } 7850145028aSToby Isaac } 7860145028aSToby Isaac ierr = PetscFree2(orients,orientPoints);CHKERRQ(ierr); 7870145028aSToby Isaac ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr); 7880145028aSToby Isaac ierr = PetscFree4(co,cellPoints,counts,dummyWeights);CHKERRQ(ierr); 7890145028aSToby Isaac } 7900145028aSToby Isaac ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 7910145028aSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 7920145028aSToby Isaac PetscFunctionReturn(0); 7930145028aSToby Isaac } 7940145028aSToby Isaac 79551a454edSToby Isaac static PetscErrorCode DMFieldInitialize_DS(DMField field) 79651a454edSToby Isaac { 79751a454edSToby Isaac PetscFunctionBegin; 79851a454edSToby Isaac field->ops->destroy = DMFieldDestroy_DS; 79951a454edSToby Isaac field->ops->evaluate = DMFieldEvaluate_DS; 80051a454edSToby Isaac field->ops->evaluateFE = DMFieldEvaluateFE_DS; 80151a454edSToby Isaac field->ops->getFEInvariance = DMFieldGetFEInvariance_DS; 80251a454edSToby Isaac field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS; 80351a454edSToby Isaac field->ops->view = DMFieldView_DS; 8040145028aSToby Isaac field->ops->computeFaceData = DMFieldComputeFaceData_DS; 80551a454edSToby Isaac PetscFunctionReturn(0); 80651a454edSToby Isaac } 80751a454edSToby Isaac 80851a454edSToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field) 80951a454edSToby Isaac { 81051a454edSToby Isaac DMField_DS *dsfield; 81151a454edSToby Isaac PetscErrorCode ierr; 81251a454edSToby Isaac 81351a454edSToby Isaac PetscFunctionBegin; 81451a454edSToby Isaac ierr = PetscNewLog(field,&dsfield);CHKERRQ(ierr); 81551a454edSToby Isaac field->data = dsfield; 81651a454edSToby Isaac ierr = DMFieldInitialize_DS(field);CHKERRQ(ierr); 81751a454edSToby Isaac PetscFunctionReturn(0); 81851a454edSToby Isaac } 81951a454edSToby Isaac 82051a454edSToby Isaac PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec,DMField *field) 82151a454edSToby Isaac { 82251a454edSToby Isaac DMField b; 82351a454edSToby Isaac DMField_DS *dsfield; 82451a454edSToby Isaac PetscObject disc = NULL; 82551a454edSToby Isaac PetscBool isContainer = PETSC_FALSE; 82651a454edSToby Isaac PetscClassId id = -1; 82751a454edSToby Isaac PetscInt numComponents = -1, dsNumFields; 82851a454edSToby Isaac PetscSection section; 82951a454edSToby Isaac PetscErrorCode ierr; 83051a454edSToby Isaac 83151a454edSToby Isaac PetscFunctionBegin; 83251a454edSToby Isaac ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 83351a454edSToby Isaac ierr = PetscSectionGetFieldComponents(section,fieldNum,&numComponents);CHKERRQ(ierr); 83451a454edSToby Isaac ierr = DMGetNumFields(dm,&dsNumFields);CHKERRQ(ierr); 83551a454edSToby Isaac if (dsNumFields) {ierr = DMGetField(dm,fieldNum,&disc);CHKERRQ(ierr);} 83651a454edSToby Isaac if (disc) { 83751a454edSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 83851a454edSToby Isaac isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE; 83951a454edSToby Isaac } 84051a454edSToby Isaac if (!disc || isContainer) { 84151a454edSToby Isaac MPI_Comm comm = PetscObjectComm((PetscObject) dm); 84251a454edSToby Isaac PetscInt cStart, cEnd, dim; 84351a454edSToby Isaac PetscInt localConeSize = 0, coneSize; 84451a454edSToby Isaac PetscFE fe; 84551a454edSToby Isaac PetscDualSpace Q; 84651a454edSToby Isaac PetscSpace P; 84751a454edSToby Isaac DM K; 84851a454edSToby Isaac PetscQuadrature quad, fquad; 84951a454edSToby Isaac PetscBool isSimplex; 85051a454edSToby Isaac 85151a454edSToby Isaac ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 85251a454edSToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 85351a454edSToby Isaac if (cEnd > cStart) { 85451a454edSToby Isaac ierr = DMPlexGetConeSize(dm, cStart, &localConeSize);CHKERRQ(ierr); 85551a454edSToby Isaac } 85651a454edSToby Isaac ierr = MPI_Allreduce(&localConeSize,&coneSize,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 85751a454edSToby Isaac isSimplex = (coneSize == (dim + 1)) ? PETSC_TRUE : PETSC_FALSE; 85851a454edSToby Isaac ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 85951a454edSToby Isaac ierr = PetscSpaceSetOrder(P, 1);CHKERRQ(ierr); 86051a454edSToby Isaac ierr = PetscSpaceSetNumComponents(P, numComponents);CHKERRQ(ierr); 86151a454edSToby Isaac ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 86251a454edSToby Isaac ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 86351a454edSToby Isaac ierr = PetscSpacePolynomialSetTensor(P, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr); 86451a454edSToby Isaac ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 86551a454edSToby Isaac ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 86651a454edSToby Isaac ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 86751a454edSToby Isaac ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 86851a454edSToby Isaac ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 86951a454edSToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 87051a454edSToby Isaac ierr = PetscDualSpaceSetNumComponents(Q, numComponents);CHKERRQ(ierr); 87151a454edSToby Isaac ierr = PetscDualSpaceSetOrder(Q, 1);CHKERRQ(ierr); 87251a454edSToby Isaac ierr = PetscDualSpaceLagrangeSetTensor(Q, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr); 87351a454edSToby Isaac ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 87451a454edSToby Isaac ierr = PetscFECreate(comm, &fe);CHKERRQ(ierr); 87551a454edSToby Isaac ierr = PetscFESetType(fe,PETSCFEBASIC);CHKERRQ(ierr); 87651a454edSToby Isaac ierr = PetscFESetBasisSpace(fe, P);CHKERRQ(ierr); 87751a454edSToby Isaac ierr = PetscFESetDualSpace(fe, Q);CHKERRQ(ierr); 87851a454edSToby Isaac ierr = PetscFESetNumComponents(fe, numComponents);CHKERRQ(ierr); 87951a454edSToby Isaac ierr = PetscFESetUp(fe);CHKERRQ(ierr); 88051a454edSToby Isaac ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 88151a454edSToby Isaac ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 88251a454edSToby Isaac if (isSimplex) { 88351a454edSToby Isaac ierr = PetscDTGaussJacobiQuadrature(dim, 1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr); 88451a454edSToby Isaac ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr); 88551a454edSToby Isaac } 88651a454edSToby Isaac else { 88751a454edSToby Isaac ierr = PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr); 88851a454edSToby Isaac ierr = PetscDTGaussTensorQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr); 88951a454edSToby Isaac } 89051a454edSToby Isaac ierr = PetscFESetQuadrature(fe, quad);CHKERRQ(ierr); 89151a454edSToby Isaac ierr = PetscFESetFaceQuadrature(fe, fquad);CHKERRQ(ierr); 89251a454edSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 89351a454edSToby Isaac ierr = PetscQuadratureDestroy(&fquad);CHKERRQ(ierr); 89451a454edSToby Isaac disc = (PetscObject) fe; 89551a454edSToby Isaac } else { 89651a454edSToby Isaac ierr = PetscObjectReference(disc);CHKERRQ(ierr); 89751a454edSToby Isaac } 89851a454edSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 89951a454edSToby Isaac if (id == PETSCFE_CLASSID) { 90051a454edSToby Isaac PetscFE fe = (PetscFE) disc; 90151a454edSToby Isaac 90251a454edSToby Isaac ierr = PetscFEGetNumComponents(fe,&numComponents);CHKERRQ(ierr); 90351a454edSToby Isaac } else {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented");} 90451a454edSToby Isaac ierr = DMFieldCreate(dm,numComponents,DMFIELD_VERTEX,&b);CHKERRQ(ierr); 90551a454edSToby Isaac ierr = DMFieldSetType(b,DMFIELDDS);CHKERRQ(ierr); 90651a454edSToby Isaac dsfield = (DMField_DS *) b->data; 90751a454edSToby Isaac dsfield->fieldNum = fieldNum; 90851a454edSToby Isaac ierr = DMGetDimension(dm,&dsfield->height);CHKERRQ(ierr); 909f99c8401SToby Isaac dsfield->height++; 91051a454edSToby Isaac ierr = PetscCalloc1(dsfield->height,&dsfield->disc);CHKERRQ(ierr); 91151a454edSToby Isaac dsfield->disc[0] = disc; 91251a454edSToby Isaac ierr = PetscObjectReference((PetscObject)vec);CHKERRQ(ierr); 91351a454edSToby Isaac dsfield->vec = vec; 91451a454edSToby Isaac *field = b; 91551a454edSToby Isaac 91651a454edSToby Isaac PetscFunctionReturn(0); 91751a454edSToby Isaac } 918