1280aadf6SBlaise Bourdin#include "petsc/finclude/petsc.h" 2c5e229c2SMartin Diehlprogram ex62f90 3280aadf6SBlaise Bourdin use petsc 4280aadf6SBlaise Bourdin implicit none 5280aadf6SBlaise Bourdin#include "exodusII.inc" 6280aadf6SBlaise Bourdin 7dd01b7e5SBarry Smith ! Get the Fortran kind associated with PetscInt and PetscReal so that we can use literal constants. 8280aadf6SBlaise Bourdin PetscInt :: dummyPetscInt 9280aadf6SBlaise Bourdin PetscReal :: dummyPetscreal 10280aadf6SBlaise Bourdin integer, parameter :: kPI = kind(dummyPetscInt) 11280aadf6SBlaise Bourdin integer, parameter :: kPR = kind(dummyPetscReal) 12280aadf6SBlaise Bourdin 13280aadf6SBlaise Bourdin type(tDM) :: dm, dmU, dmA, dmS, dmUA, dmUA2, pDM 14280aadf6SBlaise Bourdin type(tDM), dimension(:), pointer :: dmList 15280aadf6SBlaise Bourdin type(tVec) :: X, U, A, S, UA, UA2 16280aadf6SBlaise Bourdin type(tIS) :: isU, isA, isS, isUA 17e2739ba6SAlexis Marboeuf type(tPetscSection) :: section, rootSection, leafSection 18ce78bad3SBarry Smith PetscInt :: fieldU = 0 19ce78bad3SBarry Smith PetscInt :: fieldA = 2 20ce78bad3SBarry Smith PetscInt :: fieldS = 1 21280aadf6SBlaise Bourdin PetscInt, dimension(2) :: fieldUA = [0, 2] 22280aadf6SBlaise Bourdin character(len=PETSC_MAX_PATH_LEN) :: ifilename, ofilename, IOBuffer 23280aadf6SBlaise Bourdin integer :: exoid = -1 24280aadf6SBlaise Bourdin type(tIS) :: csIS 25280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: csID 26280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: pStartDepth, pEndDepth 27280aadf6SBlaise Bourdin PetscInt :: order = 1 28280aadf6SBlaise Bourdin PetscInt :: sdim, d, pStart, pEnd, p, numCS, set, i, j 29280aadf6SBlaise Bourdin PetscMPIInt :: rank, numProc 30280aadf6SBlaise Bourdin PetscBool :: flg 31280aadf6SBlaise Bourdin PetscErrorCode :: ierr 32*b06eb4cdSBarry Smith MPIU_Comm :: comm 33280aadf6SBlaise Bourdin type(tPetscViewer) :: viewer 34280aadf6SBlaise Bourdin 3502c639afSMartin Diehl character(len=MXSTLN) :: sJunk 36280aadf6SBlaise Bourdin PetscInt :: numstep = 3, step 37280aadf6SBlaise Bourdin PetscInt :: numNodalVar, numZonalVar 38280aadf6SBlaise Bourdin character(len=MXSTLN) :: nodalVarName(4) 39280aadf6SBlaise Bourdin character(len=MXSTLN) :: zonalVarName(6) 40280aadf6SBlaise Bourdin logical, dimension(:, :), pointer :: truthtable 41280aadf6SBlaise Bourdin 42280aadf6SBlaise Bourdin type(tIS) :: cellIS 43280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: cellID 44280aadf6SBlaise Bourdin PetscInt :: numCells, cell, closureSize 45280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: closureA, closure 46280aadf6SBlaise Bourdin 47280aadf6SBlaise Bourdin type(tPetscSection) :: sectionUA, coordSection 48280aadf6SBlaise Bourdin type(tVec) :: UALoc, coord 49280aadf6SBlaise Bourdin PetscScalar, dimension(:), pointer :: cval, xyz 50280aadf6SBlaise Bourdin PetscInt :: dofUA, offUA, c 51280aadf6SBlaise Bourdin 52280aadf6SBlaise Bourdin ! dof layout ordered by increasing height in the DAG: cell, face, edge, vertex 53280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofS2D = [0, 0, 3] 54280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP1Tri = [2, 0, 0] 55280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP1Tri = [1, 0, 0] 56280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP2Tri = [2, 2, 0] 57280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP2Tri = [1, 1, 0] 58280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP1Quad = [2, 0, 0] 59280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP1Quad = [1, 0, 0] 60280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP2Quad = [2, 2, 2] 61280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP2Quad = [1, 1, 1] 62280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofS3D = [0, 0, 0, 6] 63280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP1Tet = [3, 0, 0, 0] 64280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP1Tet = [1, 0, 0, 0] 65280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP2Tet = [3, 3, 0, 0] 66280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP2Tet = [1, 1, 0, 0] 67280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP1Hex = [3, 0, 0, 0] 68280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP1Hex = [1, 0, 0, 0] 69280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP2Hex = [3, 3, 3, 3] 70280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP2Hex = [1, 1, 1, 1] 71280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: dofU, dofA, dofS 72280aadf6SBlaise Bourdin 73e2739ba6SAlexis Marboeuf type(tPetscSF) :: migrationSF, natSF, natPointSF, natPointSFInv 74280aadf6SBlaise Bourdin PetscPartitioner :: part 75280aadf6SBlaise Bourdin 76280aadf6SBlaise Bourdin type(tVec) :: tmpVec 77280aadf6SBlaise Bourdin PetscReal :: norm 78280aadf6SBlaise Bourdin PetscReal :: time = 1.234_kPR 79ce78bad3SBarry Smith PetscInt, pointer :: remoteOffsets(:) 80280aadf6SBlaise Bourdin 81e2739ba6SAlexis Marboeuf PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr)) 82e2739ba6SAlexis Marboeuf if (ierr /= 0) then 83e2739ba6SAlexis Marboeuf print *, 'Unable to initialize PETSc' 84e2739ba6SAlexis Marboeuf stop 85e2739ba6SAlexis Marboeuf end if 86280aadf6SBlaise Bourdin 87e2739ba6SAlexis Marboeuf PetscCallA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) 88e2739ba6SAlexis Marboeuf PetscCallA(MPI_Comm_size(PETSC_COMM_WORLD, numProc, ierr)) 89dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr)) 90dcb3e689SBarry Smith PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i <input file name>') 91dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-o', ofilename, flg, ierr)) 92dcb3e689SBarry Smith PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing output file name -o <output file name>') 93dcb3e689SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-order', order, flg, ierr)) 94280aadf6SBlaise Bourdin if ((order > 2) .or. (order < 1)) then 95280aadf6SBlaise Bourdin write (IOBuffer, '("Unsupported polynomial order ", I2, " not in [1,2]")') order 96e2447916SStefano Zampini stop 97280aadf6SBlaise Bourdin end if 98280aadf6SBlaise Bourdin 99280aadf6SBlaise Bourdin ! Read the mesh in any supported format 100d8606c27SBarry Smith PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_NULL_CHARACTER, PETSC_TRUE, dm, ierr)) 101d8606c27SBarry Smith PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE, ierr)) 102d8606c27SBarry Smith PetscCallA(DMSetFromOptions(dm, ierr)) 103d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim, ierr)) 104ce78bad3SBarry Smith PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr)) 105280aadf6SBlaise Bourdin 106280aadf6SBlaise Bourdin ! Create the exodus result file 107280aadf6SBlaise Bourdin 108a5b23f4aSJose E. Roman ! enable exodus debugging information 109d8606c27SBarry Smith PetscCallA(exopts(EXVRBS + EXDEBG, ierr)) 110280aadf6SBlaise Bourdin ! Create the exodus file 111d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_WRITE, viewer, ierr)) 112280aadf6SBlaise Bourdin ! The long way would be 113280aadf6SBlaise Bourdin ! 114d8606c27SBarry Smith ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr)) 115d8606c27SBarry Smith ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr)) 116d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr)) 117d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr)) 118280aadf6SBlaise Bourdin 119280aadf6SBlaise Bourdin ! set the mesh order 120d8606c27SBarry Smith PetscCallA(PetscViewerExodusIISetOrder(viewer, order, ierr)) 121d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr)) 122280aadf6SBlaise Bourdin ! 123280aadf6SBlaise Bourdin ! Notice how the exodus file is actually NOT open at this point (exoid is -1) 124d5b43468SJose E. Roman ! Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to 125280aadf6SBlaise Bourdin ! write the geometry (the DM), which can only be done on a brand new file. 126280aadf6SBlaise Bourdin ! 127280aadf6SBlaise Bourdin 128280aadf6SBlaise Bourdin ! Save the geometry to the file, erasing all previous content 129d8606c27SBarry Smith PetscCallA(DMView(dm, viewer, ierr)) 130d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr)) 131280aadf6SBlaise Bourdin ! 132280aadf6SBlaise Bourdin ! Note how the exodus file is now open 133280aadf6SBlaise Bourdin ! 134dcb3e689SBarry Smith ! 'Format' the exodus result file, i.e. allocate space for nodal and zonal variables 135280aadf6SBlaise Bourdin select case (sdim) 136280aadf6SBlaise Bourdin case (2) 137280aadf6SBlaise Bourdin numNodalVar = 3 138dcb3e689SBarry Smith nodalVarName(1:numNodalVar) = ['U_x ', 'U_y ', 'Alpha'] 139280aadf6SBlaise Bourdin numZonalVar = 3 140dcb3e689SBarry Smith zonalVarName(1:numZonalVar) = ['Sigma_11', 'Sigma_22', 'Sigma_12'] 141280aadf6SBlaise Bourdin case (3) 142280aadf6SBlaise Bourdin numNodalVar = 4 143dcb3e689SBarry Smith nodalVarName(1:numNodalVar) = ['U_x ', 'U_y ', 'U_z ', 'Alpha'] 144280aadf6SBlaise Bourdin numZonalVar = 6 145dcb3e689SBarry Smith zonalVarName(1:numZonalVar) = ['Sigma_11', 'Sigma_22', 'Sigma_33', 'Sigma_23', 'Sigma_13', 'Sigma_12'] 146280aadf6SBlaise Bourdin case default 147280aadf6SBlaise Bourdin write (IOBuffer, '("No layout for dimension ",I2)') sdim 148280aadf6SBlaise Bourdin end select 149d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIGetId(viewer, exoid, ierr)) 150dcb3e689SBarry Smith PetscCallA(expvp(exoid, 'E', numZonalVar, ierr)) 151dcb3e689SBarry Smith PetscCallA(expvan(exoid, 'E', numZonalVar, zonalVarName, ierr)) 152dcb3e689SBarry Smith PetscCallA(expvp(exoid, 'N', numNodalVar, ierr)) 153dcb3e689SBarry Smith PetscCallA(expvan(exoid, 'N', numNodalVar, nodalVarName, ierr)) 154d8606c27SBarry Smith PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK, numCS, PETSC_NULL_REAL, sjunk, ierr)) 155280aadf6SBlaise Bourdin 156280aadf6SBlaise Bourdin ! An exodusII truth table specifies which fields are saved at which time step 157280aadf6SBlaise Bourdin ! It speeds up I/O but reserving space for fields in the file ahead of time. 158280aadf6SBlaise Bourdin allocate (truthtable(numCS, numZonalVar)) 159280aadf6SBlaise Bourdin truthtable = .true. 160d8606c27SBarry Smith PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr)) 161280aadf6SBlaise Bourdin deallocate (truthtable) 162280aadf6SBlaise Bourdin 163280aadf6SBlaise Bourdin ! Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 164280aadf6SBlaise Bourdin do step = 1, numstep 16502c639afSMartin Diehl PetscCallA(exptim(exoid, step, real(step, kind=kPR), ierr)) 166280aadf6SBlaise Bourdin end do 167280aadf6SBlaise Bourdin 168d8606c27SBarry Smith PetscCallA(DMSetUseNatural(dm, PETSC_TRUE, ierr)) 169d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm, part, ierr)) 170d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part, ierr)) 171d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm, 0_kPI, migrationSF, pdm, ierr)) 172280aadf6SBlaise Bourdin 173280aadf6SBlaise Bourdin if (numProc > 1) then 174d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm, migrationSF, ierr)) 175d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF, ierr)) 176e2739ba6SAlexis Marboeuf else 177e2739ba6SAlexis Marboeuf pdm = dm 178280aadf6SBlaise Bourdin end if 179ce78bad3SBarry Smith PetscCallA(DMViewFromOptions(pdm, PETSC_NULL_OBJECT, '-dm_view', ierr)) 180280aadf6SBlaise Bourdin 181e2739ba6SAlexis Marboeuf PetscCallA(PetscObjectGetComm(pdm, comm, ierr)) 182d8606c27SBarry Smith PetscCallA(PetscSectionCreate(comm, section, ierr)) 183d8606c27SBarry Smith PetscCallA(PetscSectionSetNumFields(section, 3_kPI, ierr)) 184dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, 'U', ierr)) 185dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, 'Alpha', ierr)) 186dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldS, 'Sigma', ierr)) 187e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetChart(pdm, pStart, pEnd, ierr)) 188d8606c27SBarry Smith PetscCallA(PetscSectionSetChart(section, pStart, pEnd, ierr)) 189280aadf6SBlaise Bourdin 190280aadf6SBlaise Bourdin allocate (pStartDepth(sdim + 1)) 191280aadf6SBlaise Bourdin allocate (pEndDepth(sdim + 1)) 192280aadf6SBlaise Bourdin do d = 1, sdim + 1 193e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetDepthStratum(pdm, d - 1, pStartDepth(d), pEndDepth(d), ierr)) 194280aadf6SBlaise Bourdin end do 195280aadf6SBlaise Bourdin 196280aadf6SBlaise Bourdin ! Vector field U, Scalar field Alpha, Tensor field Sigma 197ccfd86f1SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim, ierr)) 198ccfd86f1SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI, ierr)) 199ccfd86f1SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim + 1)/2, ierr)) 200280aadf6SBlaise Bourdin 201280aadf6SBlaise Bourdin ! Going through cell sets then cells, and setting up storage for the sections 202dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(pdm, 'Cell Sets', numCS, ierr)) 203dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(pdm, 'Cell Sets', csIS, ierr)) 204ce78bad3SBarry Smith PetscCallA(ISGetIndices(csIS, csID, ierr)) 205280aadf6SBlaise Bourdin do set = 1, numCS 206dcb3e689SBarry Smith PetscCallA(DMGetStratumSize(pdm, 'Cell Sets', csID(set), numCells, ierr)) 207dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(pdm, 'Cell Sets', csID(set), cellIS, ierr)) 208280aadf6SBlaise Bourdin if (numCells > 0) then 209280aadf6SBlaise Bourdin select case (sdim) 210280aadf6SBlaise Bourdin case (2) 211280aadf6SBlaise Bourdin dofs => dofS2D 212280aadf6SBlaise Bourdin case (3) 213280aadf6SBlaise Bourdin dofs => dofS3D 214280aadf6SBlaise Bourdin case default 215280aadf6SBlaise Bourdin write (IOBuffer, '("No layout for dimension ",I2)') sdim 216e2447916SStefano Zampini stop 217280aadf6SBlaise Bourdin end select ! sdim 218280aadf6SBlaise Bourdin 219280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 220280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 221ce78bad3SBarry Smith PetscCallA(ISGetIndices(cellIS, cellID, ierr)) 222280aadf6SBlaise Bourdin nullify (closureA) 223ce78bad3SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(pdm, cellID(1), PETSC_TRUE, PETSC_NULL_INTEGER, closureA, ierr)) 224280aadf6SBlaise Bourdin select case (size(closureA)/2) 225280aadf6SBlaise Bourdin case (7) ! Tri 226280aadf6SBlaise Bourdin if (order == 1) then 227280aadf6SBlaise Bourdin dofU => dofUP1Tri 228280aadf6SBlaise Bourdin dofA => dofAP1Tri 229280aadf6SBlaise Bourdin else 230280aadf6SBlaise Bourdin dofU => dofUP2Tri 231280aadf6SBlaise Bourdin dofA => dofAP2Tri 232280aadf6SBlaise Bourdin end if 233280aadf6SBlaise Bourdin case (9) ! Quad 234280aadf6SBlaise Bourdin if (order == 1) then 235280aadf6SBlaise Bourdin dofU => dofUP1Quad 236280aadf6SBlaise Bourdin dofA => dofAP1Quad 237280aadf6SBlaise Bourdin else 238280aadf6SBlaise Bourdin dofU => dofUP2Quad 239280aadf6SBlaise Bourdin dofA => dofAP2Quad 240280aadf6SBlaise Bourdin end if 241280aadf6SBlaise Bourdin case (15) ! Tet 242280aadf6SBlaise Bourdin if (order == 1) then 243280aadf6SBlaise Bourdin dofU => dofUP1Tet 244280aadf6SBlaise Bourdin dofA => dofAP1Tet 245280aadf6SBlaise Bourdin else 246280aadf6SBlaise Bourdin dofU => dofUP2Tet 247280aadf6SBlaise Bourdin dofA => dofAP2Tet 248280aadf6SBlaise Bourdin end if 249280aadf6SBlaise Bourdin case (27) ! Hex 250280aadf6SBlaise Bourdin if (order == 1) then 251280aadf6SBlaise Bourdin dofU => dofUP1Hex 252280aadf6SBlaise Bourdin dofA => dofAP1Hex 253280aadf6SBlaise Bourdin else 254280aadf6SBlaise Bourdin dofU => dofUP2Hex 255280aadf6SBlaise Bourdin dofA => dofAP2Hex 256280aadf6SBlaise Bourdin end if 257280aadf6SBlaise Bourdin case default 258280aadf6SBlaise Bourdin write (IOBuffer, '("Unknown element with closure size ",I2)') size(closureA)/2 259e2447916SStefano Zampini stop 260280aadf6SBlaise Bourdin end select 261ce78bad3SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(pdm, cellID(1), PETSC_TRUE, PETSC_NULL_INTEGER, closureA, ierr)) 262280aadf6SBlaise Bourdin do cell = 1, numCells! 263280aadf6SBlaise Bourdin nullify (closure) 264ce78bad3SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(pdm, cellID(cell), PETSC_TRUE, PETSC_NULL_INTEGER, closure, ierr)) 265280aadf6SBlaise Bourdin do p = 1, size(closure), 2 266280aadf6SBlaise Bourdin ! find the depth of p 267280aadf6SBlaise Bourdin do d = 1, sdim + 1 268280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 269d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d) + dofA(d) + dofS(d), ierr)) 270d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d), ierr)) 271d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d), ierr)) 272d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d), ierr)) 273280aadf6SBlaise Bourdin end if ! closure(p) 274280aadf6SBlaise Bourdin end do ! d 275280aadf6SBlaise Bourdin end do ! p 276ce78bad3SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(pdm, cellID(cell), PETSC_TRUE, PETSC_NULL_INTEGER, closure, ierr)) 277280aadf6SBlaise Bourdin end do ! cell 278ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(cellIS, cellID, ierr)) 279d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS, ierr)) 280280aadf6SBlaise Bourdin end if ! numCells 281280aadf6SBlaise Bourdin end do ! set 282ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(csIS, csID, ierr)) 283d8606c27SBarry Smith PetscCallA(ISDestroy(csIS, ierr)) 284d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section, ierr)) 285e2739ba6SAlexis Marboeuf PetscCallA(DMSetLocalSection(pdm, section, ierr)) 286ce78bad3SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_OBJECT, '-dm_section_view', ierr)) 287d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section, ierr)) 288280aadf6SBlaise Bourdin 289e2739ba6SAlexis Marboeuf ! Creating section on the sequential DM + creating the GlobalToNatural SF 290e2739ba6SAlexis Marboeuf if (numProc > 1) then 291e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetMigrationSF(pdm, natPointSF, ierr)) 292e2739ba6SAlexis Marboeuf PetscCallA(DMGetLocalSection(pdm, rootSection, ierr)) 293e2739ba6SAlexis Marboeuf PetscCallA(PetscSFCreateInverseSF(natPointSF, natPointSFInv, ierr)) 294e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionCreate(PETSC_COMM_WORLD, leafSection, ierr)) 295ce78bad3SBarry Smith PetscCallA(PetscSFDistributeSection(natPointSFInv, rootSection, remoteOffsets, leafSection, ierr)) 296ce78bad3SBarry Smith PetscCallA(PetscSFDestroyRemoteOffsets(remoteOffsets, ierr)) 297e2739ba6SAlexis Marboeuf PetscCallA(DMSetLocalSection(dm, leafSection, ierr)) 298e2739ba6SAlexis Marboeuf PetscCallA(DMPlexCreateGlobalToNaturalSF(pdm, leafSection, natPointSF, natSF, ierr)) 299e2739ba6SAlexis Marboeuf PetscCallA(PetscSFDestroy(natPointSFInv, ierr)) 300e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionDestroy(leafSection, ierr)) 301e2739ba6SAlexis Marboeuf PetscCallA(DMSetNaturalSF(pdm, natSF, ierr)) 302e2739ba6SAlexis Marboeuf PetscCallA(PetscObjectDereference(natSF, ierr)) 303e2739ba6SAlexis Marboeuf end if 304e2739ba6SAlexis Marboeuf 305280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 306ce78bad3SBarry Smith PetscCallA(DMCreateSubDM(pdm, 1_kPI, [fieldU], isU, dmU, ierr)) 307ce78bad3SBarry Smith PetscCallA(DMCreateSubDM(pdm, 1_kPI, [fieldA], isA, dmA, ierr)) 308ce78bad3SBarry Smith PetscCallA(DMCreateSubDM(pdm, 1_kPI, [fieldS], isS, dmS, ierr)) 309e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA, ierr)) 310280aadf6SBlaise Bourdin 311280aadf6SBlaise Bourdin !Create the exodus result file 312280aadf6SBlaise Bourdin allocate (dmList(2)) 313ccfd86f1SBarry Smith dmList(1) = dmU 314ccfd86f1SBarry Smith dmList(2) = dmA 315ce78bad3SBarry Smith PetscCallA(DMCreateSuperDM(dmList, 2_kPI, PETSC_NULL_IS_POINTER, dmUA2, ierr)) 316280aadf6SBlaise Bourdin deallocate (dmList) 317280aadf6SBlaise Bourdin 318e2739ba6SAlexis Marboeuf PetscCallA(DMGetGlobalVector(pdm, X, ierr)) 319d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U, ierr)) 320d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A, ierr)) 321d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S, ierr)) 322d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA, ierr)) 323d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2, ierr)) 324280aadf6SBlaise Bourdin 325dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(U, 'U', ierr)) 326dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(A, 'Alpha', ierr)) 327dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(S, 'Sigma', ierr)) 328dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA, 'UAlpha', ierr)) 329dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA2, 'UAlpha2', ierr)) 330d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR, ierr)) 331280aadf6SBlaise Bourdin 332280aadf6SBlaise 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 */ 333d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA, ierr)) 334d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc, ierr)) 335ce78bad3SBarry Smith PetscCallA(VecGetArray(UALoc, cval, ierr)) 336d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection, ierr)) 337d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord, ierr)) 338d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd, ierr)) 339280aadf6SBlaise Bourdin 340280aadf6SBlaise Bourdin do p = pStart, pEnd - 1 341d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA, ierr)) 342280aadf6SBlaise Bourdin if (dofUA > 0) then 343d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA, ierr)) 344ce78bad3SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, PETSC_NULL_INTEGER, xyz, ierr)) 345280aadf6SBlaise Bourdin closureSize = size(xyz) 346280aadf6SBlaise Bourdin do i = 1, sdim 347280aadf6SBlaise Bourdin do j = 0, closureSize - 1, sdim 348280aadf6SBlaise Bourdin cval(offUA + i) = cval(offUA + i) + xyz(j/sdim + i) 349280aadf6SBlaise Bourdin end do 350ccfd86f1SBarry Smith cval(offUA + i) = cval(offUA + i)*sdim/closureSize 351280aadf6SBlaise Bourdin cval(offUA + sdim + 1) = cval(offUA + sdim + 1) + cval(offUA + i)**2 352280aadf6SBlaise Bourdin end do 353ce78bad3SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, PETSC_NULL_INTEGER, xyz, ierr)) 354280aadf6SBlaise Bourdin end if 355280aadf6SBlaise Bourdin end do 356280aadf6SBlaise Bourdin 357ce78bad3SBarry Smith PetscCallA(VecRestoreArray(UALoc, cval, ierr)) 358d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA, ierr)) 359d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA, ierr)) 360d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc, ierr)) 361280aadf6SBlaise Bourdin 362280aadf6SBlaise Bourdin !Update X 363d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA, ierr)) 364280aadf6SBlaise Bourdin ! Restrict to U and Alpha 365d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U, ierr)) 366d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A, ierr)) 367ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OBJECT, '-ua_vec_view', ierr)) 368ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OBJECT, '-u_vec_view', ierr)) 369ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OBJECT, '-a_vec_view', ierr)) 370280aadf6SBlaise Bourdin ! restrict to UA2 371d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2, ierr)) 372ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OBJECT, '-ua2_vec_view', ierr)) 373280aadf6SBlaise Bourdin 374e2739ba6SAlexis Marboeuf ! Getting Natural Vec 375d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr)) 376d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr)) 377280aadf6SBlaise Bourdin 378d8606c27SBarry Smith PetscCallA(VecView(U, viewer, ierr)) 379d8606c27SBarry Smith PetscCallA(VecView(A, viewer, ierr)) 380280aadf6SBlaise Bourdin 381280aadf6SBlaise Bourdin ! Saving U and Alpha in one shot. 382280aadf6SBlaise Bourdin ! For this, we need to cheat and change the Vec's name 383280aadf6SBlaise Bourdin ! Note that in the end we write variables one component at a time, 384280aadf6SBlaise Bourdin ! so that there is no real value in doing this 385d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA, 1_kPI, time, ierr)) 386d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec, ierr)) 387d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec, ierr)) 388dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U', ierr)) 389d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer, ierr)) 390e2739ba6SAlexis Marboeuf 391280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 392d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR, ierr)) 393d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer, ierr)) 394d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec, ierr)) 395d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm, ierr)) 396280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 397280aadf6SBlaise Bourdin write (IOBuffer, '("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 398280aadf6SBlaise Bourdin end if 399d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec, ierr)) 400280aadf6SBlaise Bourdin 401e2739ba6SAlexis Marboeuf ! same thing with the UA2 Vec obtained from the superDM 402d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec, ierr)) 403d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec, ierr)) 404dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U', ierr)) 405d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2, 2_kPI, time, ierr)) 406d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer, ierr)) 407e2739ba6SAlexis Marboeuf 408280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 409d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR, ierr)) 410d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer, ierr)) 411d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec, ierr)) 412d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm, ierr)) 413280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 414280aadf6SBlaise Bourdin write (IOBuffer, '("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 415280aadf6SBlaise Bourdin end if 416d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec, ierr)) 417280aadf6SBlaise Bourdin 418280aadf6SBlaise Bourdin ! Building and saving Sigma 419280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 420280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 421280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 422d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection, ierr)) 423d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord, ierr)) 424dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, 'Cell Sets', csIS, ierr)) 425dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dmS, 'Cell Sets', numCS, ierr)) 426ce78bad3SBarry Smith PetscCallA(ISGetIndices(csIS, csID, ierr)) 427280aadf6SBlaise Bourdin 428280aadf6SBlaise Bourdin do set = 1, numCS 429dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dmS, 'Cell Sets', csID(set), cellIS, ierr)) 430ce78bad3SBarry Smith PetscCallA(ISGetIndices(cellIS, cellID, ierr)) 431d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells, ierr)) 432280aadf6SBlaise Bourdin do cell = 1, numCells 433ce78bad3SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), PETSC_NULL_INTEGER, cval, ierr)) 434ce78bad3SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), PETSC_NULL_INTEGER, xyz, ierr)) 435280aadf6SBlaise Bourdin cval(1) = rank 436280aadf6SBlaise Bourdin cval(2) = csID(set) 437280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 438280aadf6SBlaise Bourdin do c = 1, size(xyz), sdim 439280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 440280aadf6SBlaise Bourdin end do 441280aadf6SBlaise Bourdin cval(3) = cval(3)*sdim/size(xyz) 442d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES, ierr)) 443ce78bad3SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), PETSC_NULL_INTEGER, cval, ierr)) 444ce78bad3SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), PETSC_NULL_INTEGER, xyz, ierr)) 445280aadf6SBlaise Bourdin end do 446ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(cellIS, cellID, ierr)) 447d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS, ierr)) 448280aadf6SBlaise Bourdin end do 449ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(csIS, csID, ierr)) 450d8606c27SBarry Smith PetscCallA(ISDestroy(csIS, ierr)) 451ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OBJECT, '-s_vec_view', ierr)) 452280aadf6SBlaise Bourdin 453280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 454d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS, 0_kPI, time, ierr)) 455d8606c27SBarry Smith PetscCallA(VecView(S, viewer, ierr)) 456280aadf6SBlaise Bourdin 457280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 458d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec, ierr)) 459d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR, ierr)) 460dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'Sigma', ierr)) 461d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer, ierr)) 462d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec, ierr)) 463d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY, norm, ierr)) 464280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 465280aadf6SBlaise Bourdin write (IOBuffer, '("Sigma ||Vin - Vout|| = ",ES12.5)') norm 466280aadf6SBlaise Bourdin end if 467d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec, ierr)) 468280aadf6SBlaise Bourdin 469d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2, ierr)) 470d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA, ierr)) 471d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S, ierr)) 472d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A, ierr)) 473d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U, ierr)) 474e2739ba6SAlexis Marboeuf PetscCallA(DMRestoreGlobalVector(pdm, X, ierr)) 475d8606c27SBarry Smith PetscCallA(DMDestroy(dmU, ierr)) 476d8606c27SBarry Smith PetscCallA(ISDestroy(isU, ierr)) 477d8606c27SBarry Smith PetscCallA(DMDestroy(dmA, ierr)) 478d8606c27SBarry Smith PetscCallA(ISDestroy(isA, ierr)) 479d8606c27SBarry Smith PetscCallA(DMDestroy(dmS, ierr)) 480d8606c27SBarry Smith PetscCallA(ISDestroy(isS, ierr)) 481d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA, ierr)) 482d8606c27SBarry Smith PetscCallA(ISDestroy(isUA, ierr)) 483d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2, ierr)) 484e2739ba6SAlexis Marboeuf PetscCallA(DMDestroy(pdm, ierr)) 485e2739ba6SAlexis Marboeuf if (numProc > 1) then 486d8606c27SBarry Smith PetscCallA(DMDestroy(dm, ierr)) 487e2739ba6SAlexis Marboeuf end if 488280aadf6SBlaise Bourdin 489280aadf6SBlaise Bourdin deallocate (pStartDepth) 490280aadf6SBlaise Bourdin deallocate (pEndDepth) 491280aadf6SBlaise Bourdin 492d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer, ierr)) 493d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 4942e8d78feSBlaise Bourdinend program ex62f90 495280aadf6SBlaise Bourdin 496280aadf6SBlaise Bourdin! /*TEST 497280aadf6SBlaise Bourdin! 498280aadf6SBlaise Bourdin! build: 499280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 500280aadf6SBlaise Bourdin! # 2D seq 501e2739ba6SAlexis Marboeuf! test: 502e2739ba6SAlexis Marboeuf! suffix: 0 5032e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1 504e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 505e2739ba6SAlexis Marboeuf! test: 506e2739ba6SAlexis Marboeuf! suffix: 1 5072e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1 508e2739ba6SAlexis Marboeuf! 509e2739ba6SAlexis Marboeuf! test: 510e2739ba6SAlexis Marboeuf! suffix: 2 5112e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1 512e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 513e2739ba6SAlexis Marboeuf! test: 514e2739ba6SAlexis Marboeuf! suffix: 3 5152e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2 516e2739ba6SAlexis Marboeuf! test: 517e2739ba6SAlexis Marboeuf! suffix: 4 5182e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2 519e2739ba6SAlexis Marboeuf! test: 520e2739ba6SAlexis Marboeuf! suffix: 5 5212e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2 522e2739ba6SAlexis Marboeuf! # 2D par 523e2739ba6SAlexis Marboeuf! test: 524e2739ba6SAlexis Marboeuf! suffix: 6 525e2739ba6SAlexis Marboeuf! nsize: 2 5262e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1 527e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 528e2739ba6SAlexis Marboeuf! test: 529e2739ba6SAlexis Marboeuf! suffix: 7 530e2739ba6SAlexis Marboeuf! nsize: 2 5312e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1 532e2739ba6SAlexis Marboeuf! test: 533e2739ba6SAlexis Marboeuf! suffix: 8 534e2739ba6SAlexis Marboeuf! nsize: 2 5352e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1 536e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 537e2739ba6SAlexis Marboeuf! test: 538e2739ba6SAlexis Marboeuf! suffix: 9 539e2739ba6SAlexis Marboeuf! nsize: 2 5402e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2 541e2739ba6SAlexis Marboeuf! test: 542e2739ba6SAlexis Marboeuf! suffix: 10 543e2739ba6SAlexis Marboeuf! nsize: 2 5442e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2 545e2739ba6SAlexis Marboeuf! test: 5469d2ffcd0SPierre Jolivet! # Something is now broken with parallel read/write for whatever shape H is 547e2739ba6SAlexis Marboeuf! TODO: broken 548e2739ba6SAlexis Marboeuf! suffix: 11 549e2739ba6SAlexis Marboeuf! nsize: 2 5502e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2 551e2739ba6SAlexis Marboeuf! #3d seq 552e2739ba6SAlexis Marboeuf! test: 553e2739ba6SAlexis Marboeuf! suffix: 12 5542e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1 555e2739ba6SAlexis Marboeuf! test: 556e2739ba6SAlexis Marboeuf! suffix: 13 5572e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1 558e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 559e2739ba6SAlexis Marboeuf! test: 560e2739ba6SAlexis Marboeuf! suffix: 14 5612e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2 562e2739ba6SAlexis Marboeuf! test: 563e2739ba6SAlexis Marboeuf! suffix: 15 5642e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2 565e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 566e2739ba6SAlexis Marboeuf! #3d par 567e2739ba6SAlexis Marboeuf! test: 568e2739ba6SAlexis Marboeuf! suffix: 16 569e2739ba6SAlexis Marboeuf! nsize: 2 5702e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1 571e2739ba6SAlexis Marboeuf! test: 572e2739ba6SAlexis Marboeuf! suffix: 17 573e2739ba6SAlexis Marboeuf! nsize: 2 5742e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1 575e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 576e2739ba6SAlexis Marboeuf! test: 577e2739ba6SAlexis Marboeuf! suffix: 18 578e2739ba6SAlexis Marboeuf! nsize: 2 5792e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2 580e2739ba6SAlexis Marboeuf! test: 581e2739ba6SAlexis Marboeuf! suffix: 19 582e2739ba6SAlexis Marboeuf! nsize: 2 5832e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2 584e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 585e2447916SStefano Zampini! 586280aadf6SBlaise Bourdin! TEST*/ 587