1ce78bad3SBarry Smith#include "petsc/finclude/petscdmplex.h" 2c5e229c2SMartin Diehlprogram ex62f90 3ce78bad3SBarry Smith use petscdmplex 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 17280aadf6SBlaise Bourdin type(tPetscSection) :: section 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 2802c639afSMartin Diehl integer :: i 290a5cf5d8SBlaise Bourdin PetscInt :: sdim, d, pStart, pEnd, p, numCS, set, j 30280aadf6SBlaise Bourdin PetscMPIInt :: rank, numProc 31280aadf6SBlaise Bourdin PetscBool :: flg 32280aadf6SBlaise Bourdin PetscErrorCode :: ierr 33*b06eb4cdSBarry Smith MPIU_Comm :: comm 34280aadf6SBlaise Bourdin type(tPetscViewer) :: viewer 35280aadf6SBlaise Bourdin 3602c639afSMartin Diehl character(len=MXSTLN) :: sJunk 37280aadf6SBlaise Bourdin PetscInt :: numstep = 3, step 3802c639afSMartin Diehl integer :: numNodalVar, numZonalVar 3949c89c76SBlaise Bourdin character(len=MXNAME), dimension(4) :: nodalVarName2D = ["U_x ", & 4049c89c76SBlaise Bourdin "U_y ", & 4149c89c76SBlaise Bourdin "Alpha", & 4249c89c76SBlaise Bourdin "Beta "] 4349c89c76SBlaise Bourdin character(len=MXNAME), dimension(3) :: zonalVarName2D = ["Sigma_11", & 4449c89c76SBlaise Bourdin "Sigma_12", & 4549c89c76SBlaise Bourdin "Sigma_22"] 4649c89c76SBlaise Bourdin character(len=MXNAME), dimension(5) :: nodalVarName3D = ["U_x ", & 4749c89c76SBlaise Bourdin "U_y ", & 4849c89c76SBlaise Bourdin "U_z ", & 4949c89c76SBlaise Bourdin "Alpha", & 5049c89c76SBlaise Bourdin "Beta "] 5149c89c76SBlaise Bourdin character(len=MXNAME), dimension(6) :: zonalVarName3D = ["Sigma_11", & 5249c89c76SBlaise Bourdin "Sigma_22", & 5349c89c76SBlaise Bourdin "Sigma_33", & 5449c89c76SBlaise Bourdin "Sigma_23", & 5549c89c76SBlaise Bourdin "Sigma_13", & 5649c89c76SBlaise Bourdin "Sigma_12"] 57280aadf6SBlaise Bourdin logical, dimension(:, :), pointer :: truthtable 58280aadf6SBlaise Bourdin 59280aadf6SBlaise Bourdin type(tIS) :: cellIS 60280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: cellID 61280aadf6SBlaise Bourdin PetscInt :: numCells, cell, closureSize 62280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: closureA, closure 63280aadf6SBlaise Bourdin 64280aadf6SBlaise Bourdin type(tPetscSection) :: sectionUA, coordSection 65280aadf6SBlaise Bourdin type(tVec) :: UALoc, coord 66280aadf6SBlaise Bourdin PetscScalar, dimension(:), pointer :: cval, xyz 67280aadf6SBlaise Bourdin PetscInt :: dofUA, offUA, c 68280aadf6SBlaise Bourdin 69280aadf6SBlaise Bourdin ! dof layout ordered by increasing height in the DAG: cell, face, edge, vertex 70280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofS2D = [0, 0, 3] 71280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP1Tri = [2, 0, 0] 72280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP1Tri = [1, 0, 0] 73280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP2Tri = [2, 2, 0] 74280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP2Tri = [1, 1, 0] 75280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP1Quad = [2, 0, 0] 76280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP1Quad = [1, 0, 0] 77280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofUP2Quad = [2, 2, 2] 78280aadf6SBlaise Bourdin PetscInt, dimension(3), target :: dofAP2Quad = [1, 1, 1] 79280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofS3D = [0, 0, 0, 6] 80280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP1Tet = [3, 0, 0, 0] 81280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP1Tet = [1, 0, 0, 0] 82280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP2Tet = [3, 3, 0, 0] 83280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP2Tet = [1, 1, 0, 0] 84280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP1Hex = [3, 0, 0, 0] 85280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP1Hex = [1, 0, 0, 0] 86280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofUP2Hex = [3, 3, 3, 3] 87280aadf6SBlaise Bourdin PetscInt, dimension(4), target :: dofAP2Hex = [1, 1, 1, 1] 88280aadf6SBlaise Bourdin PetscInt, dimension(:), pointer :: dofU, dofA, dofS 89280aadf6SBlaise Bourdin 90280aadf6SBlaise Bourdin type(tPetscSF) :: migrationSF 91280aadf6SBlaise Bourdin PetscPartitioner :: part 92f51a5268SBarry Smith PetscLayout :: layout 93280aadf6SBlaise Bourdin 94280aadf6SBlaise Bourdin type(tVec) :: tmpVec 95280aadf6SBlaise Bourdin PetscReal :: norm 96280aadf6SBlaise Bourdin PetscReal :: time = 1.234_kPR 97280aadf6SBlaise Bourdin 98e2739ba6SAlexis Marboeuf PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr)) 99e2739ba6SAlexis Marboeuf if (ierr /= 0) then 100e2739ba6SAlexis Marboeuf print *, 'Unable to initialize PETSc' 101e2739ba6SAlexis Marboeuf stop 102e2739ba6SAlexis Marboeuf end if 103280aadf6SBlaise Bourdin 104d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) 105d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, numProc, ierr)) 106dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr)) 107dcb3e689SBarry Smith PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i <input file name>') 108dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-o', ofilename, flg, ierr)) 109dcb3e689SBarry Smith PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing output file name -o <output file name>') 110dcb3e689SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-order', order, flg, ierr)) 111280aadf6SBlaise Bourdin if ((order > 2) .or. (order < 1)) then 112280aadf6SBlaise Bourdin write (IOBuffer, '("Unsupported polynomial order ", I2, " not in [1,2]")') order 113de401611SBlaise Bourdin SETERRA(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, IOBuffer) 114280aadf6SBlaise Bourdin end if 115280aadf6SBlaise Bourdin 116280aadf6SBlaise Bourdin ! Read the mesh in any supported format 117d8606c27SBarry Smith PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_NULL_CHARACTER, PETSC_TRUE, dm, ierr)) 118ccfd86f1SBarry Smith PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE, ierr)) 119ccfd86f1SBarry Smith PetscCallA(DMSetFromOptions(dm, ierr)) 120d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim, ierr)) 121ccfd86f1SBarry Smith PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr)) 122280aadf6SBlaise Bourdin 123280aadf6SBlaise Bourdin ! Create the exodus result file 124280aadf6SBlaise Bourdin 125a5b23f4aSJose E. Roman ! enable exodus debugging information 126d8606c27SBarry Smith PetscCallA(exopts(EXVRBS + EXDEBG, ierr)) 127280aadf6SBlaise Bourdin ! Create the exodus file 128d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_WRITE, viewer, ierr)) 129280aadf6SBlaise Bourdin ! The long way would be 130280aadf6SBlaise Bourdin ! 131d8606c27SBarry Smith ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr)) 132d8606c27SBarry Smith ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr)) 133d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr)) 134d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr)) 135280aadf6SBlaise Bourdin 136280aadf6SBlaise Bourdin ! set the mesh order 137d8606c27SBarry Smith PetscCallA(PetscViewerExodusIISetOrder(viewer, order, ierr)) 138d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr)) 139280aadf6SBlaise Bourdin ! 140280aadf6SBlaise Bourdin ! Notice how the exodus file is actually NOT open at this point (exoid is -1) 141d5b43468SJose E. Roman ! Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to 142280aadf6SBlaise Bourdin ! write the geometry (the DM), which can only be done on a brand new file. 143280aadf6SBlaise Bourdin ! 144280aadf6SBlaise Bourdin 145280aadf6SBlaise Bourdin ! Save the geometry to the file, erasing all previous content 146d8606c27SBarry Smith PetscCallA(DMView(dm, viewer, ierr)) 147d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr)) 148280aadf6SBlaise Bourdin ! 149280aadf6SBlaise Bourdin ! Note how the exodus file is now open 150280aadf6SBlaise Bourdin ! 151dcb3e689SBarry Smith ! 'Format' the exodus result file, i.e. allocate space for nodal and zonal variables 152280aadf6SBlaise Bourdin select case (sdim) 153280aadf6SBlaise Bourdin case (2) 154280aadf6SBlaise Bourdin numNodalVar = 4 15549c89c76SBlaise Bourdin numZonalVar = 3 15649c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariable(viewer, numZonalVar, ierr)) 15749c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariable(viewer, numNodalVar, ierr)) 15849c89c76SBlaise Bourdin do i = 1, numZonalVar 15949c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariableName(viewer, i - 1, zonalVarName2D(i), ierr)) 16049c89c76SBlaise Bourdin end do 16149c89c76SBlaise Bourdin do i = 1, numNodalVar 16249c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariableName(viewer, i - 1, nodalVarName2D(i), ierr)) 16349c89c76SBlaise Bourdin end do 16449c89c76SBlaise Bourdin case (3) 16549c89c76SBlaise Bourdin numNodalVar = 5 166280aadf6SBlaise Bourdin numZonalVar = 6 16749c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariable(viewer, numZonalVar, ierr)) 16849c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariable(viewer, numNodalVar, ierr)) 16949c89c76SBlaise Bourdin do i = 1, numZonalVar 17049c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariableName(viewer, i - 1, zonalVarName3D(i), ierr)) 17149c89c76SBlaise Bourdin end do 17249c89c76SBlaise Bourdin do i = 1, numNodalVar 17349c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariableName(viewer, i - 1, nodalVarName3D(i), ierr)) 17449c89c76SBlaise Bourdin end do 175280aadf6SBlaise Bourdin case default 176280aadf6SBlaise Bourdin write (IOBuffer, '("No layout for dimension ",I2)') sdim 177280aadf6SBlaise Bourdin end select 17849c89c76SBlaise Bourdin PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr)) 179280aadf6SBlaise Bourdin 180caff39ffSPierre Jolivet ! An ExodusII truth table specifies which fields are saved at which time step 181280aadf6SBlaise Bourdin ! It speeds up I/O but reserving space for fields in the file ahead of time. 18249c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIIGetId(viewer, exoid, ierr)) 18349c89c76SBlaise Bourdin PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK, numCS, PETSC_NULL_REAL, sjunk, ierr)) 184280aadf6SBlaise Bourdin allocate (truthtable(numCS, numZonalVar)) 185280aadf6SBlaise Bourdin truthtable = .true. 186d8606c27SBarry Smith PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr)) 187280aadf6SBlaise Bourdin deallocate (truthtable) 188280aadf6SBlaise Bourdin 189280aadf6SBlaise Bourdin ! Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 190280aadf6SBlaise Bourdin do step = 1, numstep 19102c639afSMartin Diehl PetscCallA(exptim(exoid, step, real(step, kind=kPR), ierr)) 192280aadf6SBlaise Bourdin end do 193280aadf6SBlaise Bourdin 194d8606c27SBarry Smith PetscCallA(PetscObjectGetComm(dm, comm, ierr)) 195d8606c27SBarry Smith PetscCallA(PetscSectionCreate(comm, section, ierr)) 196d8606c27SBarry Smith PetscCallA(PetscSectionSetNumFields(section, 3_kPI, ierr)) 197dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, 'U', ierr)) 198dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, 'Alpha', ierr)) 199dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldS, 'Sigma', ierr)) 200d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr)) 201d8606c27SBarry Smith PetscCallA(PetscSectionSetChart(section, pStart, pEnd, ierr)) 202280aadf6SBlaise Bourdin 203280aadf6SBlaise Bourdin allocate (pStartDepth(sdim + 1)) 204280aadf6SBlaise Bourdin allocate (pEndDepth(sdim + 1)) 205280aadf6SBlaise Bourdin do d = 1, sdim + 1 206d8606c27SBarry Smith PetscCallA(DMPlexGetDepthStratum(dm, d - 1, pStartDepth(d), pEndDepth(d), ierr)) 207280aadf6SBlaise Bourdin end do 208280aadf6SBlaise Bourdin 209280aadf6SBlaise Bourdin ! Vector field U, Scalar field Alpha, Tensor field Sigma 210ccfd86f1SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim, ierr)) 211ccfd86f1SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI, ierr)) 212ccfd86f1SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim + 1)/2, ierr)) 213280aadf6SBlaise Bourdin 214280aadf6SBlaise Bourdin ! Going through cell sets then cells, and setting up storage for the sections 215dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dm, 'Cell Sets', numCS, ierr)) 216dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dm, 'Cell Sets', csIS, ierr)) 217ce78bad3SBarry Smith PetscCallA(ISGetIndices(csIS, csID, ierr)) 218280aadf6SBlaise Bourdin do set = 1, numCS 2198f2079feSBlaise Bourdin !!! This is pointless but will avoid a warning with gfortran and -Werror=maybe-uninitialized 2208f2079feSBlaise Bourdin nullify (dofA) 2218f2079feSBlaise Bourdin nullify (dofU) 2228f2079feSBlaise Bourdin nullify (dofS) 223dcb3e689SBarry Smith PetscCallA(DMGetStratumSize(dm, 'Cell Sets', csID(set), numCells, ierr)) 224dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dm, 'Cell Sets', csID(set), cellIS, ierr)) 225280aadf6SBlaise Bourdin if (numCells > 0) then 226280aadf6SBlaise Bourdin select case (sdim) 227280aadf6SBlaise Bourdin case (2) 228280aadf6SBlaise Bourdin dofs => dofS2D 229280aadf6SBlaise Bourdin case (3) 230280aadf6SBlaise Bourdin dofs => dofS3D 231280aadf6SBlaise Bourdin case default 232280aadf6SBlaise Bourdin write (IOBuffer, '("No layout for dimension ",I2)') sdim 233de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, IOBuffer) 234280aadf6SBlaise Bourdin end select ! sdim 235280aadf6SBlaise Bourdin 236280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 237280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 238ce78bad3SBarry Smith PetscCallA(ISGetIndices(cellIS, cellID, ierr)) 239280aadf6SBlaise Bourdin nullify (closureA) 240ce78bad3SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(1), PETSC_TRUE, PETSC_NULL_INTEGER, closureA, ierr)) 241280aadf6SBlaise Bourdin select case (size(closureA)/2) 242280aadf6SBlaise Bourdin case (7) ! Tri 243280aadf6SBlaise Bourdin if (order == 1) then 244280aadf6SBlaise Bourdin dofU => dofUP1Tri 245280aadf6SBlaise Bourdin dofA => dofAP1Tri 246280aadf6SBlaise Bourdin else 247280aadf6SBlaise Bourdin dofU => dofUP2Tri 248280aadf6SBlaise Bourdin dofA => dofAP2Tri 249280aadf6SBlaise Bourdin end if 250280aadf6SBlaise Bourdin case (9) ! Quad 251280aadf6SBlaise Bourdin if (order == 1) then 252280aadf6SBlaise Bourdin dofU => dofUP1Quad 253280aadf6SBlaise Bourdin dofA => dofAP1Quad 254280aadf6SBlaise Bourdin else 255280aadf6SBlaise Bourdin dofU => dofUP2Quad 256280aadf6SBlaise Bourdin dofA => dofAP2Quad 257280aadf6SBlaise Bourdin end if 258280aadf6SBlaise Bourdin case (15) ! Tet 259280aadf6SBlaise Bourdin if (order == 1) then 260280aadf6SBlaise Bourdin dofU => dofUP1Tet 261280aadf6SBlaise Bourdin dofA => dofAP1Tet 262280aadf6SBlaise Bourdin else 263280aadf6SBlaise Bourdin dofU => dofUP2Tet 264280aadf6SBlaise Bourdin dofA => dofAP2Tet 265280aadf6SBlaise Bourdin end if 266280aadf6SBlaise Bourdin case (27) ! Hex 267280aadf6SBlaise Bourdin if (order == 1) then 268280aadf6SBlaise Bourdin dofU => dofUP1Hex 269280aadf6SBlaise Bourdin dofA => dofAP1Hex 270280aadf6SBlaise Bourdin else 271280aadf6SBlaise Bourdin dofU => dofUP2Hex 272280aadf6SBlaise Bourdin dofA => dofAP2Hex 273280aadf6SBlaise Bourdin end if 274280aadf6SBlaise Bourdin case default 275280aadf6SBlaise Bourdin write (IOBuffer, '("Unknown element with closure size ",I2)') size(closureA)/2 276de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, IOBuffer) 277280aadf6SBlaise Bourdin end select 278ce78bad3SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE, PETSC_NULL_INTEGER, closureA, ierr)) 279280aadf6SBlaise Bourdin do cell = 1, numCells! 280280aadf6SBlaise Bourdin nullify (closure) 281ce78bad3SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, PETSC_NULL_INTEGER, closure, ierr)) 282280aadf6SBlaise Bourdin do p = 1, size(closure), 2 283280aadf6SBlaise Bourdin ! find the depth of p 284280aadf6SBlaise Bourdin do d = 1, sdim + 1 285280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 286d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d) + dofA(d) + dofS(d), ierr)) 287d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d), ierr)) 288d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d), ierr)) 289d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d), ierr)) 290280aadf6SBlaise Bourdin end if ! closure(p) 291280aadf6SBlaise Bourdin end do ! d 292280aadf6SBlaise Bourdin end do ! p 293ce78bad3SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, PETSC_NULL_INTEGER, closure, ierr)) 294280aadf6SBlaise Bourdin end do ! cell 295ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(cellIS, cellID, ierr)) 296d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS, ierr)) 297280aadf6SBlaise Bourdin end if ! numCells 298280aadf6SBlaise Bourdin end do ! set 299ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(csIS, csID, ierr)) 300d8606c27SBarry Smith PetscCallA(ISDestroy(csIS, ierr)) 301d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section, ierr)) 302d8606c27SBarry Smith PetscCallA(DMSetLocalSection(dm, section, ierr)) 303ce78bad3SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_OBJECT, '-dm_section_view', ierr)) 304f51a5268SBarry Smith PetscCallA(PetscSectionGetValueLayout(PETSC_COMM_WORLD, section, layout, ierr)) 305f51a5268SBarry Smith PetscCallA(PetscLayoutDestroy(layout, ierr)) 306d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section, ierr)) 307280aadf6SBlaise Bourdin 308d8606c27SBarry Smith PetscCallA(DMSetUseNatural(dm, PETSC_TRUE, ierr)) 309d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm, part, ierr)) 310d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part, ierr)) 311d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm, 0_kPI, migrationSF, pdm, ierr)) 312280aadf6SBlaise Bourdin 313280aadf6SBlaise Bourdin if (numProc > 1) then 314d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm, migrationSF, ierr)) 315d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF, ierr)) 316e2739ba6SAlexis Marboeuf else 317e2739ba6SAlexis Marboeuf pdm = dm 318280aadf6SBlaise Bourdin end if 319ce78bad3SBarry Smith PetscCallA(DMViewFromOptions(pdm, PETSC_NULL_OBJECT, '-dm_view', ierr)) 320280aadf6SBlaise Bourdin 321280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 322ce78bad3SBarry Smith PetscCallA(DMCreateSubDM(pdm, 1_kPI, [fieldU], isU, dmU, ierr)) 323ce78bad3SBarry Smith PetscCallA(DMCreateSubDM(pdm, 1_kPI, [fieldA], isA, dmA, ierr)) 324ce78bad3SBarry Smith PetscCallA(DMCreateSubDM(pdm, 1_kPI, [fieldS], isS, dmS, ierr)) 325e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA, ierr)) 326280aadf6SBlaise Bourdin 327280aadf6SBlaise Bourdin !Create the exodus result file 328280aadf6SBlaise Bourdin allocate (dmList(2)) 329ccfd86f1SBarry Smith dmList(1) = dmU 330ccfd86f1SBarry Smith dmList(2) = dmA 331ce78bad3SBarry Smith PetscCallA(DMCreateSuperDM(dmList, 2_kPI, PETSC_NULL_IS_POINTER, dmUA2, ierr)) 332280aadf6SBlaise Bourdin deallocate (dmList) 333280aadf6SBlaise Bourdin 334e2739ba6SAlexis Marboeuf PetscCallA(DMGetGlobalVector(pdm, X, ierr)) 335d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U, ierr)) 336d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A, ierr)) 337d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S, ierr)) 338d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA, ierr)) 339d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2, ierr)) 340280aadf6SBlaise Bourdin 341dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(U, 'U', ierr)) 342dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(A, 'Alpha', ierr)) 343dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(S, 'Sigma', ierr)) 344dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA, 'UAlpha', ierr)) 345dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA2, 'UAlpha2', ierr)) 346d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR, ierr)) 347280aadf6SBlaise Bourdin 348280aadf6SBlaise 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 */ 349d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA, ierr)) 350d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc, ierr)) 351ce78bad3SBarry Smith PetscCallA(VecGetArray(UALoc, cval, ierr)) 352d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection, ierr)) 353d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord, ierr)) 354d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd, ierr)) 355280aadf6SBlaise Bourdin 356280aadf6SBlaise Bourdin do p = pStart, pEnd - 1 357d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA, ierr)) 358280aadf6SBlaise Bourdin if (dofUA > 0) then 359d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA, ierr)) 360ce78bad3SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, PETSC_NULL_INTEGER, xyz, ierr)) 361280aadf6SBlaise Bourdin closureSize = size(xyz) 362280aadf6SBlaise Bourdin do i = 1, sdim 363280aadf6SBlaise Bourdin do j = 0, closureSize - 1, sdim 364280aadf6SBlaise Bourdin cval(offUA + i) = cval(offUA + i) + xyz(j/sdim + i) 365280aadf6SBlaise Bourdin end do 366ccfd86f1SBarry Smith cval(offUA + i) = cval(offUA + i)*sdim/closureSize 367280aadf6SBlaise Bourdin cval(offUA + sdim + 1) = cval(offUA + sdim + 1) + cval(offUA + i)**2 368280aadf6SBlaise Bourdin end do 369ce78bad3SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, PETSC_NULL_INTEGER, xyz, ierr)) 370280aadf6SBlaise Bourdin end if 371280aadf6SBlaise Bourdin end do 372280aadf6SBlaise Bourdin 373ce78bad3SBarry Smith PetscCallA(VecRestoreArray(UALoc, cval, ierr)) 374d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA, ierr)) 375d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA, ierr)) 376d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc, ierr)) 377280aadf6SBlaise Bourdin 378280aadf6SBlaise Bourdin !Update X 379d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA, ierr)) 380280aadf6SBlaise Bourdin ! Restrict to U and Alpha 381d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U, ierr)) 382d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A, ierr)) 383ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OBJECT, '-ua_vec_view', ierr)) 384ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OBJECT, '-u_vec_view', ierr)) 385ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OBJECT, '-a_vec_view', ierr)) 386280aadf6SBlaise Bourdin ! restrict to UA2 387d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2, ierr)) 388ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OBJECT, '-ua2_vec_view', ierr)) 389280aadf6SBlaise Bourdin 390e2739ba6SAlexis Marboeuf ! Getting Natural Vec 391d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr)) 392d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr)) 393280aadf6SBlaise Bourdin 394d8606c27SBarry Smith PetscCallA(VecView(U, viewer, ierr)) 395d8606c27SBarry Smith PetscCallA(VecView(A, viewer, ierr)) 396280aadf6SBlaise Bourdin 397280aadf6SBlaise Bourdin !Saving U and Alpha in one shot. 398280aadf6SBlaise Bourdin !For this, we need to cheat and change the Vec's name 399280aadf6SBlaise Bourdin !Note that in the end we write variables one component at a time, 400280aadf6SBlaise Bourdin !so that there is no real value in doing this 401d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA, 1_kPI, time, ierr)) 402d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec, ierr)) 403d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec, ierr)) 404dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U', ierr)) 405d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer, ierr)) 406280aadf6SBlaise Bourdin 407280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 408d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR, ierr)) 409d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer, ierr)) 410d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec, ierr)) 411d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm, ierr)) 412280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 413280aadf6SBlaise Bourdin write (IOBuffer, '("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 414280aadf6SBlaise Bourdin end if 415d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec, ierr)) 416280aadf6SBlaise Bourdin 417280aadf6SBlaise Bourdin ! same thing with the UA2 Vec obtained from the superDM 418d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec, ierr)) 419d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec, ierr)) 420dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U', ierr)) 421d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2, 2_kPI, time, ierr)) 422d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer, ierr)) 423280aadf6SBlaise Bourdin 424280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 425d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR, ierr)) 426d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer, ierr)) 427d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec, ierr)) 428d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm, ierr)) 429280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 430280aadf6SBlaise Bourdin write (IOBuffer, '("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 431280aadf6SBlaise Bourdin end if 432d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec, ierr)) 433280aadf6SBlaise Bourdin 434280aadf6SBlaise Bourdin ! Building and saving Sigma 435280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 436280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 437280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 438d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection, ierr)) 439d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord, ierr)) 440dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, 'Cell Sets', csIS, ierr)) 441dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dmS, 'Cell Sets', numCS, ierr)) 442ce78bad3SBarry Smith PetscCallA(ISGetIndices(csIS, csID, ierr)) 443280aadf6SBlaise Bourdin 444280aadf6SBlaise Bourdin do set = 1, numCS 445dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dmS, 'Cell Sets', csID(set), cellIS, ierr)) 446ce78bad3SBarry Smith PetscCallA(ISGetIndices(cellIS, cellID, ierr)) 447d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells, ierr)) 448280aadf6SBlaise Bourdin do cell = 1, numCells 449ce78bad3SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), PETSC_NULL_INTEGER, cval, ierr)) 450ce78bad3SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), PETSC_NULL_INTEGER, xyz, ierr)) 451280aadf6SBlaise Bourdin cval(1) = rank 452280aadf6SBlaise Bourdin cval(2) = csID(set) 453280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 454280aadf6SBlaise Bourdin do c = 1, size(xyz), sdim 455280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 456280aadf6SBlaise Bourdin end do 457280aadf6SBlaise Bourdin cval(3) = cval(3)*sdim/size(xyz) 458d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES, ierr)) 459ce78bad3SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), PETSC_NULL_INTEGER, cval, ierr)) 460ce78bad3SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), PETSC_NULL_INTEGER, xyz, ierr)) 461280aadf6SBlaise Bourdin end do 462ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(cellIS, cellID, ierr)) 463d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS, ierr)) 464280aadf6SBlaise Bourdin end do 465ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(csIS, csID, ierr)) 466d8606c27SBarry Smith PetscCallA(ISDestroy(csIS, ierr)) 467ce78bad3SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OBJECT, '-s_vec_view', ierr)) 468280aadf6SBlaise Bourdin 469280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 470d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS, 0_kPI, time, ierr)) 471d8606c27SBarry Smith PetscCallA(VecView(S, viewer, ierr)) 472280aadf6SBlaise Bourdin 473280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 474d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec, ierr)) 475d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR, ierr)) 476dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'Sigma', ierr)) 477d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer, ierr)) 478d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec, ierr)) 479d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY, norm, ierr)) 480280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 481280aadf6SBlaise Bourdin write (IOBuffer, '("Sigma ||Vin - Vout|| = ",ES12.5)') norm 482280aadf6SBlaise Bourdin end if 483d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec, ierr)) 484280aadf6SBlaise Bourdin 485d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2, ierr)) 486d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA, ierr)) 487d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S, ierr)) 488d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A, ierr)) 489d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U, ierr)) 490e2739ba6SAlexis Marboeuf PetscCallA(DMRestoreGlobalVector(pdm, X, ierr)) 491d8606c27SBarry Smith PetscCallA(DMDestroy(dmU, ierr)) 492d8606c27SBarry Smith PetscCallA(ISDestroy(isU, ierr)) 493d8606c27SBarry Smith PetscCallA(DMDestroy(dmA, ierr)) 494d8606c27SBarry Smith PetscCallA(ISDestroy(isA, ierr)) 495d8606c27SBarry Smith PetscCallA(DMDestroy(dmS, ierr)) 496d8606c27SBarry Smith PetscCallA(ISDestroy(isS, ierr)) 497d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA, ierr)) 498d8606c27SBarry Smith PetscCallA(ISDestroy(isUA, ierr)) 499d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2, ierr)) 500e2739ba6SAlexis Marboeuf PetscCallA(DMDestroy(pdm, ierr)) 501e2739ba6SAlexis Marboeuf if (numProc > 1) then 502d8606c27SBarry Smith PetscCallA(DMDestroy(dm, ierr)) 503e2739ba6SAlexis Marboeuf end if 504280aadf6SBlaise Bourdin 505280aadf6SBlaise Bourdin deallocate (pStartDepth) 506280aadf6SBlaise Bourdin deallocate (pEndDepth) 507280aadf6SBlaise Bourdin 508d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer, ierr)) 509d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 5102e8d78feSBlaise Bourdinend program ex62f90 511280aadf6SBlaise Bourdin 512280aadf6SBlaise Bourdin! /*TEST 513280aadf6SBlaise Bourdin! 514280aadf6SBlaise Bourdin! build: 515280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 516280aadf6SBlaise Bourdin! # 2D seq 517e2739ba6SAlexis Marboeuf! test: 518e2739ba6SAlexis Marboeuf! suffix: 0 519e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 520e2739ba6SAlexis 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 521e2739ba6SAlexis Marboeuf! test: 522e2739ba6SAlexis Marboeuf! suffix: 1 523e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 524e2739ba6SAlexis Marboeuf! 525e2739ba6SAlexis Marboeuf! test: 526e2739ba6SAlexis Marboeuf! suffix: 2 527e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 528e2739ba6SAlexis 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 529e2739ba6SAlexis Marboeuf! test: 530e2739ba6SAlexis Marboeuf! suffix: 3 531e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 532e2739ba6SAlexis Marboeuf! test: 533e2739ba6SAlexis Marboeuf! suffix: 4 534e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 535e2739ba6SAlexis Marboeuf! test: 536e2739ba6SAlexis Marboeuf! suffix: 5 537e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 538e2739ba6SAlexis Marboeuf! # 2D par 539e2739ba6SAlexis Marboeuf! test: 540e2739ba6SAlexis Marboeuf! suffix: 6 541e2739ba6SAlexis Marboeuf! nsize: 2 542e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 543e2739ba6SAlexis 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 544e2739ba6SAlexis Marboeuf! test: 545e2739ba6SAlexis Marboeuf! suffix: 7 546e2739ba6SAlexis Marboeuf! nsize: 2 547e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 548e2739ba6SAlexis Marboeuf! test: 549e2739ba6SAlexis Marboeuf! suffix: 8 550e2739ba6SAlexis Marboeuf! nsize: 2 551e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 552e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 553e2739ba6SAlexis Marboeuf! test: 554e2739ba6SAlexis Marboeuf! suffix: 9 555e2739ba6SAlexis Marboeuf! nsize: 2 556e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 557e2739ba6SAlexis Marboeuf! test: 558e2739ba6SAlexis Marboeuf! suffix: 10 559e2739ba6SAlexis Marboeuf! nsize: 2 560e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 561e2739ba6SAlexis Marboeuf! test: 5629c89aa79SPierre Jolivet! # Something is now broken with parallel read/write for whatever shape H is 563e2739ba6SAlexis Marboeuf! TODO: broken 564e2739ba6SAlexis Marboeuf! suffix: 11 565e2739ba6SAlexis Marboeuf! nsize: 2 566e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 567e2739ba6SAlexis Marboeuf 568e2739ba6SAlexis Marboeuf! #3d seq 569e2739ba6SAlexis Marboeuf! test: 570e2739ba6SAlexis Marboeuf! suffix: 12 571e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 572e2739ba6SAlexis Marboeuf! test: 573e2739ba6SAlexis Marboeuf! suffix: 13 574e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_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: 14 578e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 579e2739ba6SAlexis Marboeuf! test: 580e2739ba6SAlexis Marboeuf! suffix: 15 581e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 582e2739ba6SAlexis 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 583e2739ba6SAlexis Marboeuf! #3d par 584e2739ba6SAlexis Marboeuf! test: 585e2739ba6SAlexis Marboeuf! suffix: 16 586e2739ba6SAlexis Marboeuf! nsize: 2 587e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 588e2739ba6SAlexis Marboeuf! test: 589e2739ba6SAlexis Marboeuf! suffix: 17 590e2739ba6SAlexis Marboeuf! nsize: 2 591e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1 592e2739ba6SAlexis 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 593e2739ba6SAlexis Marboeuf! test: 594e2739ba6SAlexis Marboeuf! suffix: 18 595e2739ba6SAlexis Marboeuf! nsize: 2 596e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 597e2739ba6SAlexis Marboeuf! test: 598e2739ba6SAlexis Marboeuf! suffix: 19 599e2739ba6SAlexis Marboeuf! nsize: 2 600e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 601e2739ba6SAlexis 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 602e2739ba6SAlexis Marboeuf 603280aadf6SBlaise Bourdin! TEST*/ 604