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