1280aadf6SBlaise Bourdinprogram ex26f90 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 80d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 81280aadf6SBlaise Bourdin 82d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 83d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,numProc,ierr)) 84d8606c27SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-i",ifilename,flg,ierr)) 85280aadf6SBlaise Bourdin if (.not. flg) then 86*de401611SBlaise Bourdin SETERRA(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing input file name -i <input file name>") 87280aadf6SBlaise Bourdin end if 88d8606c27SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-o",ofilename,flg,ierr)) 89280aadf6SBlaise Bourdin if (.not. flg) then 90*de401611SBlaise Bourdin SETERRA(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing output file name -o <output file name>") 91280aadf6SBlaise Bourdin end if 92d8606c27SBarry 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 95*de401611SBlaise 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)) 100d8606c27SBarry Smith PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr)) 101d8606c27SBarry Smith PetscCallA(DMSetFromOptions(dm,ierr)) 102d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim,ierr)) 103d8606c27SBarry 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) 1236aad120cSJose 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 ! 133280aadf6SBlaise Bourdin ! "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 137280aadf6SBlaise Bourdin nodalVarName(1:numNodalVar) = ["U_x ","U_y ","Alpha"] 138280aadf6SBlaise Bourdin numZonalVar = 3 139280aadf6SBlaise Bourdin zonalVarName(1:numZonalVar) = ["Sigma_11","Sigma_22","Sigma_12"] 140280aadf6SBlaise Bourdin case(3) 141280aadf6SBlaise Bourdin numNodalVar = 4 142280aadf6SBlaise Bourdin nodalVarName(1:numNodalVar) = ["U_x ","U_y ","U_z ","Alpha"] 143280aadf6SBlaise Bourdin numZonalVar = 6 144280aadf6SBlaise Bourdin 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)) 149d8606c27SBarry Smith PetscCallA(expvp(exoid, "E", numZonalVar,ierr)) 150d8606c27SBarry Smith PetscCallA(expvan(exoid, "E", numZonalVar, zonalVarName,ierr)) 151d8606c27SBarry Smith PetscCallA(expvp(exoid, "N", numNodalVar,ierr)) 152d8606c27SBarry 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)) 170d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, "U",ierr)) 171d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, "Alpha",ierr)) 172d8606c27SBarry 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 183d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr)) 184d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr)) 185d8606c27SBarry Smith 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 188d8606c27SBarry Smith PetscCallA(DMGetLabelSize(dm, "Cell Sets", numCS, ierr)) 189d8606c27SBarry Smith PetscCallA(DMGetLabelIdIS(dm, "Cell Sets", csIS, ierr)) 190d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID, ierr)) 191280aadf6SBlaise Bourdin do set = 1,numCS 192d8606c27SBarry Smith PetscCallA(DMGetStratumSize(dm, "Cell Sets", csID(set), numCells,ierr)) 193d8606c27SBarry Smith PetscCallA(DMGetStratumIS(dm, "Cell Sets", csID(set), cellIS,ierr)) 194280aadf6SBlaise Bourdin if (numCells > 0) then 195280aadf6SBlaise Bourdin select case(sdim) 196280aadf6SBlaise Bourdin case(2) 197280aadf6SBlaise Bourdin dofs => dofS2D 198280aadf6SBlaise Bourdin case(3) 199280aadf6SBlaise Bourdin dofs => dofS3D 200280aadf6SBlaise Bourdin case default 201280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 202*de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 203280aadf6SBlaise Bourdin end select ! sdim 204280aadf6SBlaise Bourdin 205280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 206280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 207d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 208280aadf6SBlaise Bourdin nullify(closureA) 209d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr)) 210280aadf6SBlaise Bourdin select case(size(closureA)/2) 211280aadf6SBlaise Bourdin case(7) ! Tri 212280aadf6SBlaise Bourdin if (order == 1) then 213280aadf6SBlaise Bourdin dofU => dofUP1Tri 214280aadf6SBlaise Bourdin dofA => dofAP1Tri 215280aadf6SBlaise Bourdin else 216280aadf6SBlaise Bourdin dofU => dofUP2Tri 217280aadf6SBlaise Bourdin dofA => dofAP2Tri 218280aadf6SBlaise Bourdin end if 219280aadf6SBlaise Bourdin case(9) ! Quad 220280aadf6SBlaise Bourdin if (order == 1) then 221280aadf6SBlaise Bourdin dofU => dofUP1Quad 222280aadf6SBlaise Bourdin dofA => dofAP1Quad 223280aadf6SBlaise Bourdin else 224280aadf6SBlaise Bourdin dofU => dofUP2Quad 225280aadf6SBlaise Bourdin dofA => dofAP2Quad 226280aadf6SBlaise Bourdin end if 227280aadf6SBlaise Bourdin case(15) ! Tet 228280aadf6SBlaise Bourdin if (order == 1) then 229280aadf6SBlaise Bourdin dofU => dofUP1Tet 230280aadf6SBlaise Bourdin dofA => dofAP1Tet 231280aadf6SBlaise Bourdin else 232280aadf6SBlaise Bourdin dofU => dofUP2Tet 233280aadf6SBlaise Bourdin dofA => dofAP2Tet 234280aadf6SBlaise Bourdin end if 235280aadf6SBlaise Bourdin case(27) ! Hex 236280aadf6SBlaise Bourdin if (order == 1) then 237280aadf6SBlaise Bourdin dofU => dofUP1Hex 238280aadf6SBlaise Bourdin dofA => dofAP1Hex 239280aadf6SBlaise Bourdin else 240280aadf6SBlaise Bourdin dofU => dofUP2Hex 241280aadf6SBlaise Bourdin dofA => dofAP2Hex 242280aadf6SBlaise Bourdin end if 243280aadf6SBlaise Bourdin case default 244280aadf6SBlaise Bourdin write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2 245*de401611SBlaise Bourdin SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer) 246280aadf6SBlaise Bourdin end select 247d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr)) 248280aadf6SBlaise Bourdin do cell = 1,numCells! 249280aadf6SBlaise Bourdin nullify(closure) 250d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 251280aadf6SBlaise Bourdin do p = 1,size(closure),2 252280aadf6SBlaise Bourdin ! find the depth of p 253280aadf6SBlaise Bourdin do d = 1,sdim+1 254280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 255d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr)) 256d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr)) 257d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr)) 258d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr)) 259280aadf6SBlaise Bourdin end if ! closure(p) 260280aadf6SBlaise Bourdin end do ! d 261280aadf6SBlaise Bourdin end do ! p 262d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 263280aadf6SBlaise Bourdin end do ! cell 264d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 265d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 266280aadf6SBlaise Bourdin end if ! numCells 267280aadf6SBlaise Bourdin end do ! set 268d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 269d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 270d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section,ierr)) 271d8606c27SBarry Smith PetscCallA(DMSetLocalSection(dm, section,ierr)) 272d8606c27SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, "-dm_section_view",ierr)) 273d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section,ierr)) 274280aadf6SBlaise Bourdin 275d8606c27SBarry Smith PetscCallA(DMSetUseNatural(dm,PETSC_TRUE,ierr)) 276d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm,part,ierr)) 277d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part,ierr)) 278d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr)) 279280aadf6SBlaise Bourdin 280280aadf6SBlaise Bourdin if (numProc > 1) then 281d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr)) 282d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF,ierr)) 283d8606c27SBarry Smith PetscCallA(DMDestroy(dm,ierr)) 284280aadf6SBlaise Bourdin dm = pdm 285280aadf6SBlaise Bourdin end if 286d8606c27SBarry Smith PetscCallA(DMViewFromOptions(dm,PETSC_NULL_OPTIONS,"-dm_view",ierr)) 287280aadf6SBlaise Bourdin 288280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 289d8606c27SBarry Smith PetscCallA(DMCreateSubDM(dm, 1_kPI, fieldU, isU, dmU,ierr)) 290d8606c27SBarry Smith PetscCallA(DMCreateSubDM(dm, 1_kPI, fieldA, isA, dmA,ierr)) 291d8606c27SBarry Smith PetscCallA(DMCreateSubDM(dm, 1_kPI, fieldS, isS, dmS,ierr)) 292d8606c27SBarry Smith PetscCallA(DMCreateSubDM(dm, 2_kPI, fieldUA, isUA, dmUA,ierr)) 293280aadf6SBlaise Bourdin 294280aadf6SBlaise Bourdin !Create the exodus result file 295280aadf6SBlaise Bourdin allocate(dmList(2)) 296280aadf6SBlaise Bourdin dmList(1) = dmU; 297280aadf6SBlaise Bourdin dmList(2) = dmA; 298d8606c27SBarry Smith PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr)) 299280aadf6SBlaise Bourdin deallocate(dmList) 300280aadf6SBlaise Bourdin 301d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dm, X,ierr)) 302d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U,ierr)) 303d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A,ierr)) 304d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S,ierr)) 305d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA,ierr)) 306d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr)) 307280aadf6SBlaise Bourdin 308d8606c27SBarry Smith PetscCallA(PetscObjectSetName(U, "U",ierr)) 309d8606c27SBarry Smith PetscCallA(PetscObjectSetName(A, "Alpha",ierr)) 310d8606c27SBarry Smith PetscCallA(PetscObjectSetName(S, "Sigma",ierr)) 311d8606c27SBarry Smith PetscCallA(PetscObjectSetName(UA, "UAlpha",ierr)) 312d8606c27SBarry Smith PetscCallA(PetscObjectSetName(UA2, "UAlpha2",ierr)) 313d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR,ierr)) 314280aadf6SBlaise Bourdin 315280aadf6SBlaise 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 */ 316d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr)) 317d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr)) 318d8606c27SBarry Smith PetscCallA(VecGetArrayF90(UALoc, cval,ierr)) 319d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr)) 320d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr)) 321d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr)) 322280aadf6SBlaise Bourdin 323280aadf6SBlaise Bourdin do p = pStart,pEnd-1 324d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr)) 325280aadf6SBlaise Bourdin if (dofUA > 0) then 326d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr)) 327d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr)) 328280aadf6SBlaise Bourdin closureSize = size(xyz) 329280aadf6SBlaise Bourdin do i = 1,sdim 330280aadf6SBlaise Bourdin do j = 0, closureSize-1,sdim 331280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i) 332280aadf6SBlaise Bourdin end do 333280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) * sdim / closureSize; 334280aadf6SBlaise Bourdin cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2 335280aadf6SBlaise Bourdin end do 336d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr)) 337280aadf6SBlaise Bourdin end if 338280aadf6SBlaise Bourdin end do 339280aadf6SBlaise Bourdin 340d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr)) 341d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 342d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 343d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr)) 344280aadf6SBlaise Bourdin 345280aadf6SBlaise Bourdin !Update X 346d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr)) 347280aadf6SBlaise Bourdin ! Restrict to U and Alpha 348d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr)) 349d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr)) 350d8606c27SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, "-ua_vec_view",ierr)) 351d8606c27SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, "-u_vec_view",ierr)) 352d8606c27SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, "-a_vec_view",ierr)) 353280aadf6SBlaise Bourdin ! restrict to UA2 354d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr)) 355d8606c27SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, "-ua2_vec_view",ierr)) 356280aadf6SBlaise Bourdin 357280aadf6SBlaise Bourdin ! Writing nodal variables to ExodusII file 358d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU,0_kPI,time,ierr)) 359d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA,0_kPI,time,ierr)) 360280aadf6SBlaise Bourdin 361d8606c27SBarry Smith PetscCallA(VecView(U, viewer,ierr)) 362d8606c27SBarry Smith PetscCallA(VecView(A, viewer,ierr)) 363280aadf6SBlaise Bourdin 364280aadf6SBlaise Bourdin ! Saving U and Alpha in one shot. 365280aadf6SBlaise Bourdin ! For this, we need to cheat and change the Vec's name 366280aadf6SBlaise Bourdin ! Note that in the end we write variables one component at a time, 367280aadf6SBlaise Bourdin ! so that there is no real value in doing this 368d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr)) 369d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr)) 370d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec,ierr)) 371d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "U",ierr)) 372d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 373280aadf6SBlaise Bourdin 374280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 375d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 376d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer,ierr)) 377d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr)) 378d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr)) 379280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 380280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 381280aadf6SBlaise Bourdin end if 382d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr)) 383280aadf6SBlaise Bourdin 384280aadf6SBlaise Bourdin ! same thing with the UA2 Vec obtained from the superDM 385d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr)) 386d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec,ierr)) 387d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "U",ierr)) 388d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr)) 389d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 390280aadf6SBlaise Bourdin 391280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 392d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 393d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 394d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr)) 395d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr)) 396280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 397280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 398280aadf6SBlaise Bourdin end if 399d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr)) 400280aadf6SBlaise Bourdin 401280aadf6SBlaise Bourdin ! Building and saving Sigma 402280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 403280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 404280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 405d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr)) 406d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr)) 407d8606c27SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, "Cell Sets", csIS,ierr)) 408d8606c27SBarry Smith PetscCallA(DMGetLabelSize(dmS, "Cell Sets",numCS,ierr)) 409d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID,ierr)) 410280aadf6SBlaise Bourdin 411280aadf6SBlaise Bourdin do set = 1, numCS 412d8606c27SBarry Smith PetscCallA(DMGetStratumIS(dmS, "Cell Sets", csID(set), cellIS,ierr)) 413d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 414d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells,ierr)) 415280aadf6SBlaise Bourdin do cell = 1,numCells 416d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 417d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 418280aadf6SBlaise Bourdin cval(1) = rank 419280aadf6SBlaise Bourdin cval(2) = csID(set) 420280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 421280aadf6SBlaise Bourdin do c = 1, size(xyz),sdim 422280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 423280aadf6SBlaise Bourdin end do 424280aadf6SBlaise Bourdin cval(3) = cval(3) * sdim / size(xyz) 425d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr)) 426d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 427d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 428280aadf6SBlaise Bourdin end do 429d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 430d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 431280aadf6SBlaise Bourdin end do 432d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 433d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 434d8606c27SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, "-s_vec_view",ierr)) 435280aadf6SBlaise Bourdin 436280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 437d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr)) 438d8606c27SBarry Smith PetscCallA(VecView(S,viewer,ierr)) 439280aadf6SBlaise Bourdin 440280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 441d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr)) 442d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 443d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "Sigma",ierr)) 444d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 445d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr)) 446d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr)) 447280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 448280aadf6SBlaise Bourdin write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm 449280aadf6SBlaise Bourdin end if 450d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr)) 451280aadf6SBlaise Bourdin 452d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr)) 453d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr)) 454d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S,ierr)) 455d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A,ierr)) 456d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U,ierr)) 457d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dm, X,ierr)) 458d8606c27SBarry Smith PetscCallA(DMDestroy(dmU,ierr)) 459d8606c27SBarry Smith PetscCallA(ISDestroy(isU,ierr)) 460d8606c27SBarry Smith PetscCallA(DMDestroy(dmA,ierr)) 461d8606c27SBarry Smith PetscCallA(ISDestroy(isA,ierr)) 462d8606c27SBarry Smith PetscCallA(DMDestroy(dmS,ierr)) 463d8606c27SBarry Smith PetscCallA(ISDestroy(isS,ierr)) 464d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA,ierr)) 465d8606c27SBarry Smith PetscCallA(ISDestroy(isUA,ierr)) 466d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2,ierr)) 467d8606c27SBarry Smith PetscCallA(DMDestroy(dm,ierr)) 468280aadf6SBlaise Bourdin 469280aadf6SBlaise Bourdin deallocate(pStartDepth) 470280aadf6SBlaise Bourdin deallocate(pEndDepth) 471280aadf6SBlaise Bourdin 472d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer,ierr)) 473d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 474280aadf6SBlaise Bourdinend program ex26f90 475280aadf6SBlaise Bourdin 476280aadf6SBlaise Bourdin! /*TEST 477280aadf6SBlaise Bourdin! 478280aadf6SBlaise Bourdin! build: 479280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 480280aadf6SBlaise Bourdin! # 2D seq 481445019dbSMatthew G. Knepley! # test: 482445019dbSMatthew G. Knepley! # suffix: 0 483445019dbSMatthew G. Knepley! # 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 484445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 485445019dbSMatthew G. Knepley! # test: 486445019dbSMatthew G. Knepley! # suffix: 1 487445019dbSMatthew G. Knepley! # 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 488445019dbSMatthew G. Knepley! # 489445019dbSMatthew G. Knepley! # test: 490445019dbSMatthew G. Knepley! # suffix: 2 491445019dbSMatthew G. Knepley! # 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 492445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 493445019dbSMatthew G. Knepley! # test: 494445019dbSMatthew G. Knepley! # suffix: 3 495445019dbSMatthew G. Knepley! # 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 496445019dbSMatthew G. Knepley! # test: 497445019dbSMatthew G. Knepley! # suffix: 4 498445019dbSMatthew G. Knepley! # 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 499445019dbSMatthew G. Knepley! # test: 500445019dbSMatthew G. Knepley! # suffix: 5 501445019dbSMatthew G. Knepley! # 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 502445019dbSMatthew G. Knepley! # # 2D par 503445019dbSMatthew G. Knepley! # test: 504445019dbSMatthew G. Knepley! # suffix: 6 505445019dbSMatthew G. Knepley! # nsize: 2 506445019dbSMatthew G. Knepley! # 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 507445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 508445019dbSMatthew G. Knepley! # test: 509445019dbSMatthew G. Knepley! # suffix: 7 510445019dbSMatthew G. Knepley! # nsize: 2 511445019dbSMatthew G. Knepley! # 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 512445019dbSMatthew G. Knepley! # test: 513445019dbSMatthew G. Knepley! # suffix: 8 514445019dbSMatthew G. Knepley! # nsize: 2 515445019dbSMatthew G. Knepley! # 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 516445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 517445019dbSMatthew G. Knepley! # test: 518445019dbSMatthew G. Knepley! # suffix: 9 519445019dbSMatthew G. Knepley! # nsize: 2 520445019dbSMatthew G. Knepley! # 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 521445019dbSMatthew G. Knepley! # test: 522445019dbSMatthew G. Knepley! # suffix: 10 523445019dbSMatthew G. Knepley! # nsize: 2 524445019dbSMatthew G. Knepley! # 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 525445019dbSMatthew G. Knepley! # test: 526445019dbSMatthew G. Knepley! # # Something is now broken with parallel read/write for wahtever shape H is 527445019dbSMatthew G. Knepley! # TODO: broken 528445019dbSMatthew G. Knepley! # suffix: 11 529445019dbSMatthew G. Knepley! # nsize: 2 530445019dbSMatthew G. Knepley! # 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 531445019dbSMatthew G. Knepley! # 532445019dbSMatthew G. Knepley! # #3d seq 533445019dbSMatthew G. Knepley! # test: 534445019dbSMatthew G. Knepley! # suffix: 12 535445019dbSMatthew G. Knepley! # 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 536445019dbSMatthew G. Knepley! # test: 537445019dbSMatthew G. Knepley! # suffix: 13 538445019dbSMatthew G. Knepley! # 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 539445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 540445019dbSMatthew G. Knepley! # test: 541445019dbSMatthew G. Knepley! # suffix: 14 542445019dbSMatthew G. Knepley! # 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 543445019dbSMatthew G. Knepley! # test: 544445019dbSMatthew G. Knepley! # suffix: 15 545445019dbSMatthew G. Knepley! # 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 546445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 547445019dbSMatthew G. Knepley! # #3d par 548445019dbSMatthew G. Knepley! # test: 549445019dbSMatthew G. Knepley! # suffix: 16 550445019dbSMatthew G. Knepley! # nsize: 2 551445019dbSMatthew G. Knepley! # 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 552445019dbSMatthew G. Knepley! # test: 553445019dbSMatthew G. Knepley! # suffix: 17 554445019dbSMatthew G. Knepley! # nsize: 2 555445019dbSMatthew G. Knepley! # 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 556445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 557445019dbSMatthew G. Knepley! # test: 558445019dbSMatthew G. Knepley! # suffix: 18 559445019dbSMatthew G. Knepley! # nsize: 2 560445019dbSMatthew G. Knepley! # 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 561445019dbSMatthew G. Knepley! # test: 562445019dbSMatthew G. Knepley! # suffix: 19 563445019dbSMatthew G. Knepley! # nsize: 2 564445019dbSMatthew G. Knepley! # 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 565445019dbSMatthew G. Knepley! # #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints 566280aadf6SBlaise Bourdin! TEST*/ 567