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 7280aadf6SBlaise Bourdin ! 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 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 73280aadf6SBlaise Bourdin type(tPetscSF) :: migrationSF 74280aadf6SBlaise Bourdin PetscPartitioner :: part 75280aadf6SBlaise Bourdin 76280aadf6SBlaise Bourdin type(tVec) :: tmpVec 77280aadf6SBlaise Bourdin PetscReal :: norm 78280aadf6SBlaise Bourdin PetscReal :: time = 1.234_kPR 79280aadf6SBlaise Bourdin 80e2739ba6SAlexis Marboeuf PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER,ierr)) 81e2739ba6SAlexis Marboeuf if (ierr /= 0) then 82e2739ba6SAlexis Marboeuf print*,'Unable to initialize PETSc' 83e2739ba6SAlexis Marboeuf stop 84e2739ba6SAlexis Marboeuf endif 85280aadf6SBlaise Bourdin 86d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 87d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,numProc,ierr)) 88*dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-i',ifilename,flg,ierr)) 89*dcb3e689SBarry Smith PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing input file name -i <input file name>') 90*dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-o',ofilename,flg,ierr)) 91*dcb3e689SBarry Smith PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing output file name -o <output file name>') 92*dcb3e689SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-order',order,flg,ierr)) 93280aadf6SBlaise Bourdin if ((order > 2) .or. (order < 1)) then 94280aadf6SBlaise Bourdin write(IOBuffer,'("Unsupported polynomial order ", I2, " not in [1,2]")') order 95de401611SBlaise Bourdin SETERRA(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 96280aadf6SBlaise Bourdin end if 97280aadf6SBlaise Bourdin 98280aadf6SBlaise Bourdin ! Read the mesh in any supported format 99d8606c27SBarry Smith PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr)) 100e2739ba6SAlexis Marboeuf PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr)); 101e2739ba6SAlexis Marboeuf PetscCallA(DMSetFromOptions(dm,ierr)); 102d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim,ierr)) 103*dcb3e689SBarry Smith PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,'-dm_view',ierr)); 104280aadf6SBlaise Bourdin 105280aadf6SBlaise Bourdin ! Create the exodus result file 106280aadf6SBlaise Bourdin 107a5b23f4aSJose E. Roman ! enable exodus debugging information 108d8606c27SBarry Smith PetscCallA(exopts(EXVRBS+EXDEBG,ierr)) 109280aadf6SBlaise Bourdin ! Create the exodus file 110d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr)) 111280aadf6SBlaise Bourdin ! The long way would be 112280aadf6SBlaise Bourdin ! 113d8606c27SBarry Smith ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr)) 114d8606c27SBarry Smith ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr)) 115d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr)) 116d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr)) 117280aadf6SBlaise Bourdin 118280aadf6SBlaise Bourdin ! set the mesh order 119d8606c27SBarry Smith PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr)) 120d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr)) 121280aadf6SBlaise Bourdin ! 122280aadf6SBlaise Bourdin ! Notice how the exodus file is actually NOT open at this point (exoid is -1) 123d5b43468SJose E. Roman ! Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to 124280aadf6SBlaise Bourdin ! write the geometry (the DM), which can only be done on a brand new file. 125280aadf6SBlaise Bourdin ! 126280aadf6SBlaise Bourdin 127280aadf6SBlaise Bourdin ! Save the geometry to the file, erasing all previous content 128d8606c27SBarry Smith PetscCallA(DMView(dm,viewer,ierr)) 129d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr)) 130280aadf6SBlaise Bourdin ! 131280aadf6SBlaise Bourdin ! Note how the exodus file is now open 132280aadf6SBlaise Bourdin ! 133*dcb3e689SBarry Smith ! 'Format' the exodus result file, i.e. allocate space for nodal and zonal variables 134280aadf6SBlaise Bourdin select case(sdim) 135280aadf6SBlaise Bourdin case(2) 136280aadf6SBlaise Bourdin numNodalVar = 3 137*dcb3e689SBarry Smith nodalVarName(1:numNodalVar) = ['U_x ','U_y ','Alpha'] 138280aadf6SBlaise Bourdin numZonalVar = 3 139*dcb3e689SBarry Smith zonalVarName(1:numZonalVar) = ['Sigma_11','Sigma_22','Sigma_12'] 140280aadf6SBlaise Bourdin case(3) 141280aadf6SBlaise Bourdin numNodalVar = 4 142*dcb3e689SBarry Smith nodalVarName(1:numNodalVar) = ['U_x ','U_y ','U_z ','Alpha'] 143280aadf6SBlaise Bourdin numZonalVar = 6 144*dcb3e689SBarry Smith zonalVarName(1:numZonalVar) = ['Sigma_11','Sigma_22','Sigma_33','Sigma_23','Sigma_13','Sigma_12'] 145280aadf6SBlaise Bourdin case default 146280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 147280aadf6SBlaise Bourdin end select 148d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIGetId(viewer,exoid,ierr)) 149*dcb3e689SBarry Smith PetscCallA(expvp(exoid, 'E', numZonalVar,ierr)) 150*dcb3e689SBarry Smith PetscCallA(expvan(exoid, 'E', numZonalVar, zonalVarName,ierr)) 151*dcb3e689SBarry Smith PetscCallA(expvp(exoid, 'N', numNodalVar,ierr)) 152*dcb3e689SBarry Smith PetscCallA(expvan(exoid, 'N', numNodalVar, nodalVarName,ierr)) 153d8606c27SBarry Smith PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK,numCS,PETSC_NULL_REAL,sjunk,ierr)) 154280aadf6SBlaise Bourdin 155280aadf6SBlaise Bourdin ! An exodusII truth table specifies which fields are saved at which time step 156280aadf6SBlaise Bourdin ! It speeds up I/O but reserving space for fields in the file ahead of time. 157280aadf6SBlaise Bourdin allocate(truthtable(numCS,numZonalVar)) 158280aadf6SBlaise Bourdin truthtable = .true. 159d8606c27SBarry Smith PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr)) 160280aadf6SBlaise Bourdin deallocate(truthtable) 161280aadf6SBlaise Bourdin 162280aadf6SBlaise Bourdin ! Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 163280aadf6SBlaise Bourdin do step = 1,numstep 164d8606c27SBarry Smith PetscCallA(exptim(exoid,step,Real(step,kind=kPR),ierr)) 165280aadf6SBlaise Bourdin end do 166280aadf6SBlaise Bourdin 167d8606c27SBarry Smith PetscCallA(PetscObjectGetComm(dm,comm,ierr)) 168d8606c27SBarry Smith PetscCallA(PetscSectionCreate(comm, section,ierr)) 169d8606c27SBarry Smith PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr)) 170*dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, 'U',ierr)) 171*dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, 'Alpha',ierr)) 172*dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldS, 'Sigma',ierr)) 173d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dm, pStart, pEnd,ierr)) 174d8606c27SBarry Smith PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr)) 175280aadf6SBlaise Bourdin 176280aadf6SBlaise Bourdin allocate(pStartDepth(sdim+1)) 177280aadf6SBlaise Bourdin allocate(pEndDepth(sdim+1)) 178280aadf6SBlaise Bourdin do d = 1, sdim+1 179d8606c27SBarry Smith PetscCallA(DMPlexGetDepthStratum(dm, d-1, pStartDepth(d), pEndDepth(d),ierr)) 180280aadf6SBlaise Bourdin end do 181280aadf6SBlaise Bourdin 182280aadf6SBlaise Bourdin ! Vector field U, Scalar field Alpha, Tensor field Sigma 183e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr)); 184e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr)); 185e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr)); 186280aadf6SBlaise Bourdin 187280aadf6SBlaise Bourdin ! Going through cell sets then cells, and setting up storage for the sections 188*dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dm, 'Cell Sets', numCS, ierr)) 189*dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dm, 'Cell Sets', csIS, ierr)) 190d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID, ierr)) 191280aadf6SBlaise Bourdin do set = 1,numCS 1928f2079feSBlaise Bourdin !!! This is pointless but will avoid a warning with gfortran and -Werror=maybe-uninitialized 1938f2079feSBlaise Bourdin nullify(dofA) 1948f2079feSBlaise Bourdin nullify(dofU) 1958f2079feSBlaise Bourdin nullify(dofS) 196*dcb3e689SBarry Smith PetscCallA(DMGetStratumSize(dm, 'Cell Sets', csID(set), numCells,ierr)) 197*dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dm, 'Cell Sets', csID(set), cellIS,ierr)) 198280aadf6SBlaise Bourdin if (numCells > 0) then 199280aadf6SBlaise Bourdin select case(sdim) 200280aadf6SBlaise Bourdin case(2) 201280aadf6SBlaise Bourdin dofs => dofS2D 202280aadf6SBlaise Bourdin case(3) 203280aadf6SBlaise Bourdin dofs => dofS3D 204280aadf6SBlaise Bourdin case default 205280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 206de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 207280aadf6SBlaise Bourdin end select ! sdim 208280aadf6SBlaise Bourdin 209280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 210280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 211d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 212280aadf6SBlaise Bourdin nullify(closureA) 213d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr)) 214280aadf6SBlaise Bourdin select case(size(closureA)/2) 215280aadf6SBlaise Bourdin case(7) ! Tri 216280aadf6SBlaise Bourdin if (order == 1) then 217280aadf6SBlaise Bourdin dofU => dofUP1Tri 218280aadf6SBlaise Bourdin dofA => dofAP1Tri 219280aadf6SBlaise Bourdin else 220280aadf6SBlaise Bourdin dofU => dofUP2Tri 221280aadf6SBlaise Bourdin dofA => dofAP2Tri 222280aadf6SBlaise Bourdin end if 223280aadf6SBlaise Bourdin case(9) ! Quad 224280aadf6SBlaise Bourdin if (order == 1) then 225280aadf6SBlaise Bourdin dofU => dofUP1Quad 226280aadf6SBlaise Bourdin dofA => dofAP1Quad 227280aadf6SBlaise Bourdin else 228280aadf6SBlaise Bourdin dofU => dofUP2Quad 229280aadf6SBlaise Bourdin dofA => dofAP2Quad 230280aadf6SBlaise Bourdin end if 231280aadf6SBlaise Bourdin case(15) ! Tet 232280aadf6SBlaise Bourdin if (order == 1) then 233280aadf6SBlaise Bourdin dofU => dofUP1Tet 234280aadf6SBlaise Bourdin dofA => dofAP1Tet 235280aadf6SBlaise Bourdin else 236280aadf6SBlaise Bourdin dofU => dofUP2Tet 237280aadf6SBlaise Bourdin dofA => dofAP2Tet 238280aadf6SBlaise Bourdin end if 239280aadf6SBlaise Bourdin case(27) ! Hex 240280aadf6SBlaise Bourdin if (order == 1) then 241280aadf6SBlaise Bourdin dofU => dofUP1Hex 242280aadf6SBlaise Bourdin dofA => dofAP1Hex 243280aadf6SBlaise Bourdin else 244280aadf6SBlaise Bourdin dofU => dofUP2Hex 245280aadf6SBlaise Bourdin dofA => dofAP2Hex 246280aadf6SBlaise Bourdin end if 247280aadf6SBlaise Bourdin case default 248280aadf6SBlaise Bourdin write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2 249de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer) 250280aadf6SBlaise Bourdin end select 251d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr)) 252280aadf6SBlaise Bourdin do cell = 1,numCells! 253280aadf6SBlaise Bourdin nullify(closure) 254d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 255280aadf6SBlaise Bourdin do p = 1,size(closure),2 256280aadf6SBlaise Bourdin ! find the depth of p 257280aadf6SBlaise Bourdin do d = 1,sdim+1 258280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 259d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr)) 260d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr)) 261d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr)) 262d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr)) 263280aadf6SBlaise Bourdin end if ! closure(p) 264280aadf6SBlaise Bourdin end do ! d 265280aadf6SBlaise Bourdin end do ! p 266d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 267280aadf6SBlaise Bourdin end do ! cell 268d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 269d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 270280aadf6SBlaise Bourdin end if ! numCells 271280aadf6SBlaise Bourdin end do ! set 272d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 273d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 274d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section,ierr)) 275d8606c27SBarry Smith PetscCallA(DMSetLocalSection(dm, section,ierr)) 276*dcb3e689SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, '-dm_section_view',ierr)) 277d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section,ierr)) 278280aadf6SBlaise Bourdin 279d8606c27SBarry Smith PetscCallA(DMSetUseNatural(dm,PETSC_TRUE,ierr)) 280d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm,part,ierr)) 281d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part,ierr)) 282d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr)) 283280aadf6SBlaise Bourdin 284280aadf6SBlaise Bourdin if (numProc > 1) then 285d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr)) 286d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF,ierr)) 287e2739ba6SAlexis Marboeuf else 288e2739ba6SAlexis Marboeuf pdm = dm 289280aadf6SBlaise Bourdin end if 290*dcb3e689SBarry Smith PetscCallA(DMViewFromOptions(pdm,PETSC_NULL_OPTIONS,'-dm_view',ierr)) 291280aadf6SBlaise Bourdin 292280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 293e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldU, isU, dmU,ierr)) 294e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldA, isA, dmA,ierr)) 295e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldS, isS, dmS,ierr)) 296e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA,ierr)) 297280aadf6SBlaise Bourdin 298280aadf6SBlaise Bourdin !Create the exodus result file 299280aadf6SBlaise Bourdin allocate(dmList(2)) 300280aadf6SBlaise Bourdin dmList(1) = dmU; 301280aadf6SBlaise Bourdin dmList(2) = dmA; 302d8606c27SBarry Smith PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr)) 303280aadf6SBlaise Bourdin deallocate(dmList) 304280aadf6SBlaise Bourdin 305e2739ba6SAlexis Marboeuf PetscCallA(DMGetGlobalVector(pdm, X,ierr)) 306d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U,ierr)) 307d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A,ierr)) 308d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S,ierr)) 309d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA,ierr)) 310d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr)) 311280aadf6SBlaise Bourdin 312*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(U, 'U',ierr)) 313*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(A, 'Alpha',ierr)) 314*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(S, 'Sigma',ierr)) 315*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA, 'UAlpha',ierr)) 316*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(UA2, 'UAlpha2',ierr)) 317d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR,ierr)) 318280aadf6SBlaise Bourdin 319280aadf6SBlaise 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 */ 320d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr)) 321d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr)) 322d8606c27SBarry Smith PetscCallA(VecGetArrayF90(UALoc, cval,ierr)) 323d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr)) 324d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr)) 325d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr)) 326280aadf6SBlaise Bourdin 327280aadf6SBlaise Bourdin do p = pStart,pEnd-1 328d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr)) 329280aadf6SBlaise Bourdin if (dofUA > 0) then 330d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr)) 331d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr)) 332280aadf6SBlaise Bourdin closureSize = size(xyz) 333280aadf6SBlaise Bourdin do i = 1,sdim 334280aadf6SBlaise Bourdin do j = 0, closureSize-1,sdim 335280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i) 336280aadf6SBlaise Bourdin end do 337280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) * sdim / closureSize; 338280aadf6SBlaise Bourdin cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2 339280aadf6SBlaise Bourdin end do 340d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr)) 341280aadf6SBlaise Bourdin end if 342280aadf6SBlaise Bourdin end do 343280aadf6SBlaise Bourdin 344d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr)) 345d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 346d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 347d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr)) 348280aadf6SBlaise Bourdin 349280aadf6SBlaise Bourdin !Update X 350d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr)) 351280aadf6SBlaise Bourdin ! Restrict to U and Alpha 352d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr)) 353d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr)) 354*dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, '-ua_vec_view',ierr)) 355*dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, '-u_vec_view',ierr)) 356*dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, '-a_vec_view',ierr)) 357280aadf6SBlaise Bourdin ! restrict to UA2 358d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr)) 359*dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, '-ua2_vec_view',ierr)) 360280aadf6SBlaise Bourdin 361e2739ba6SAlexis Marboeuf ! Getting Natural Vec 362d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr)) 363d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr)) 364280aadf6SBlaise Bourdin 365d8606c27SBarry Smith PetscCallA(VecView(U, viewer,ierr)) 366d8606c27SBarry Smith PetscCallA(VecView(A, viewer,ierr)) 367280aadf6SBlaise Bourdin 368280aadf6SBlaise Bourdin !Saving U and Alpha in one shot. 369280aadf6SBlaise Bourdin !For this, we need to cheat and change the Vec's name 370280aadf6SBlaise Bourdin !Note that in the end we write variables one component at a time, 371280aadf6SBlaise Bourdin !so that there is no real value in doing this 372d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr)) 373d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr)) 374d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec,ierr)) 375*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr)) 376d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 377280aadf6SBlaise Bourdin 378280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 379d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 380d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer,ierr)) 381d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr)) 382d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr)) 383280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 384280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 385280aadf6SBlaise Bourdin end if 386d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr)) 387280aadf6SBlaise Bourdin 388280aadf6SBlaise Bourdin ! same thing with the UA2 Vec obtained from the superDM 389d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr)) 390d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec,ierr)) 391*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr)) 392d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr)) 393d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 394280aadf6SBlaise Bourdin 395280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 396d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 397d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 398d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr)) 399d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr)) 400280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 401280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 402280aadf6SBlaise Bourdin end if 403d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr)) 404280aadf6SBlaise Bourdin 405280aadf6SBlaise Bourdin ! Building and saving Sigma 406280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 407280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 408280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 409d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr)) 410d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr)) 411*dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, 'Cell Sets', csIS,ierr)) 412*dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dmS, 'Cell Sets',numCS,ierr)) 413d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID,ierr)) 414280aadf6SBlaise Bourdin 415280aadf6SBlaise Bourdin do set = 1, numCS 416*dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dmS, 'Cell Sets', csID(set), cellIS,ierr)) 417d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 418d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells,ierr)) 419280aadf6SBlaise Bourdin do cell = 1,numCells 420d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 421d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 422280aadf6SBlaise Bourdin cval(1) = rank 423280aadf6SBlaise Bourdin cval(2) = csID(set) 424280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 425280aadf6SBlaise Bourdin do c = 1, size(xyz),sdim 426280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 427280aadf6SBlaise Bourdin end do 428280aadf6SBlaise Bourdin cval(3) = cval(3) * sdim / size(xyz) 429d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr)) 430d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 431d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 432280aadf6SBlaise Bourdin end do 433d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 434d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 435280aadf6SBlaise Bourdin end do 436d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 437d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 438*dcb3e689SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, '-s_vec_view',ierr)) 439280aadf6SBlaise Bourdin 440280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 441d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr)) 442d8606c27SBarry Smith PetscCallA(VecView(S,viewer,ierr)) 443280aadf6SBlaise Bourdin 444280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 445d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr)) 446d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 447*dcb3e689SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, 'Sigma',ierr)) 448d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 449d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr)) 450d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr)) 451280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 452280aadf6SBlaise Bourdin write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm 453280aadf6SBlaise Bourdin end if 454d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr)) 455280aadf6SBlaise Bourdin 456d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr)) 457d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr)) 458d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S,ierr)) 459d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A,ierr)) 460d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U,ierr)) 461e2739ba6SAlexis Marboeuf PetscCallA(DMRestoreGlobalVector(pdm, X,ierr)) 462d8606c27SBarry Smith PetscCallA(DMDestroy(dmU,ierr)) 463d8606c27SBarry Smith PetscCallA(ISDestroy(isU,ierr)) 464d8606c27SBarry Smith PetscCallA(DMDestroy(dmA,ierr)) 465d8606c27SBarry Smith PetscCallA(ISDestroy(isA,ierr)) 466d8606c27SBarry Smith PetscCallA(DMDestroy(dmS,ierr)) 467d8606c27SBarry Smith PetscCallA(ISDestroy(isS,ierr)) 468d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA,ierr)) 469d8606c27SBarry Smith PetscCallA(ISDestroy(isUA,ierr)) 470d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2,ierr)) 471e2739ba6SAlexis Marboeuf PetscCallA(DMDestroy(pdm,ierr)) 472e2739ba6SAlexis Marboeuf if (numProc > 1) then 473d8606c27SBarry Smith PetscCallA(DMDestroy(dm,ierr)) 474e2739ba6SAlexis Marboeuf end if 475280aadf6SBlaise Bourdin 476280aadf6SBlaise Bourdin deallocate(pStartDepth) 477280aadf6SBlaise Bourdin deallocate(pEndDepth) 478280aadf6SBlaise Bourdin 479d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer,ierr)) 480d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 4812e8d78feSBlaise Bourdinend program ex62f90 482280aadf6SBlaise Bourdin 483280aadf6SBlaise Bourdin! /*TEST 484280aadf6SBlaise Bourdin! 485280aadf6SBlaise Bourdin! build: 486280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 487280aadf6SBlaise Bourdin! # 2D seq 488e2739ba6SAlexis Marboeuf! test: 489e2739ba6SAlexis Marboeuf! suffix: 0 490e2739ba6SAlexis 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 491e2739ba6SAlexis 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 492e2739ba6SAlexis Marboeuf! test: 493e2739ba6SAlexis Marboeuf! suffix: 1 494e2739ba6SAlexis 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 495e2739ba6SAlexis Marboeuf! 496e2739ba6SAlexis Marboeuf! test: 497e2739ba6SAlexis Marboeuf! suffix: 2 498e2739ba6SAlexis 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 499e2739ba6SAlexis 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 500e2739ba6SAlexis Marboeuf! test: 501e2739ba6SAlexis Marboeuf! suffix: 3 502e2739ba6SAlexis 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 503e2739ba6SAlexis Marboeuf! test: 504e2739ba6SAlexis Marboeuf! suffix: 4 505e2739ba6SAlexis 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 506e2739ba6SAlexis Marboeuf! test: 507e2739ba6SAlexis Marboeuf! suffix: 5 508e2739ba6SAlexis 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 509e2739ba6SAlexis Marboeuf! # 2D par 510e2739ba6SAlexis Marboeuf! test: 511e2739ba6SAlexis Marboeuf! suffix: 6 512e2739ba6SAlexis Marboeuf! nsize: 2 513e2739ba6SAlexis 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 514e2739ba6SAlexis 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 515e2739ba6SAlexis Marboeuf! test: 516e2739ba6SAlexis Marboeuf! suffix: 7 517e2739ba6SAlexis Marboeuf! nsize: 2 518e2739ba6SAlexis 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 519e2739ba6SAlexis Marboeuf! test: 520e2739ba6SAlexis Marboeuf! suffix: 8 521e2739ba6SAlexis Marboeuf! nsize: 2 522e2739ba6SAlexis 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 523e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 524e2739ba6SAlexis Marboeuf! test: 525e2739ba6SAlexis Marboeuf! suffix: 9 526e2739ba6SAlexis Marboeuf! nsize: 2 527e2739ba6SAlexis 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 528e2739ba6SAlexis Marboeuf! test: 529e2739ba6SAlexis Marboeuf! suffix: 10 530e2739ba6SAlexis Marboeuf! nsize: 2 531e2739ba6SAlexis 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 532e2739ba6SAlexis Marboeuf! test: 533e2739ba6SAlexis Marboeuf! # Something is now broken with parallel read/write for wahtever shape H is 534e2739ba6SAlexis Marboeuf! TODO: broken 535e2739ba6SAlexis Marboeuf! suffix: 11 536e2739ba6SAlexis Marboeuf! nsize: 2 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 539e2739ba6SAlexis Marboeuf! #3d seq 540e2739ba6SAlexis Marboeuf! test: 541e2739ba6SAlexis Marboeuf! suffix: 12 542e2739ba6SAlexis 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 543e2739ba6SAlexis Marboeuf! test: 544e2739ba6SAlexis Marboeuf! suffix: 13 545e2739ba6SAlexis 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 546e2739ba6SAlexis 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 547e2739ba6SAlexis Marboeuf! test: 548e2739ba6SAlexis Marboeuf! suffix: 14 549e2739ba6SAlexis 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 550e2739ba6SAlexis Marboeuf! test: 551e2739ba6SAlexis Marboeuf! suffix: 15 552e2739ba6SAlexis 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 553e2739ba6SAlexis 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 554e2739ba6SAlexis Marboeuf! #3d par 555e2739ba6SAlexis Marboeuf! test: 556e2739ba6SAlexis Marboeuf! suffix: 16 557e2739ba6SAlexis Marboeuf! nsize: 2 558e2739ba6SAlexis 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 559e2739ba6SAlexis Marboeuf! test: 560e2739ba6SAlexis Marboeuf! suffix: 17 561e2739ba6SAlexis Marboeuf! nsize: 2 562e2739ba6SAlexis 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 563e2739ba6SAlexis 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 564e2739ba6SAlexis Marboeuf! test: 565e2739ba6SAlexis Marboeuf! suffix: 18 566e2739ba6SAlexis Marboeuf! nsize: 2 567e2739ba6SAlexis 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 568e2739ba6SAlexis Marboeuf! test: 569e2739ba6SAlexis Marboeuf! suffix: 19 570e2739ba6SAlexis Marboeuf! nsize: 2 571e2739ba6SAlexis 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 572e2739ba6SAlexis 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 573e2739ba6SAlexis Marboeuf 574280aadf6SBlaise Bourdin! TEST*/ 575