12e8d78feSBlaise Bourdinprogram ex62f90 2280aadf6SBlaise Bourdin#include "petsc/finclude/petsc.h" 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 17280aadf6SBlaise Bourdin type(tPetscSection) :: section 18280aadf6SBlaise Bourdin PetscInt,dimension(1) :: fieldU = [0] 19280aadf6SBlaise Bourdin PetscInt,dimension(1) :: fieldA = [2] 20280aadf6SBlaise Bourdin PetscInt,dimension(1) :: 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 32280aadf6SBlaise Bourdin MPI_Comm :: comm 33280aadf6SBlaise Bourdin type(tPetscViewer) :: viewer 34280aadf6SBlaise Bourdin 35280aadf6SBlaise Bourdin Character(len=MXSTLN) :: sJunk 36280aadf6SBlaise Bourdin PetscInt :: numstep = 3, step 37280aadf6SBlaise Bourdin PetscInt :: numNodalVar,numZonalVar 38*49c89c76SBlaise Bourdin character(len=MXNAME),dimension(4) :: nodalVarName2D = ["U_x ", & 39*49c89c76SBlaise Bourdin "U_y ", & 40*49c89c76SBlaise Bourdin "Alpha", & 41*49c89c76SBlaise Bourdin "Beta "] 42*49c89c76SBlaise Bourdin character(len=MXNAME),dimension(3) :: zonalVarName2D = ["Sigma_11", & 43*49c89c76SBlaise Bourdin "Sigma_12", & 44*49c89c76SBlaise Bourdin "Sigma_22"] 45*49c89c76SBlaise Bourdin character(len=MXNAME),dimension(5) :: nodalVarName3D = ["U_x ", & 46*49c89c76SBlaise Bourdin "U_y ", & 47*49c89c76SBlaise Bourdin "U_z ", & 48*49c89c76SBlaise Bourdin "Alpha", & 49*49c89c76SBlaise Bourdin "Beta "] 50*49c89c76SBlaise Bourdin character(len=MXNAME),dimension(6) :: zonalVarName3D = ["Sigma_11", & 51*49c89c76SBlaise Bourdin "Sigma_22", & 52*49c89c76SBlaise Bourdin "Sigma_33", & 53*49c89c76SBlaise Bourdin "Sigma_23", & 54*49c89c76SBlaise Bourdin "Sigma_13", & 55*49c89c76SBlaise Bourdin "Sigma_12"] 56280aadf6SBlaise Bourdin logical,dimension(:,:),pointer :: truthtable 57280aadf6SBlaise Bourdin 58280aadf6SBlaise Bourdin type(tIS) :: cellIS 59280aadf6SBlaise Bourdin PetscInt,dimension(:),pointer :: cellID 60280aadf6SBlaise Bourdin PetscInt :: numCells, cell, closureSize 61280aadf6SBlaise Bourdin PetscInt,dimension(:),pointer :: closureA,closure 62280aadf6SBlaise Bourdin 63280aadf6SBlaise Bourdin type(tPetscSection) :: sectionUA,coordSection 64280aadf6SBlaise Bourdin type(tVec) :: UALoc,coord 65280aadf6SBlaise Bourdin PetscScalar,dimension(:),pointer :: cval,xyz 66280aadf6SBlaise Bourdin PetscInt :: dofUA,offUA,c 67280aadf6SBlaise Bourdin 68280aadf6SBlaise Bourdin ! dof layout ordered by increasing height in the DAG: cell, face, edge, vertex 69280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofS2D = [0, 0, 3] 70280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofUP1Tri = [2, 0, 0] 71280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofAP1Tri = [1, 0, 0] 72280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofUP2Tri = [2, 2, 0] 73280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofAP2Tri = [1, 1, 0] 74280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofUP1Quad = [2, 0, 0] 75280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofAP1Quad = [1, 0, 0] 76280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofUP2Quad = [2, 2, 2] 77280aadf6SBlaise Bourdin PetscInt,dimension(3),target :: dofAP2Quad = [1, 1, 1] 78280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofS3D = [0, 0, 0, 6] 79280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofUP1Tet = [3, 0, 0, 0] 80280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofAP1Tet = [1, 0, 0, 0] 81280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofUP2Tet = [3, 3, 0, 0] 82280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofAP2Tet = [1, 1, 0, 0] 83280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofUP1Hex = [3, 0, 0, 0] 84280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofAP1Hex = [1, 0, 0, 0] 85280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofUP2Hex = [3, 3, 3, 3] 86280aadf6SBlaise Bourdin PetscInt,dimension(4),target :: dofAP2Hex = [1, 1, 1, 1] 87280aadf6SBlaise Bourdin PetscInt,dimension(:),pointer :: dofU,dofA,dofS 88280aadf6SBlaise Bourdin 89280aadf6SBlaise Bourdin type(tPetscSF) :: migrationSF 90280aadf6SBlaise Bourdin PetscPartitioner :: part 91f51a5268SBarry Smith PetscLayout :: layout 92280aadf6SBlaise Bourdin 93280aadf6SBlaise Bourdin type(tVec) :: tmpVec 94280aadf6SBlaise Bourdin PetscReal :: norm 95280aadf6SBlaise Bourdin PetscReal :: time = 1.234_kPR 96280aadf6SBlaise Bourdin 97e2739ba6SAlexis Marboeuf PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER,ierr)) 98e2739ba6SAlexis Marboeuf if (ierr /= 0) then 99e2739ba6SAlexis Marboeuf print*,'Unable to initialize PETSc' 100e2739ba6SAlexis Marboeuf stop 101e2739ba6SAlexis Marboeuf endif 102280aadf6SBlaise Bourdin 103d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 104d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,numProc,ierr)) 105dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-i',ifilename,flg,ierr)) 106dcb3e689SBarry Smith PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing input file name -i <input file name>') 107dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-o',ofilename,flg,ierr)) 108dcb3e689SBarry Smith PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing output file name -o <output file name>') 109dcb3e689SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-order',order,flg,ierr)) 110280aadf6SBlaise Bourdin if ((order > 2) .or. (order < 1)) then 111280aadf6SBlaise Bourdin write(IOBuffer,'("Unsupported polynomial order ", I2, " not in [1,2]")') order 112de401611SBlaise Bourdin SETERRA(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 113280aadf6SBlaise Bourdin end if 114280aadf6SBlaise Bourdin 115280aadf6SBlaise Bourdin ! Read the mesh in any supported format 116d8606c27SBarry Smith PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr)) 117e2739ba6SAlexis Marboeuf PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr)); 118e2739ba6SAlexis Marboeuf PetscCallA(DMSetFromOptions(dm,ierr)); 119d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim,ierr)) 120dcb3e689SBarry Smith PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,'-dm_view',ierr)); 121280aadf6SBlaise Bourdin 122280aadf6SBlaise Bourdin ! Create the exodus result file 123280aadf6SBlaise Bourdin 124a5b23f4aSJose E. Roman ! enable exodus debugging information 125d8606c27SBarry Smith PetscCallA(exopts(EXVRBS+EXDEBG,ierr)) 126280aadf6SBlaise Bourdin ! Create the exodus file 127d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr)) 128280aadf6SBlaise Bourdin ! The long way would be 129280aadf6SBlaise Bourdin ! 130d8606c27SBarry Smith ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr)) 131d8606c27SBarry Smith ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr)) 132d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr)) 133d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr)) 134280aadf6SBlaise Bourdin 135280aadf6SBlaise Bourdin ! set the mesh order 136d8606c27SBarry Smith PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr)) 137d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr)) 138280aadf6SBlaise Bourdin ! 139280aadf6SBlaise Bourdin ! Notice how the exodus file is actually NOT open at this point (exoid is -1) 140d5b43468SJose E. Roman ! Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to 141280aadf6SBlaise Bourdin ! write the geometry (the DM), which can only be done on a brand new file. 142280aadf6SBlaise Bourdin ! 143280aadf6SBlaise Bourdin 144280aadf6SBlaise Bourdin ! Save the geometry to the file, erasing all previous content 145d8606c27SBarry Smith PetscCallA(DMView(dm,viewer,ierr)) 146d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr)) 147280aadf6SBlaise Bourdin ! 148280aadf6SBlaise Bourdin ! Note how the exodus file is now open 149280aadf6SBlaise Bourdin ! 150dcb3e689SBarry Smith ! 'Format' the exodus result file, i.e. allocate space for nodal and zonal variables 151280aadf6SBlaise Bourdin select case(sdim) 152280aadf6SBlaise Bourdin case(2) 153280aadf6SBlaise Bourdin numNodalVar = 4 154*49c89c76SBlaise Bourdin numZonalVar = 3 155*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariable(viewer, numZonalVar, ierr)) 156*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariable(viewer, numNodalVar, ierr)) 157*49c89c76SBlaise Bourdin do i = 1, numZonalVar 158*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariableName(viewer, i-1, zonalVarName2D(i), ierr)) 159*49c89c76SBlaise Bourdin end do 160*49c89c76SBlaise Bourdin do i = 1, numNodalVar 161*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariableName(viewer, i-1, nodalVarName2D(i), ierr)) 162*49c89c76SBlaise Bourdin end do 163*49c89c76SBlaise Bourdin case(3) 164*49c89c76SBlaise Bourdin numNodalVar = 5 165280aadf6SBlaise Bourdin numZonalVar = 6 166*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariable(viewer, numZonalVar, ierr)) 167*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariable(viewer, numNodalVar, ierr)) 168*49c89c76SBlaise Bourdin do i = 1, numZonalVar 169*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetZonalVariableName(viewer, i-1, zonalVarName3D(i), ierr)) 170*49c89c76SBlaise Bourdin end do 171*49c89c76SBlaise Bourdin do i = 1, numNodalVar 172*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIISetNodalVariableName(viewer, i-1, nodalVarName3D(i), ierr)) 173*49c89c76SBlaise Bourdin end do 174280aadf6SBlaise Bourdin case default 175280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 176280aadf6SBlaise Bourdin end select 177*49c89c76SBlaise Bourdin PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD,ierr)) 178280aadf6SBlaise Bourdin 179280aadf6SBlaise Bourdin ! An exodusII truth table specifies which fields are saved at which time step 180280aadf6SBlaise Bourdin ! It speeds up I/O but reserving space for fields in the file ahead of time. 181*49c89c76SBlaise Bourdin PetscCallA(PetscViewerExodusIIGetId(viewer,exoid,ierr)) 182*49c89c76SBlaise Bourdin PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK,numCS,PETSC_NULL_REAL,sjunk,ierr)) 183280aadf6SBlaise Bourdin allocate(truthtable(numCS,numZonalVar)) 184280aadf6SBlaise Bourdin truthtable = .true. 185d8606c27SBarry Smith PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr)) 186280aadf6SBlaise Bourdin deallocate(truthtable) 187280aadf6SBlaise Bourdin 188280aadf6SBlaise Bourdin ! Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 189280aadf6SBlaise Bourdin do step = 1,numstep 190d8606c27SBarry Smith PetscCallA(exptim(exoid,step,Real(step,kind=kPR),ierr)) 191280aadf6SBlaise Bourdin end do 192280aadf6SBlaise Bourdin 193d8606c27SBarry Smith PetscCallA(PetscObjectGetComm(dm,comm,ierr)) 194d8606c27SBarry Smith PetscCallA(PetscSectionCreate(comm, section,ierr)) 195d8606c27SBarry Smith PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr)) 196dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, 'U',ierr)) 197dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, 'Alpha',ierr)) 198dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldS, 'Sigma',ierr)) 199d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dm, pStart, pEnd,ierr)) 200d8606c27SBarry Smith PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr)) 201280aadf6SBlaise Bourdin 202280aadf6SBlaise Bourdin allocate(pStartDepth(sdim+1)) 203280aadf6SBlaise Bourdin allocate(pEndDepth(sdim+1)) 204280aadf6SBlaise Bourdin do d = 1, sdim+1 205d8606c27SBarry Smith PetscCallA(DMPlexGetDepthStratum(dm, d-1, pStartDepth(d), pEndDepth(d),ierr)) 206280aadf6SBlaise Bourdin end do 207280aadf6SBlaise Bourdin 208280aadf6SBlaise Bourdin ! Vector field U, Scalar field Alpha, Tensor field Sigma 209e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr)); 210e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr)); 211e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr)); 212280aadf6SBlaise Bourdin 213280aadf6SBlaise Bourdin ! Going through cell sets then cells, and setting up storage for the sections 214dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dm, 'Cell Sets', numCS, ierr)) 215dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dm, 'Cell Sets', csIS, ierr)) 216d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID, ierr)) 217280aadf6SBlaise Bourdin do set = 1,numCS 2188f2079feSBlaise Bourdin !!! This is pointless but will avoid a warning with gfortran and -Werror=maybe-uninitialized 2198f2079feSBlaise Bourdin nullify(dofA) 2208f2079feSBlaise Bourdin nullify(dofU) 2218f2079feSBlaise Bourdin nullify(dofS) 222dcb3e689SBarry Smith PetscCallA(DMGetStratumSize(dm, 'Cell Sets', csID(set), numCells,ierr)) 223dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dm, 'Cell Sets', csID(set), cellIS,ierr)) 224280aadf6SBlaise Bourdin if (numCells > 0) then 225280aadf6SBlaise Bourdin select case(sdim) 226280aadf6SBlaise Bourdin case(2) 227280aadf6SBlaise Bourdin dofs => dofS2D 228280aadf6SBlaise Bourdin case(3) 229280aadf6SBlaise Bourdin dofs => dofS3D 230280aadf6SBlaise Bourdin case default 231280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 232de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 233280aadf6SBlaise Bourdin end select ! sdim 234280aadf6SBlaise Bourdin 235280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 236280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 237d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 238280aadf6SBlaise Bourdin nullify(closureA) 239d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr)) 240280aadf6SBlaise Bourdin select case(size(closureA)/2) 241280aadf6SBlaise Bourdin case(7) ! Tri 242280aadf6SBlaise Bourdin if (order == 1) then 243280aadf6SBlaise Bourdin dofU => dofUP1Tri 244280aadf6SBlaise Bourdin dofA => dofAP1Tri 245280aadf6SBlaise Bourdin else 246280aadf6SBlaise Bourdin dofU => dofUP2Tri 247280aadf6SBlaise Bourdin dofA => dofAP2Tri 248280aadf6SBlaise Bourdin end if 249280aadf6SBlaise Bourdin case(9) ! Quad 250280aadf6SBlaise Bourdin if (order == 1) then 251280aadf6SBlaise Bourdin dofU => dofUP1Quad 252280aadf6SBlaise Bourdin dofA => dofAP1Quad 253280aadf6SBlaise Bourdin else 254280aadf6SBlaise Bourdin dofU => dofUP2Quad 255280aadf6SBlaise Bourdin dofA => dofAP2Quad 256280aadf6SBlaise Bourdin end if 257280aadf6SBlaise Bourdin case(15) ! Tet 258280aadf6SBlaise Bourdin if (order == 1) then 259280aadf6SBlaise Bourdin dofU => dofUP1Tet 260280aadf6SBlaise Bourdin dofA => dofAP1Tet 261280aadf6SBlaise Bourdin else 262280aadf6SBlaise Bourdin dofU => dofUP2Tet 263280aadf6SBlaise Bourdin dofA => dofAP2Tet 264280aadf6SBlaise Bourdin end if 265280aadf6SBlaise Bourdin case(27) ! Hex 266280aadf6SBlaise Bourdin if (order == 1) then 267280aadf6SBlaise Bourdin dofU => dofUP1Hex 268280aadf6SBlaise Bourdin dofA => dofAP1Hex 269280aadf6SBlaise Bourdin else 270280aadf6SBlaise Bourdin dofU => dofUP2Hex 271280aadf6SBlaise Bourdin dofA => dofAP2Hex 272280aadf6SBlaise Bourdin end if 273280aadf6SBlaise Bourdin case default 274280aadf6SBlaise Bourdin write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2 275de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer) 276280aadf6SBlaise Bourdin end select 277d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr)) 278280aadf6SBlaise Bourdin do cell = 1,numCells! 279280aadf6SBlaise Bourdin nullify(closure) 280d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 281280aadf6SBlaise Bourdin do p = 1,size(closure),2 282280aadf6SBlaise Bourdin ! find the depth of p 283280aadf6SBlaise Bourdin do d = 1,sdim+1 284280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 285d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr)) 286d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr)) 287d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr)) 288d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr)) 289280aadf6SBlaise Bourdin end if ! closure(p) 290280aadf6SBlaise Bourdin end do ! d 291280aadf6SBlaise Bourdin end do ! p 292d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 293280aadf6SBlaise Bourdin end do ! cell 294d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 295d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 296280aadf6SBlaise Bourdin end if ! numCells 297280aadf6SBlaise Bourdin end do ! set 298d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 299d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 300d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section,ierr)) 301d8606c27SBarry Smith PetscCallA(DMSetLocalSection(dm, section,ierr)) 302dcb3e689SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, '-dm_section_view',ierr)) 303f51a5268SBarry Smith PetscCallA(PetscSectionGetValueLayout(PETSC_COMM_WORLD, section, layout, ierr)) 304f51a5268SBarry Smith PetscCallA(PetscLayoutDestroy(layout,ierr)) 305d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section,ierr)) 306280aadf6SBlaise Bourdin 307d8606c27SBarry Smith PetscCallA(DMSetUseNatural(dm,PETSC_TRUE,ierr)) 308d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm,part,ierr)) 309d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part,ierr)) 310d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr)) 311280aadf6SBlaise Bourdin 312280aadf6SBlaise Bourdin if (numProc > 1) then 313d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr)) 314d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF,ierr)) 315e2739ba6SAlexis Marboeuf else 316e2739ba6SAlexis Marboeuf pdm = dm 317280aadf6SBlaise Bourdin end if 318dcb3e689SBarry Smith PetscCallA(DMViewFromOptions(pdm,PETSC_NULL_OPTIONS,'-dm_view',ierr)) 319280aadf6SBlaise Bourdin 320280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 321e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldU, isU, dmU,ierr)) 322e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldA, isA, dmA,ierr)) 323e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldS, isS, dmS,ierr)) 324e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA,ierr)) 325280aadf6SBlaise Bourdin 326280aadf6SBlaise Bourdin !Create the exodus result file 327280aadf6SBlaise Bourdin allocate(dmList(2)) 328280aadf6SBlaise Bourdin dmList(1) = dmU; 329280aadf6SBlaise Bourdin dmList(2) = dmA; 330d8606c27SBarry Smith PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr)) 331280aadf6SBlaise Bourdin deallocate(dmList) 332280aadf6SBlaise Bourdin 333e2739ba6SAlexis Marboeuf PetscCallA(DMGetGlobalVector(pdm, X,ierr)) 334d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U,ierr)) 335d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A,ierr)) 336d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S,ierr)) 337d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA,ierr)) 338d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr)) 339280aadf6SBlaise Bourdin 340dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(U, 'U',ierr)) 341dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(A, 'Alpha',ierr)) 342dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(S, 'Sigma',ierr)) 343dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA, 'UAlpha',ierr)) 344dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA2, 'UAlpha2',ierr)) 345d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR,ierr)) 346280aadf6SBlaise Bourdin 347280aadf6SBlaise 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 */ 348d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr)) 349d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr)) 350d8606c27SBarry Smith PetscCallA(VecGetArrayF90(UALoc, cval,ierr)) 351d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr)) 352d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr)) 353d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr)) 354280aadf6SBlaise Bourdin 355280aadf6SBlaise Bourdin do p = pStart,pEnd-1 356d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr)) 357280aadf6SBlaise Bourdin if (dofUA > 0) then 358d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr)) 359d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr)) 360280aadf6SBlaise Bourdin closureSize = size(xyz) 361280aadf6SBlaise Bourdin do i = 1,sdim 362280aadf6SBlaise Bourdin do j = 0, closureSize-1,sdim 363280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i) 364280aadf6SBlaise Bourdin end do 365280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) * sdim / closureSize; 366280aadf6SBlaise Bourdin cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2 367280aadf6SBlaise Bourdin end do 368d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr)) 369280aadf6SBlaise Bourdin end if 370280aadf6SBlaise Bourdin end do 371280aadf6SBlaise Bourdin 372d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr)) 373d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 374d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 375d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr)) 376280aadf6SBlaise Bourdin 377280aadf6SBlaise Bourdin !Update X 378d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr)) 379280aadf6SBlaise Bourdin ! Restrict to U and Alpha 380d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr)) 381d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr)) 382dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, '-ua_vec_view',ierr)) 383dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, '-u_vec_view',ierr)) 384dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, '-a_vec_view',ierr)) 385280aadf6SBlaise Bourdin ! restrict to UA2 386d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr)) 387dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, '-ua2_vec_view',ierr)) 388280aadf6SBlaise Bourdin 389e2739ba6SAlexis Marboeuf ! Getting Natural Vec 390d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr)) 391d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr)) 392280aadf6SBlaise Bourdin 393d8606c27SBarry Smith PetscCallA(VecView(U, viewer,ierr)) 394d8606c27SBarry Smith PetscCallA(VecView(A, viewer,ierr)) 395280aadf6SBlaise Bourdin 396280aadf6SBlaise Bourdin !Saving U and Alpha in one shot. 397280aadf6SBlaise Bourdin !For this, we need to cheat and change the Vec's name 398280aadf6SBlaise Bourdin !Note that in the end we write variables one component at a time, 399280aadf6SBlaise Bourdin !so that there is no real value in doing this 400d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr)) 401d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr)) 402d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec,ierr)) 403dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr)) 404d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 405280aadf6SBlaise Bourdin 406280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 407d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 408d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer,ierr)) 409d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr)) 410d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr)) 411280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 412280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 413280aadf6SBlaise Bourdin end if 414d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr)) 415280aadf6SBlaise Bourdin 416280aadf6SBlaise Bourdin ! same thing with the UA2 Vec obtained from the superDM 417d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr)) 418d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec,ierr)) 419dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr)) 420d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr)) 421d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 422280aadf6SBlaise Bourdin 423280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 424d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 425d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 426d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr)) 427d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr)) 428280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 429280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 430280aadf6SBlaise Bourdin end if 431d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr)) 432280aadf6SBlaise Bourdin 433280aadf6SBlaise Bourdin ! Building and saving Sigma 434280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 435280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 436280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 437d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr)) 438d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr)) 439dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, 'Cell Sets', csIS,ierr)) 440dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dmS, 'Cell Sets',numCS,ierr)) 441d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID,ierr)) 442280aadf6SBlaise Bourdin 443280aadf6SBlaise Bourdin do set = 1, numCS 444dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dmS, 'Cell Sets', csID(set), cellIS,ierr)) 445d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 446d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells,ierr)) 447280aadf6SBlaise Bourdin do cell = 1,numCells 448d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 449d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 450280aadf6SBlaise Bourdin cval(1) = rank 451280aadf6SBlaise Bourdin cval(2) = csID(set) 452280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 453280aadf6SBlaise Bourdin do c = 1, size(xyz),sdim 454280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 455280aadf6SBlaise Bourdin end do 456280aadf6SBlaise Bourdin cval(3) = cval(3) * sdim / size(xyz) 457d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr)) 458d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 459d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 460280aadf6SBlaise Bourdin end do 461d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 462d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 463280aadf6SBlaise Bourdin end do 464d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 465d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 466dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, '-s_vec_view',ierr)) 467280aadf6SBlaise Bourdin 468280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 469d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr)) 470d8606c27SBarry Smith PetscCallA(VecView(S,viewer,ierr)) 471280aadf6SBlaise Bourdin 472280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 473d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr)) 474d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 475dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'Sigma',ierr)) 476d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 477d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr)) 478d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr)) 479280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 480280aadf6SBlaise Bourdin write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm 481280aadf6SBlaise Bourdin end if 482d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr)) 483280aadf6SBlaise Bourdin 484d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr)) 485d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr)) 486d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S,ierr)) 487d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A,ierr)) 488d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U,ierr)) 489e2739ba6SAlexis Marboeuf PetscCallA(DMRestoreGlobalVector(pdm, X,ierr)) 490d8606c27SBarry Smith PetscCallA(DMDestroy(dmU,ierr)) 491d8606c27SBarry Smith PetscCallA(ISDestroy(isU,ierr)) 492d8606c27SBarry Smith PetscCallA(DMDestroy(dmA,ierr)) 493d8606c27SBarry Smith PetscCallA(ISDestroy(isA,ierr)) 494d8606c27SBarry Smith PetscCallA(DMDestroy(dmS,ierr)) 495d8606c27SBarry Smith PetscCallA(ISDestroy(isS,ierr)) 496d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA,ierr)) 497d8606c27SBarry Smith PetscCallA(ISDestroy(isUA,ierr)) 498d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2,ierr)) 499e2739ba6SAlexis Marboeuf PetscCallA(DMDestroy(pdm,ierr)) 500e2739ba6SAlexis Marboeuf if (numProc > 1) then 501d8606c27SBarry Smith PetscCallA(DMDestroy(dm,ierr)) 502e2739ba6SAlexis Marboeuf end if 503280aadf6SBlaise Bourdin 504280aadf6SBlaise Bourdin deallocate(pStartDepth) 505280aadf6SBlaise Bourdin deallocate(pEndDepth) 506280aadf6SBlaise Bourdin 507d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer,ierr)) 508d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 5092e8d78feSBlaise Bourdinend program ex62f90 510280aadf6SBlaise Bourdin 511280aadf6SBlaise Bourdin! /*TEST 512280aadf6SBlaise Bourdin! 513280aadf6SBlaise Bourdin! build: 514280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 515280aadf6SBlaise Bourdin! # 2D seq 516e2739ba6SAlexis Marboeuf! test: 517e2739ba6SAlexis Marboeuf! suffix: 0 518e2739ba6SAlexis 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 519e2739ba6SAlexis 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 520e2739ba6SAlexis Marboeuf! test: 521e2739ba6SAlexis Marboeuf! suffix: 1 522e2739ba6SAlexis 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 523e2739ba6SAlexis Marboeuf! 524e2739ba6SAlexis Marboeuf! test: 525e2739ba6SAlexis Marboeuf! suffix: 2 526e2739ba6SAlexis 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 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: 3 530e2739ba6SAlexis 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 531e2739ba6SAlexis Marboeuf! test: 532e2739ba6SAlexis Marboeuf! suffix: 4 533e2739ba6SAlexis 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 534e2739ba6SAlexis Marboeuf! test: 535e2739ba6SAlexis Marboeuf! suffix: 5 536e2739ba6SAlexis 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 537e2739ba6SAlexis Marboeuf! # 2D par 538e2739ba6SAlexis Marboeuf! test: 539e2739ba6SAlexis Marboeuf! suffix: 6 540e2739ba6SAlexis Marboeuf! nsize: 2 541e2739ba6SAlexis 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 542e2739ba6SAlexis 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 543e2739ba6SAlexis Marboeuf! test: 544e2739ba6SAlexis Marboeuf! suffix: 7 545e2739ba6SAlexis Marboeuf! nsize: 2 546e2739ba6SAlexis 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 547e2739ba6SAlexis Marboeuf! test: 548e2739ba6SAlexis Marboeuf! suffix: 8 549e2739ba6SAlexis Marboeuf! nsize: 2 550e2739ba6SAlexis 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 551e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 552e2739ba6SAlexis Marboeuf! test: 553e2739ba6SAlexis Marboeuf! suffix: 9 554e2739ba6SAlexis Marboeuf! nsize: 2 555e2739ba6SAlexis 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 556e2739ba6SAlexis Marboeuf! test: 557e2739ba6SAlexis Marboeuf! suffix: 10 558e2739ba6SAlexis Marboeuf! nsize: 2 559e2739ba6SAlexis 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 560e2739ba6SAlexis Marboeuf! test: 5619c89aa79SPierre Jolivet! # Something is now broken with parallel read/write for whatever shape H is 562e2739ba6SAlexis Marboeuf! TODO: broken 563e2739ba6SAlexis Marboeuf! suffix: 11 564e2739ba6SAlexis Marboeuf! nsize: 2 565e2739ba6SAlexis 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 566e2739ba6SAlexis Marboeuf 567e2739ba6SAlexis Marboeuf! #3d seq 568e2739ba6SAlexis Marboeuf! test: 569e2739ba6SAlexis Marboeuf! suffix: 12 570e2739ba6SAlexis 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 571e2739ba6SAlexis Marboeuf! test: 572e2739ba6SAlexis Marboeuf! suffix: 13 573e2739ba6SAlexis 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 574e2739ba6SAlexis 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 575e2739ba6SAlexis Marboeuf! test: 576e2739ba6SAlexis Marboeuf! suffix: 14 577e2739ba6SAlexis 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 578e2739ba6SAlexis Marboeuf! test: 579e2739ba6SAlexis Marboeuf! suffix: 15 580e2739ba6SAlexis 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 581e2739ba6SAlexis 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 582e2739ba6SAlexis Marboeuf! #3d par 583e2739ba6SAlexis Marboeuf! test: 584e2739ba6SAlexis Marboeuf! suffix: 16 585e2739ba6SAlexis Marboeuf! nsize: 2 586e2739ba6SAlexis 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 587e2739ba6SAlexis Marboeuf! test: 588e2739ba6SAlexis Marboeuf! suffix: 17 589e2739ba6SAlexis Marboeuf! nsize: 2 590e2739ba6SAlexis 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 591e2739ba6SAlexis 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 592e2739ba6SAlexis Marboeuf! test: 593e2739ba6SAlexis Marboeuf! suffix: 18 594e2739ba6SAlexis Marboeuf! nsize: 2 595e2739ba6SAlexis 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 596e2739ba6SAlexis Marboeuf! test: 597e2739ba6SAlexis Marboeuf! suffix: 19 598e2739ba6SAlexis Marboeuf! nsize: 2 599e2739ba6SAlexis 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 600e2739ba6SAlexis 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 601e2739ba6SAlexis Marboeuf 602280aadf6SBlaise Bourdin! TEST*/ 603