1*a8db3e61SBlaise Bourdin static char help[] = "Test FEM layout and GlobalToNaturalSF\n\n"; 2*a8db3e61SBlaise Bourdin 3*a8db3e61SBlaise Bourdin /* 4*a8db3e61SBlaise Bourdin In order to see the vectors which are being tested, use 5*a8db3e61SBlaise Bourdin 6*a8db3e61SBlaise Bourdin -ua_vec_view -s_vec_view 7*a8db3e61SBlaise Bourdin */ 8*a8db3e61SBlaise Bourdin 9*a8db3e61SBlaise Bourdin #include <petsc.h> 10*a8db3e61SBlaise Bourdin #include <exodusII.h> 11*a8db3e61SBlaise Bourdin 12*a8db3e61SBlaise Bourdin #include <petsc/private/dmpleximpl.h> 13*a8db3e61SBlaise Bourdin 14*a8db3e61SBlaise Bourdin int main(int argc, char **argv) { 15*a8db3e61SBlaise Bourdin DM dm, pdm, dmU, dmA, dmS, dmUA, dmUA2, *dmList; 16*a8db3e61SBlaise Bourdin DM seqdmU, seqdmA, seqdmS, seqdmUA, seqdmUA2, *seqdmList; 17*a8db3e61SBlaise Bourdin Vec X, U, A, S, UA, UA2; 18*a8db3e61SBlaise Bourdin Vec seqX, seqU, seqA, seqS, seqUA, seqUA2; 19*a8db3e61SBlaise Bourdin IS isU, isA, isS, isUA; 20*a8db3e61SBlaise Bourdin IS seqisU, seqisA, seqisS, seqisUA; 21*a8db3e61SBlaise Bourdin PetscSection section; 22*a8db3e61SBlaise Bourdin const PetscInt fieldU = 0; 23*a8db3e61SBlaise Bourdin const PetscInt fieldA = 2; 24*a8db3e61SBlaise Bourdin const PetscInt fieldS = 1; 25*a8db3e61SBlaise Bourdin const PetscInt fieldUA[2] = {0, 2}; 26*a8db3e61SBlaise Bourdin char ifilename[PETSC_MAX_PATH_LEN]; 27*a8db3e61SBlaise Bourdin IS csIS; 28*a8db3e61SBlaise Bourdin const PetscInt *csID; 29*a8db3e61SBlaise Bourdin PetscInt *pStartDepth, *pEndDepth; 30*a8db3e61SBlaise Bourdin PetscInt order = 1; 31*a8db3e61SBlaise Bourdin PetscInt sdim, d, pStart, pEnd, p, numCS, set; 32*a8db3e61SBlaise Bourdin PetscMPIInt rank, size; 33*a8db3e61SBlaise Bourdin PetscSF migrationSF; 34*a8db3e61SBlaise Bourdin 35*a8db3e61SBlaise Bourdin PetscCall(PetscInitialize(&argc, &argv,NULL, help)); 36*a8db3e61SBlaise Bourdin PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 37*a8db3e61SBlaise Bourdin PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 38*a8db3e61SBlaise Bourdin PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex64"); 39*a8db3e61SBlaise Bourdin PetscCall(PetscOptionsString("-i", "Filename to read", "ex64", ifilename, ifilename, sizeof(ifilename), NULL)); 40*a8db3e61SBlaise Bourdin PetscOptionsEnd(); 41*a8db3e61SBlaise Bourdin 42*a8db3e61SBlaise Bourdin /* Read the mesh from a file in any supported format */ 43*a8db3e61SBlaise Bourdin PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, "ex64_plex", PETSC_TRUE, &dm)); 44*a8db3e61SBlaise Bourdin PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 45*a8db3e61SBlaise Bourdin PetscCall(DMSetFromOptions(dm)); 46*a8db3e61SBlaise Bourdin PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 47*a8db3e61SBlaise Bourdin PetscCall(DMGetDimension(dm, &sdim)); 48*a8db3e61SBlaise Bourdin 49*a8db3e61SBlaise Bourdin /* Create the main section containing all fields */ 50*a8db3e61SBlaise Bourdin PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion)); 51*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetNumFields(section, 3)); 52*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldName(section, fieldU, "U")); 53*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldName(section, fieldA, "Alpha")); 54*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldName(section, fieldS, "Sigma")); 55*a8db3e61SBlaise Bourdin PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 56*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetChart(section, pStart, pEnd)); 57*a8db3e61SBlaise Bourdin PetscCall(PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth)); 58*a8db3e61SBlaise Bourdin for (d = 0; d <= sdim; ++d) {PetscCall(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]));} 59*a8db3e61SBlaise Bourdin /* Vector field U, Scalar field Alpha, Tensor field Sigma */ 60*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldComponents(section, fieldU, sdim)); 61*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldComponents(section, fieldA, 1)); 62*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2)); 63*a8db3e61SBlaise Bourdin 64*a8db3e61SBlaise Bourdin /* Going through cell sets then cells, and setting up storage for the sections */ 65*a8db3e61SBlaise Bourdin PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS)); 66*a8db3e61SBlaise Bourdin PetscCall(DMGetLabelIdIS(dm, "Cell Sets", &csIS)); 67*a8db3e61SBlaise Bourdin if (csIS) {PetscCall(ISGetIndices(csIS, &csID));} 68*a8db3e61SBlaise Bourdin for (set = 0; set < numCS; set++) { 69*a8db3e61SBlaise Bourdin IS cellIS; 70*a8db3e61SBlaise Bourdin const PetscInt *cellID; 71*a8db3e61SBlaise Bourdin PetscInt numCells, cell, closureSize, *closureA = NULL; 72*a8db3e61SBlaise Bourdin 73*a8db3e61SBlaise Bourdin PetscCall(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells)); 74*a8db3e61SBlaise Bourdin PetscCall(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS)); 75*a8db3e61SBlaise Bourdin if (numCells > 0) { 76*a8db3e61SBlaise Bourdin /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */ 77*a8db3e61SBlaise Bourdin PetscInt dofUP1Tri[] = {2, 0, 0}; 78*a8db3e61SBlaise Bourdin PetscInt dofAP1Tri[] = {1, 0, 0}; 79*a8db3e61SBlaise Bourdin PetscInt dofUP2Tri[] = {2, 2, 0}; 80*a8db3e61SBlaise Bourdin PetscInt dofAP2Tri[] = {1, 1, 0}; 81*a8db3e61SBlaise Bourdin PetscInt dofUP1Quad[] = {2, 0, 0}; 82*a8db3e61SBlaise Bourdin PetscInt dofAP1Quad[] = {1, 0, 0}; 83*a8db3e61SBlaise Bourdin PetscInt dofUP2Quad[] = {2, 2, 2}; 84*a8db3e61SBlaise Bourdin PetscInt dofAP2Quad[] = {1, 1, 1}; 85*a8db3e61SBlaise Bourdin PetscInt dofS2D[] = {0, 0, 3}; 86*a8db3e61SBlaise Bourdin PetscInt dofUP1Tet[] = {3, 0, 0, 0}; 87*a8db3e61SBlaise Bourdin PetscInt dofAP1Tet[] = {1, 0, 0, 0}; 88*a8db3e61SBlaise Bourdin PetscInt dofUP2Tet[] = {3, 3, 0, 0}; 89*a8db3e61SBlaise Bourdin PetscInt dofAP2Tet[] = {1, 1, 0, 0}; 90*a8db3e61SBlaise Bourdin PetscInt dofUP1Hex[] = {3, 0, 0, 0}; 91*a8db3e61SBlaise Bourdin PetscInt dofAP1Hex[] = {1, 0, 0, 0}; 92*a8db3e61SBlaise Bourdin PetscInt dofUP2Hex[] = {3, 3, 3, 3}; 93*a8db3e61SBlaise Bourdin PetscInt dofAP2Hex[] = {1, 1, 1, 1}; 94*a8db3e61SBlaise Bourdin PetscInt dofS3D[] = {0, 0, 0, 6}; 95*a8db3e61SBlaise Bourdin PetscInt *dofU, *dofA, *dofS; 96*a8db3e61SBlaise Bourdin 97*a8db3e61SBlaise Bourdin switch (sdim) { 98*a8db3e61SBlaise Bourdin case 2: dofS = dofS2D;break; 99*a8db3e61SBlaise Bourdin case 3: dofS = dofS3D;break; 100*a8db3e61SBlaise Bourdin default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim); 101*a8db3e61SBlaise Bourdin } 102*a8db3e61SBlaise Bourdin 103*a8db3e61SBlaise Bourdin /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 104*a8db3e61SBlaise Bourdin It will not be enough to identify more exotic elements like pyramid or prisms... */ 105*a8db3e61SBlaise Bourdin PetscCall(ISGetIndices(cellIS, &cellID)); 106*a8db3e61SBlaise Bourdin PetscCall(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA)); 107*a8db3e61SBlaise Bourdin switch (closureSize) { 108*a8db3e61SBlaise Bourdin case 7: /* Tri */ 109*a8db3e61SBlaise Bourdin if (order == 1) { 110*a8db3e61SBlaise Bourdin dofU = dofUP1Tri; 111*a8db3e61SBlaise Bourdin dofA = dofAP1Tri; 112*a8db3e61SBlaise Bourdin } else { 113*a8db3e61SBlaise Bourdin dofU = dofUP2Tri; 114*a8db3e61SBlaise Bourdin dofA = dofAP2Tri; 115*a8db3e61SBlaise Bourdin } 116*a8db3e61SBlaise Bourdin break; 117*a8db3e61SBlaise Bourdin case 9: /* Quad */ 118*a8db3e61SBlaise Bourdin if (order == 1) { 119*a8db3e61SBlaise Bourdin dofU = dofUP1Quad; 120*a8db3e61SBlaise Bourdin dofA = dofAP1Quad; 121*a8db3e61SBlaise Bourdin } else { 122*a8db3e61SBlaise Bourdin dofU = dofUP2Quad; 123*a8db3e61SBlaise Bourdin dofA = dofAP2Quad; 124*a8db3e61SBlaise Bourdin } 125*a8db3e61SBlaise Bourdin break; 126*a8db3e61SBlaise Bourdin case 15: /* Tet */ 127*a8db3e61SBlaise Bourdin if (order == 1) { 128*a8db3e61SBlaise Bourdin dofU = dofUP1Tet; 129*a8db3e61SBlaise Bourdin dofA = dofAP1Tet; 130*a8db3e61SBlaise Bourdin } else { 131*a8db3e61SBlaise Bourdin dofU = dofUP2Tet; 132*a8db3e61SBlaise Bourdin dofA = dofAP2Tet; 133*a8db3e61SBlaise Bourdin } 134*a8db3e61SBlaise Bourdin break; 135*a8db3e61SBlaise Bourdin case 27: /* Hex */ 136*a8db3e61SBlaise Bourdin if (order == 1) { 137*a8db3e61SBlaise Bourdin dofU = dofUP1Hex; 138*a8db3e61SBlaise Bourdin dofA = dofAP1Hex; 139*a8db3e61SBlaise Bourdin } else { 140*a8db3e61SBlaise Bourdin dofU = dofUP2Hex; 141*a8db3e61SBlaise Bourdin dofA = dofAP2Hex; 142*a8db3e61SBlaise Bourdin } 143*a8db3e61SBlaise Bourdin break; 144*a8db3e61SBlaise Bourdin default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %" PetscInt_FMT, closureSize); 145*a8db3e61SBlaise Bourdin } 146*a8db3e61SBlaise Bourdin PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA)); 147*a8db3e61SBlaise Bourdin 148*a8db3e61SBlaise Bourdin for (cell = 0; cell < numCells; cell++) { 149*a8db3e61SBlaise Bourdin PetscInt *closure = NULL; 150*a8db3e61SBlaise Bourdin 151*a8db3e61SBlaise Bourdin PetscCall(DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure)); 152*a8db3e61SBlaise Bourdin for (p = 0; p < closureSize; ++p) { 153*a8db3e61SBlaise Bourdin /* Find depth of p */ 154*a8db3e61SBlaise Bourdin for (d = 0; d <= sdim; ++d) { 155*a8db3e61SBlaise Bourdin if ((closure[2*p] >= pStartDepth[d]) && (closure[2*p] < pEndDepth[d])) { 156*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d])); 157*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d])); 158*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d])); 159*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d])); 160*a8db3e61SBlaise Bourdin } 161*a8db3e61SBlaise Bourdin } 162*a8db3e61SBlaise Bourdin } 163*a8db3e61SBlaise Bourdin PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure)); 164*a8db3e61SBlaise Bourdin } 165*a8db3e61SBlaise Bourdin PetscCall(ISRestoreIndices(cellIS, &cellID)); 166*a8db3e61SBlaise Bourdin PetscCall(ISDestroy(&cellIS)); 167*a8db3e61SBlaise Bourdin } 168*a8db3e61SBlaise Bourdin } 169*a8db3e61SBlaise Bourdin if (csIS) {PetscCall(ISRestoreIndices(csIS, &csID));} 170*a8db3e61SBlaise Bourdin PetscCall(ISDestroy(&csIS)); 171*a8db3e61SBlaise Bourdin PetscCall(PetscSectionSetUp(section)); 172*a8db3e61SBlaise Bourdin PetscCall(DMSetLocalSection(dm, section)); 173*a8db3e61SBlaise Bourdin PetscCall(PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view")); 174*a8db3e61SBlaise Bourdin PetscCall(PetscSectionDestroy(§ion)); 175*a8db3e61SBlaise Bourdin 176*a8db3e61SBlaise Bourdin /* Get DM and IS for each field of dm */ 177*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(dm, 1, &fieldU, &seqisU, &seqdmU)); 178*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(dm, 1, &fieldA, &seqisA, &seqdmA)); 179*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(dm, 1, &fieldS, &seqisS, &seqdmS)); 180*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(dm, 2, fieldUA, &seqisUA, &seqdmUA)); 181*a8db3e61SBlaise Bourdin 182*a8db3e61SBlaise Bourdin PetscCall(PetscMalloc1(2,&seqdmList)); 183*a8db3e61SBlaise Bourdin seqdmList[0] = seqdmU; 184*a8db3e61SBlaise Bourdin seqdmList[1] = seqdmA; 185*a8db3e61SBlaise Bourdin 186*a8db3e61SBlaise Bourdin PetscCall(DMCreateSuperDM(seqdmList,2,NULL,&seqdmUA2)); 187*a8db3e61SBlaise Bourdin PetscCall(PetscFree(seqdmList)); 188*a8db3e61SBlaise Bourdin 189*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dm, &seqX)); 190*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmU, &seqU)); 191*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmA, &seqA)); 192*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmS, &seqS)); 193*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmUA, &seqUA)); 194*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmUA2, &seqUA2)); 195*a8db3e61SBlaise Bourdin 196*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) seqX, "seqX")); 197*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) seqU, "seqU")); 198*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) seqA, "seqAlpha")); 199*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) seqS, "seqSigma")); 200*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) seqUA, "seqUAlpha")); 201*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) seqUA2, "seqUAlpha2")); 202*a8db3e61SBlaise Bourdin PetscCall(VecSet(seqX, -111.)); 203*a8db3e61SBlaise Bourdin 204*a8db3e61SBlaise Bourdin /* Setting u to [x,y,z] and alpha to x^2+y^2+z^2 by writing in UAlpha then restricting to U and Alpha */ 205*a8db3e61SBlaise Bourdin { 206*a8db3e61SBlaise Bourdin PetscSection sectionUA; 207*a8db3e61SBlaise Bourdin Vec UALoc; 208*a8db3e61SBlaise Bourdin PetscSection coordSection; 209*a8db3e61SBlaise Bourdin Vec coord; 210*a8db3e61SBlaise Bourdin PetscScalar *cval, *xyz; 211*a8db3e61SBlaise Bourdin PetscInt clSize, i, j; 212*a8db3e61SBlaise Bourdin 213*a8db3e61SBlaise Bourdin PetscCall(DMGetLocalSection(seqdmUA, §ionUA)); 214*a8db3e61SBlaise Bourdin PetscCall(DMGetLocalVector(seqdmUA, &UALoc)); 215*a8db3e61SBlaise Bourdin PetscCall(VecGetArray(UALoc, &cval)); 216*a8db3e61SBlaise Bourdin PetscCall(DMGetCoordinateSection(seqdmUA, &coordSection)); 217*a8db3e61SBlaise Bourdin PetscCall(DMGetCoordinatesLocal(seqdmUA, &coord)); 218*a8db3e61SBlaise Bourdin PetscCall(DMPlexGetChart(seqdmUA, &pStart, &pEnd)); 219*a8db3e61SBlaise Bourdin 220*a8db3e61SBlaise Bourdin for (p = pStart; p < pEnd; ++p) { 221*a8db3e61SBlaise Bourdin PetscInt dofUA, offUA; 222*a8db3e61SBlaise Bourdin 223*a8db3e61SBlaise Bourdin PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA)); 224*a8db3e61SBlaise Bourdin if (dofUA > 0) { 225*a8db3e61SBlaise Bourdin xyz=NULL; 226*a8db3e61SBlaise Bourdin PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA)); 227*a8db3e61SBlaise Bourdin PetscCall(DMPlexVecGetClosure(seqdmUA, coordSection, coord, p, &clSize, &xyz)); 228*a8db3e61SBlaise Bourdin cval[offUA+sdim] = 0.; 229*a8db3e61SBlaise Bourdin for (i = 0; i < sdim; ++i) { 230*a8db3e61SBlaise Bourdin cval[offUA+i] = 0; 231*a8db3e61SBlaise Bourdin for (j = 0; j < clSize/sdim; ++j) { 232*a8db3e61SBlaise Bourdin cval[offUA+i] += xyz[j*sdim+i]; 233*a8db3e61SBlaise Bourdin } 234*a8db3e61SBlaise Bourdin cval[offUA+i] = cval[offUA+i] * sdim / clSize; 235*a8db3e61SBlaise Bourdin cval[offUA+sdim] += PetscSqr(cval[offUA+i]); 236*a8db3e61SBlaise Bourdin } 237*a8db3e61SBlaise Bourdin PetscCall(DMPlexVecRestoreClosure(seqdmUA, coordSection, coord, p, &clSize, &xyz)); 238*a8db3e61SBlaise Bourdin } 239*a8db3e61SBlaise Bourdin } 240*a8db3e61SBlaise Bourdin PetscCall(VecRestoreArray(UALoc, &cval)); 241*a8db3e61SBlaise Bourdin PetscCall(DMLocalToGlobalBegin(seqdmUA, UALoc, INSERT_VALUES, seqUA)); 242*a8db3e61SBlaise Bourdin PetscCall(DMLocalToGlobalEnd(seqdmUA, UALoc, INSERT_VALUES, seqUA)); 243*a8db3e61SBlaise Bourdin PetscCall(DMRestoreLocalVector(seqdmUA, &UALoc)); 244*a8db3e61SBlaise Bourdin 245*a8db3e61SBlaise Bourdin /* Update X */ 246*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(seqX, seqisUA, SCATTER_FORWARD, seqUA)); 247*a8db3e61SBlaise Bourdin PetscCall(VecViewFromOptions(seqUA, NULL, "-ua_vec_view")); 248*a8db3e61SBlaise Bourdin 249*a8db3e61SBlaise Bourdin /* Restrict to U and Alpha */ 250*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(seqX, seqisU, SCATTER_REVERSE, seqU)); 251*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(seqX, seqisA, SCATTER_REVERSE, seqA)); 252*a8db3e61SBlaise Bourdin 253*a8db3e61SBlaise Bourdin /* restrict to UA2 */ 254*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(seqX, seqisUA, SCATTER_REVERSE, seqUA2)); 255*a8db3e61SBlaise Bourdin PetscCall(VecViewFromOptions(seqUA2, NULL, "-ua2_vec_view")); 256*a8db3e61SBlaise Bourdin } 257*a8db3e61SBlaise Bourdin 258*a8db3e61SBlaise Bourdin { 259*a8db3e61SBlaise Bourdin PetscInt ovlp = 0; 260*a8db3e61SBlaise Bourdin PetscPartitioner part; 261*a8db3e61SBlaise Bourdin 262*a8db3e61SBlaise Bourdin PetscCall(DMSetUseNatural(dm,PETSC_TRUE)); 263*a8db3e61SBlaise Bourdin PetscCall(DMPlexGetPartitioner(dm,&part)); 264*a8db3e61SBlaise Bourdin PetscCall(PetscPartitionerSetFromOptions(part)); 265*a8db3e61SBlaise Bourdin PetscCall(DMPlexDistribute(dm,ovlp,&migrationSF,&pdm)); 266*a8db3e61SBlaise Bourdin if (!pdm) pdm = dm; 267*a8db3e61SBlaise Bourdin if (migrationSF) { 268*a8db3e61SBlaise Bourdin PetscCall(DMPlexSetMigrationSF(pdm, migrationSF)); 269*a8db3e61SBlaise Bourdin PetscCall(PetscSFDestroy(&migrationSF)); 270*a8db3e61SBlaise Bourdin } 271*a8db3e61SBlaise Bourdin PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view")); 272*a8db3e61SBlaise Bourdin } 273*a8db3e61SBlaise Bourdin 274*a8db3e61SBlaise Bourdin /* Get DM and IS for each field of dm */ 275*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(pdm, 1, &fieldU, &isU, &dmU)); 276*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(pdm, 1, &fieldA, &isA, &dmA)); 277*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(pdm, 1, &fieldS, &isS, &dmS)); 278*a8db3e61SBlaise Bourdin PetscCall(DMCreateSubDM(pdm, 2, fieldUA, &isUA, &dmUA)); 279*a8db3e61SBlaise Bourdin 280*a8db3e61SBlaise Bourdin PetscCall(PetscMalloc1(2,&dmList)); 281*a8db3e61SBlaise Bourdin dmList[0] = dmU; 282*a8db3e61SBlaise Bourdin dmList[1] = dmA; 283*a8db3e61SBlaise Bourdin 284*a8db3e61SBlaise Bourdin PetscCall(DMCreateSuperDM(dmList,2,NULL,&dmUA2)); 285*a8db3e61SBlaise Bourdin PetscCall(PetscFree(dmList)); 286*a8db3e61SBlaise Bourdin 287*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(pdm, &X)); 288*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dmU, &U)); 289*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dmA, &A)); 290*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dmS, &S)); 291*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dmUA, &UA)); 292*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dmUA2, &UA2)); 293*a8db3e61SBlaise Bourdin 294*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) X, "X")); 295*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) U, "U")); 296*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) A, "Alpha")); 297*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) S, "Sigma")); 298*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) UA, "UAlpha")); 299*a8db3e61SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject) UA2, "UAlpha2")); 300*a8db3e61SBlaise Bourdin PetscCall(VecSet(X, -111.)); 301*a8db3e61SBlaise Bourdin 302*a8db3e61SBlaise Bourdin /* Setting u to [x,y,z] and alpha to x^2+y^2+z^2 by writing in UAlpha then restricting to U and Alpha */ 303*a8db3e61SBlaise Bourdin { 304*a8db3e61SBlaise Bourdin PetscSection sectionUA; 305*a8db3e61SBlaise Bourdin Vec UALoc; 306*a8db3e61SBlaise Bourdin PetscSection coordSection; 307*a8db3e61SBlaise Bourdin Vec coord; 308*a8db3e61SBlaise Bourdin PetscScalar *cval, *xyz; 309*a8db3e61SBlaise Bourdin PetscInt clSize, i, j; 310*a8db3e61SBlaise Bourdin 311*a8db3e61SBlaise Bourdin PetscCall(DMGetLocalSection(dmUA, §ionUA)); 312*a8db3e61SBlaise Bourdin PetscCall(DMGetLocalVector(dmUA, &UALoc)); 313*a8db3e61SBlaise Bourdin PetscCall(VecGetArray(UALoc, &cval)); 314*a8db3e61SBlaise Bourdin PetscCall(DMGetCoordinateSection(dmUA, &coordSection)); 315*a8db3e61SBlaise Bourdin PetscCall(DMGetCoordinatesLocal(dmUA, &coord)); 316*a8db3e61SBlaise Bourdin PetscCall(DMPlexGetChart(dmUA, &pStart, &pEnd)); 317*a8db3e61SBlaise Bourdin 318*a8db3e61SBlaise Bourdin for (p = pStart; p < pEnd; ++p) { 319*a8db3e61SBlaise Bourdin PetscInt dofUA, offUA; 320*a8db3e61SBlaise Bourdin 321*a8db3e61SBlaise Bourdin PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA)); 322*a8db3e61SBlaise Bourdin if (dofUA > 0) { 323*a8db3e61SBlaise Bourdin xyz=NULL; 324*a8db3e61SBlaise Bourdin PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA)); 325*a8db3e61SBlaise Bourdin PetscCall(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz)); 326*a8db3e61SBlaise Bourdin cval[offUA+sdim] = 0.; 327*a8db3e61SBlaise Bourdin for (i = 0; i < sdim; ++i) { 328*a8db3e61SBlaise Bourdin cval[offUA+i] = 0; 329*a8db3e61SBlaise Bourdin for (j = 0; j < clSize/sdim; ++j) { 330*a8db3e61SBlaise Bourdin cval[offUA+i] += xyz[j*sdim+i]; 331*a8db3e61SBlaise Bourdin } 332*a8db3e61SBlaise Bourdin cval[offUA+i] = cval[offUA+i] * sdim / clSize; 333*a8db3e61SBlaise Bourdin cval[offUA+sdim] += PetscSqr(cval[offUA+i]); 334*a8db3e61SBlaise Bourdin } 335*a8db3e61SBlaise Bourdin PetscCall(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz)); 336*a8db3e61SBlaise Bourdin } 337*a8db3e61SBlaise Bourdin } 338*a8db3e61SBlaise Bourdin PetscCall(VecRestoreArray(UALoc, &cval)); 339*a8db3e61SBlaise Bourdin PetscCall(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA)); 340*a8db3e61SBlaise Bourdin PetscCall(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA)); 341*a8db3e61SBlaise Bourdin PetscCall(DMRestoreLocalVector(dmUA, &UALoc)); 342*a8db3e61SBlaise Bourdin 343*a8db3e61SBlaise Bourdin /* Update X */ 344*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(X, isUA, SCATTER_FORWARD, UA)); 345*a8db3e61SBlaise Bourdin PetscCall(VecViewFromOptions(UA, NULL, "-ua_vec_view")); 346*a8db3e61SBlaise Bourdin 347*a8db3e61SBlaise Bourdin /* Restrict to U and Alpha */ 348*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(X, isU, SCATTER_REVERSE, U)); 349*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(X, isA, SCATTER_REVERSE, A)); 350*a8db3e61SBlaise Bourdin 351*a8db3e61SBlaise Bourdin /* restrict to UA2 */ 352*a8db3e61SBlaise Bourdin PetscCall(VecISCopy(X, isUA, SCATTER_REVERSE, UA2)); 353*a8db3e61SBlaise Bourdin PetscCall(VecViewFromOptions(UA2, NULL, "-ua2_vec_view")); 354*a8db3e61SBlaise Bourdin } 355*a8db3e61SBlaise Bourdin 356*a8db3e61SBlaise Bourdin /* Verification */ 357*a8db3e61SBlaise Bourdin 358*a8db3e61SBlaise Bourdin Vec Xnat, Unat, Anat, UAnat, Snat, UA2nat; 359*a8db3e61SBlaise Bourdin PetscReal norm; 360*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(dm, &Xnat)); 361*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmU, &Unat)); 362*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmA, &Anat)); 363*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmUA, &UAnat)); 364*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmS, &Snat)); 365*a8db3e61SBlaise Bourdin PetscCall(DMGetGlobalVector(seqdmUA2, &UA2nat)); 366*a8db3e61SBlaise Bourdin 367*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(pdm, X, Xnat)); 368*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(pdm, X, Xnat)); 369*a8db3e61SBlaise Bourdin PetscCall(VecAXPY(Xnat, -1.0, seqX)); 370*a8db3e61SBlaise Bourdin PetscCall(VecNorm(Xnat, NORM_INFINITY, &norm)); 371*a8db3e61SBlaise Bourdin PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "X ||Vin - Vout|| = %g", (double) norm); 372*a8db3e61SBlaise Bourdin 373*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dmU, U, Unat)); 374*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dmU, U, Unat)); 375*a8db3e61SBlaise Bourdin PetscCall(VecAXPY(Unat, -1.0, seqU)); 376*a8db3e61SBlaise Bourdin PetscCall(VecNorm(Unat, NORM_INFINITY, &norm)); 377*a8db3e61SBlaise Bourdin PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "U ||Vin - Vout|| = %g", (double) norm); 378*a8db3e61SBlaise Bourdin 379*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dmA, A, Anat)); 380*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dmA, A, Anat)); 381*a8db3e61SBlaise Bourdin PetscCall(VecAXPY(Anat, -1.0, seqA)); 382*a8db3e61SBlaise Bourdin PetscCall(VecNorm(Anat, NORM_INFINITY, &norm)); 383*a8db3e61SBlaise Bourdin PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "A ||Vin - Vout|| = %g", (double) norm); 384*a8db3e61SBlaise Bourdin 385*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dmUA, UA, UAnat)); 386*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dmUA, UA, UAnat)); 387*a8db3e61SBlaise Bourdin PetscCall(VecAXPY(UAnat, -1.0, seqUA)); 388*a8db3e61SBlaise Bourdin PetscCall(VecNorm(UAnat, NORM_INFINITY, &norm)); 389*a8db3e61SBlaise Bourdin PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UA ||Vin - Vout|| = %g", (double) norm); 390*a8db3e61SBlaise Bourdin 391*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dmS, S, Snat)); 392*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dmS, S, Snat)); 393*a8db3e61SBlaise Bourdin PetscCall(VecAXPY(Snat, -1.0, seqS)); 394*a8db3e61SBlaise Bourdin PetscCall(VecNorm(Snat, NORM_INFINITY, &norm)); 395*a8db3e61SBlaise Bourdin PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "S ||Vin - Vout|| = %g", (double) norm); 396*a8db3e61SBlaise Bourdin 397*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dmUA2, UA2, UA2nat)); 398*a8db3e61SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dmUA2, UA2, UA2nat)); 399*a8db3e61SBlaise Bourdin PetscCall(VecAXPY(UA2nat, -1.0, seqUA2)); 400*a8db3e61SBlaise Bourdin PetscCall(VecNorm(UA2nat, NORM_INFINITY, &norm)); 401*a8db3e61SBlaise Bourdin PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double) norm); 402*a8db3e61SBlaise Bourdin 403*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dm, &Xnat)); 404*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmU, &Unat)); 405*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmA, &Anat)); 406*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmUA, &UAnat)); 407*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmS, &Snat)); 408*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmUA2, &UA2nat)); 409*a8db3e61SBlaise Bourdin 410*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmUA2, &seqUA2)); 411*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmUA, &seqUA)); 412*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmS, &seqS)); 413*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmA, &seqA)); 414*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(seqdmU, &seqU)); 415*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dm, &seqX)); 416*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&seqdmU));PetscCall(ISDestroy(&seqisU)); 417*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&seqdmA));PetscCall(ISDestroy(&seqisA)); 418*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&seqdmS));PetscCall(ISDestroy(&seqisS)); 419*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&seqdmUA));PetscCall(ISDestroy(&seqisUA)); 420*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&seqdmUA2)); 421*a8db3e61SBlaise Bourdin 422*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dmUA2, &UA2)); 423*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dmUA, &UA)); 424*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dmS, &S)); 425*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dmA, &A)); 426*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(dmU, &U)); 427*a8db3e61SBlaise Bourdin PetscCall(DMRestoreGlobalVector(pdm, &X)); 428*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&dmU)); PetscCall(ISDestroy(&isU)); 429*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&dmA)); PetscCall(ISDestroy(&isA)); 430*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&dmS)); PetscCall(ISDestroy(&isS)); 431*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&dmUA));PetscCall(ISDestroy(&isUA)); 432*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&dmUA2)); 433*a8db3e61SBlaise Bourdin PetscCall(DMDestroy(&pdm)); 434*a8db3e61SBlaise Bourdin if (size > 1) PetscCall(DMDestroy(&dm)); 435*a8db3e61SBlaise Bourdin PetscCall(PetscFree2(pStartDepth, pEndDepth)); 436*a8db3e61SBlaise Bourdin PetscCall(PetscFinalize()); 437*a8db3e61SBlaise Bourdin return 0; 438*a8db3e61SBlaise Bourdin } 439*a8db3e61SBlaise Bourdin 440*a8db3e61SBlaise Bourdin /*TEST 441*a8db3e61SBlaise Bourdin 442*a8db3e61SBlaise Bourdin build: 443*a8db3e61SBlaise Bourdin requires: !complex parmetis exodusii pnetcdf 444*a8db3e61SBlaise Bourdin # 2D seq 445*a8db3e61SBlaise Bourdin test: 446*a8db3e61SBlaise Bourdin suffix: 0 447*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -dm_view -petscpartitioner_type parmetis 448*a8db3e61SBlaise Bourdin test: 449*a8db3e61SBlaise Bourdin suffix: 1 450*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -dm_view -petscpartitioner_type parmetis 451*a8db3e61SBlaise Bourdin test: 452*a8db3e61SBlaise Bourdin suffix: 2 453*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -dm_view -petscpartitioner_type parmetis 454*a8db3e61SBlaise Bourdin 455*a8db3e61SBlaise Bourdin # 2D par 456*a8db3e61SBlaise Bourdin test: 457*a8db3e61SBlaise Bourdin suffix: 3 458*a8db3e61SBlaise Bourdin nsize: 2 459*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -dm_view -petscpartitioner_type parmetis 460*a8db3e61SBlaise Bourdin test: 461*a8db3e61SBlaise Bourdin suffix: 4 462*a8db3e61SBlaise Bourdin nsize: 2 463*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -dm_view -petscpartitioner_type parmetis 464*a8db3e61SBlaise Bourdin test: 465*a8db3e61SBlaise Bourdin suffix: 5 466*a8db3e61SBlaise Bourdin nsize: 2 467*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -dm_view -petscpartitioner_type parmetis 468*a8db3e61SBlaise Bourdin 469*a8db3e61SBlaise Bourdin #3d seq 470*a8db3e61SBlaise Bourdin test: 471*a8db3e61SBlaise Bourdin suffix: 6 472*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -dm_view -petscpartitioner_type parmetis 473*a8db3e61SBlaise Bourdin test: 474*a8db3e61SBlaise Bourdin suffix: 7 475*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -dm_view -petscpartitioner_type parmetis 476*a8db3e61SBlaise Bourdin 477*a8db3e61SBlaise Bourdin #3d par 478*a8db3e61SBlaise Bourdin test: 479*a8db3e61SBlaise Bourdin suffix: 8 480*a8db3e61SBlaise Bourdin nsize: 2 481*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -dm_view -petscpartitioner_type parmetis 482*a8db3e61SBlaise Bourdin test: 483*a8db3e61SBlaise Bourdin suffix: 9 484*a8db3e61SBlaise Bourdin nsize: 2 485*a8db3e61SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -dm_view -petscpartitioner_type parmetis 486*a8db3e61SBlaise Bourdin 487*a8db3e61SBlaise Bourdin TEST*/ 488