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)) 88d8606c27SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-i",ifilename,flg,ierr)) 89280aadf6SBlaise Bourdin if (.not. flg) then 90280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing input file name -i <input file name>") 91280aadf6SBlaise Bourdin end if 92d8606c27SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-o",ofilename,flg,ierr)) 93280aadf6SBlaise Bourdin if (.not. flg) then 94280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing output file name -o <output file name>") 95280aadf6SBlaise Bourdin end if 96d8606c27SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-order",order,flg,ierr)) 97280aadf6SBlaise Bourdin if ((order > 2) .or. (order < 1)) then 98280aadf6SBlaise Bourdin write(IOBuffer,'("Unsupported polynomial order ", I2, " not in [1,2]")') order 99280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 100280aadf6SBlaise Bourdin end if 101280aadf6SBlaise Bourdin 102280aadf6SBlaise Bourdin ! Read the mesh in any supported format 103d8606c27SBarry Smith PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr)) 104e2739ba6SAlexis Marboeuf PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr)); 105e2739ba6SAlexis Marboeuf PetscCallA(DMSetFromOptions(dm,ierr)); 106d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim,ierr)) 107e2739ba6SAlexis Marboeuf PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,"-dm_view",ierr)); 108280aadf6SBlaise Bourdin 109280aadf6SBlaise Bourdin ! Create the exodus result file 110280aadf6SBlaise Bourdin 111a5b23f4aSJose E. Roman ! enable exodus debugging information 112d8606c27SBarry Smith PetscCallA(exopts(EXVRBS+EXDEBG,ierr)) 113280aadf6SBlaise Bourdin ! Create the exodus file 114d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr)) 115280aadf6SBlaise Bourdin ! The long way would be 116280aadf6SBlaise Bourdin ! 117d8606c27SBarry Smith ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr)) 118d8606c27SBarry Smith ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr)) 119d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr)) 120d8606c27SBarry Smith ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr)) 121280aadf6SBlaise Bourdin 122280aadf6SBlaise Bourdin ! set the mesh order 123d8606c27SBarry Smith PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr)) 124d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr)) 125280aadf6SBlaise Bourdin ! 126280aadf6SBlaise Bourdin ! Notice how the exodus file is actually NOT open at this point (exoid is -1) 127e2739ba6SAlexis Marboeuf ! Since we are overwritting the file (mode is FILE_MODE_WRITE), we are going to have to 128280aadf6SBlaise Bourdin ! write the geometry (the DM), which can only be done on a brand new file. 129280aadf6SBlaise Bourdin ! 130280aadf6SBlaise Bourdin 131280aadf6SBlaise Bourdin ! Save the geometry to the file, erasing all previous content 132d8606c27SBarry Smith PetscCallA(DMView(dm,viewer,ierr)) 133d8606c27SBarry Smith PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr)) 134280aadf6SBlaise Bourdin ! 135280aadf6SBlaise Bourdin ! Note how the exodus file is now open 136280aadf6SBlaise Bourdin ! 137280aadf6SBlaise Bourdin ! "Format" the exodus result file, i.e. allocate space for nodal and zonal variables 138280aadf6SBlaise Bourdin select case(sdim) 139280aadf6SBlaise Bourdin case(2) 140280aadf6SBlaise Bourdin numNodalVar = 3 141280aadf6SBlaise Bourdin nodalVarName(1:numNodalVar) = ["U_x ","U_y ","Alpha"] 142280aadf6SBlaise Bourdin numZonalVar = 3 143280aadf6SBlaise Bourdin zonalVarName(1:numZonalVar) = ["Sigma_11","Sigma_22","Sigma_12"] 144280aadf6SBlaise Bourdin case(3) 145280aadf6SBlaise Bourdin numNodalVar = 4 146280aadf6SBlaise Bourdin nodalVarName(1:numNodalVar) = ["U_x ","U_y ","U_z ","Alpha"] 147280aadf6SBlaise Bourdin numZonalVar = 6 148280aadf6SBlaise Bourdin zonalVarName(1:numZonalVar) = ["Sigma_11","Sigma_22","Sigma_33","Sigma_23","Sigma_13","Sigma_12"] 149280aadf6SBlaise Bourdin case default 150280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 151280aadf6SBlaise Bourdin end select 152d8606c27SBarry Smith PetscCallA(PetscViewerExodusIIGetId(viewer,exoid,ierr)) 153d8606c27SBarry Smith PetscCallA(expvp(exoid, "E", numZonalVar,ierr)) 154d8606c27SBarry Smith PetscCallA(expvan(exoid, "E", numZonalVar, zonalVarName,ierr)) 155d8606c27SBarry Smith PetscCallA(expvp(exoid, "N", numNodalVar,ierr)) 156d8606c27SBarry Smith PetscCallA(expvan(exoid, "N", numNodalVar, nodalVarName,ierr)) 157d8606c27SBarry Smith PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK,numCS,PETSC_NULL_REAL,sjunk,ierr)) 158280aadf6SBlaise Bourdin 159280aadf6SBlaise Bourdin ! An exodusII truth table specifies which fields are saved at which time step 160280aadf6SBlaise Bourdin ! It speeds up I/O but reserving space for fields in the file ahead of time. 161280aadf6SBlaise Bourdin allocate(truthtable(numCS,numZonalVar)) 162280aadf6SBlaise Bourdin truthtable = .true. 163d8606c27SBarry Smith PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr)) 164280aadf6SBlaise Bourdin deallocate(truthtable) 165280aadf6SBlaise Bourdin 166280aadf6SBlaise Bourdin ! Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 167280aadf6SBlaise Bourdin do step = 1,numstep 168d8606c27SBarry Smith PetscCallA(exptim(exoid,step,Real(step,kind=kPR),ierr)) 169280aadf6SBlaise Bourdin end do 170280aadf6SBlaise Bourdin 171d8606c27SBarry Smith PetscCallA(PetscObjectGetComm(dm,comm,ierr)) 172d8606c27SBarry Smith PetscCallA(PetscSectionCreate(comm, section,ierr)) 173d8606c27SBarry Smith PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr)) 174d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, "U",ierr)) 175d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, "Alpha",ierr)) 176d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldS, "Sigma",ierr)) 177d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dm, pStart, pEnd,ierr)) 178d8606c27SBarry Smith PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr)) 179280aadf6SBlaise Bourdin 180280aadf6SBlaise Bourdin allocate(pStartDepth(sdim+1)) 181280aadf6SBlaise Bourdin allocate(pEndDepth(sdim+1)) 182280aadf6SBlaise Bourdin do d = 1, sdim+1 183d8606c27SBarry Smith PetscCallA(DMPlexGetDepthStratum(dm, d-1, pStartDepth(d), pEndDepth(d),ierr)) 184280aadf6SBlaise Bourdin end do 185280aadf6SBlaise Bourdin 186280aadf6SBlaise Bourdin ! Vector field U, Scalar field Alpha, Tensor field Sigma 187e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr)); 188e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr)); 189e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr)); 190280aadf6SBlaise Bourdin 191280aadf6SBlaise Bourdin ! Going through cell sets then cells, and setting up storage for the sections 192d8606c27SBarry Smith PetscCallA(DMGetLabelSize(dm, "Cell Sets", numCS, ierr)) 193d8606c27SBarry Smith PetscCallA(DMGetLabelIdIS(dm, "Cell Sets", csIS, ierr)) 194d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID, ierr)) 195280aadf6SBlaise Bourdin do set = 1,numCS 196*8f2079feSBlaise Bourdin !!! This is pointless but will avoid a warning with gfortran and -Werror=maybe-uninitialized 197*8f2079feSBlaise Bourdin nullify(dofA) 198*8f2079feSBlaise Bourdin nullify(dofU) 199*8f2079feSBlaise Bourdin nullify(dofS) 200d8606c27SBarry Smith PetscCallA(DMGetStratumSize(dm, "Cell Sets", csID(set), numCells,ierr)) 201d8606c27SBarry Smith PetscCallA(DMGetStratumIS(dm, "Cell Sets", csID(set), cellIS,ierr)) 202280aadf6SBlaise Bourdin if (numCells > 0) then 203280aadf6SBlaise Bourdin select case(sdim) 204280aadf6SBlaise Bourdin case(2) 205280aadf6SBlaise Bourdin dofs => dofS2D 206280aadf6SBlaise Bourdin case(3) 207280aadf6SBlaise Bourdin dofs => dofS3D 208280aadf6SBlaise Bourdin case default 209280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 210280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 211280aadf6SBlaise Bourdin end select ! sdim 212280aadf6SBlaise Bourdin 213280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 214280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 215d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 216280aadf6SBlaise Bourdin nullify(closureA) 217d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr)) 218280aadf6SBlaise Bourdin select case(size(closureA)/2) 219280aadf6SBlaise Bourdin case(7) ! Tri 220280aadf6SBlaise Bourdin if (order == 1) then 221280aadf6SBlaise Bourdin dofU => dofUP1Tri 222280aadf6SBlaise Bourdin dofA => dofAP1Tri 223280aadf6SBlaise Bourdin else 224280aadf6SBlaise Bourdin dofU => dofUP2Tri 225280aadf6SBlaise Bourdin dofA => dofAP2Tri 226280aadf6SBlaise Bourdin end if 227280aadf6SBlaise Bourdin case(9) ! Quad 228280aadf6SBlaise Bourdin if (order == 1) then 229280aadf6SBlaise Bourdin dofU => dofUP1Quad 230280aadf6SBlaise Bourdin dofA => dofAP1Quad 231280aadf6SBlaise Bourdin else 232280aadf6SBlaise Bourdin dofU => dofUP2Quad 233280aadf6SBlaise Bourdin dofA => dofAP2Quad 234280aadf6SBlaise Bourdin end if 235280aadf6SBlaise Bourdin case(15) ! Tet 236280aadf6SBlaise Bourdin if (order == 1) then 237280aadf6SBlaise Bourdin dofU => dofUP1Tet 238280aadf6SBlaise Bourdin dofA => dofAP1Tet 239280aadf6SBlaise Bourdin else 240280aadf6SBlaise Bourdin dofU => dofUP2Tet 241280aadf6SBlaise Bourdin dofA => dofAP2Tet 242280aadf6SBlaise Bourdin end if 243280aadf6SBlaise Bourdin case(27) ! Hex 244280aadf6SBlaise Bourdin if (order == 1) then 245280aadf6SBlaise Bourdin dofU => dofUP1Hex 246280aadf6SBlaise Bourdin dofA => dofAP1Hex 247280aadf6SBlaise Bourdin else 248280aadf6SBlaise Bourdin dofU => dofUP2Hex 249280aadf6SBlaise Bourdin dofA => dofAP2Hex 250280aadf6SBlaise Bourdin end if 251280aadf6SBlaise Bourdin case default 252280aadf6SBlaise Bourdin write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2 253280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer) 254280aadf6SBlaise Bourdin end select 255d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr)) 256280aadf6SBlaise Bourdin do cell = 1,numCells! 257280aadf6SBlaise Bourdin nullify(closure) 258d8606c27SBarry Smith PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 259280aadf6SBlaise Bourdin do p = 1,size(closure),2 260280aadf6SBlaise Bourdin ! find the depth of p 261280aadf6SBlaise Bourdin do d = 1,sdim+1 262280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 263d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr)) 264d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr)) 265d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr)) 266d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr)) 267280aadf6SBlaise Bourdin end if ! closure(p) 268280aadf6SBlaise Bourdin end do ! d 269280aadf6SBlaise Bourdin end do ! p 270d8606c27SBarry Smith PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr)) 271280aadf6SBlaise Bourdin end do ! cell 272d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 273d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 274280aadf6SBlaise Bourdin end if ! numCells 275280aadf6SBlaise Bourdin end do ! set 276d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 277d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 278d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section,ierr)) 279d8606c27SBarry Smith PetscCallA(DMSetLocalSection(dm, section,ierr)) 280d8606c27SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, "-dm_section_view",ierr)) 281d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section,ierr)) 282280aadf6SBlaise Bourdin 283d8606c27SBarry Smith PetscCallA(DMSetUseNatural(dm,PETSC_TRUE,ierr)) 284d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm,part,ierr)) 285d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part,ierr)) 286d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr)) 287280aadf6SBlaise Bourdin 288280aadf6SBlaise Bourdin if (numProc > 1) then 289d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr)) 290d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF,ierr)) 291e2739ba6SAlexis Marboeuf else 292e2739ba6SAlexis Marboeuf pdm = dm 293280aadf6SBlaise Bourdin end if 294e2739ba6SAlexis Marboeuf PetscCallA(DMViewFromOptions(pdm,PETSC_NULL_OPTIONS,"-dm_view",ierr)) 295280aadf6SBlaise Bourdin 296280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 297e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldU, isU, dmU,ierr)) 298e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldA, isA, dmA,ierr)) 299e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldS, isS, dmS,ierr)) 300e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA,ierr)) 301280aadf6SBlaise Bourdin 302280aadf6SBlaise Bourdin !Create the exodus result file 303280aadf6SBlaise Bourdin allocate(dmList(2)) 304280aadf6SBlaise Bourdin dmList(1) = dmU; 305280aadf6SBlaise Bourdin dmList(2) = dmA; 306d8606c27SBarry Smith PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr)) 307280aadf6SBlaise Bourdin deallocate(dmList) 308280aadf6SBlaise Bourdin 309e2739ba6SAlexis Marboeuf PetscCallA(DMGetGlobalVector(pdm, X,ierr)) 310d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U,ierr)) 311d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A,ierr)) 312d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S,ierr)) 313d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA,ierr)) 314d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr)) 315280aadf6SBlaise Bourdin 316d8606c27SBarry Smith PetscCallA(PetscObjectSetName(U, "U",ierr)) 317d8606c27SBarry Smith PetscCallA(PetscObjectSetName(A, "Alpha",ierr)) 318d8606c27SBarry Smith PetscCallA(PetscObjectSetName(S, "Sigma",ierr)) 319d8606c27SBarry Smith PetscCallA(PetscObjectSetName(UA, "UAlpha",ierr)) 320d8606c27SBarry Smith PetscCallA(PetscObjectSetName(UA2, "UAlpha2",ierr)) 321d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR,ierr)) 322280aadf6SBlaise Bourdin 323280aadf6SBlaise 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 */ 324d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr)) 325d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr)) 326d8606c27SBarry Smith PetscCallA(VecGetArrayF90(UALoc, cval,ierr)) 327d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr)) 328d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr)) 329d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr)) 330280aadf6SBlaise Bourdin 331280aadf6SBlaise Bourdin do p = pStart,pEnd-1 332d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr)) 333280aadf6SBlaise Bourdin if (dofUA > 0) then 334d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr)) 335d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr)) 336280aadf6SBlaise Bourdin closureSize = size(xyz) 337280aadf6SBlaise Bourdin do i = 1,sdim 338280aadf6SBlaise Bourdin do j = 0, closureSize-1,sdim 339280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i) 340280aadf6SBlaise Bourdin end do 341280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) * sdim / closureSize; 342280aadf6SBlaise Bourdin cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2 343280aadf6SBlaise Bourdin end do 344d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr)) 345280aadf6SBlaise Bourdin end if 346280aadf6SBlaise Bourdin end do 347280aadf6SBlaise Bourdin 348d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr)) 349d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 350d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 351d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr)) 352280aadf6SBlaise Bourdin 353280aadf6SBlaise Bourdin !Update X 354d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr)) 355280aadf6SBlaise Bourdin ! Restrict to U and Alpha 356d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr)) 357d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr)) 358d8606c27SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, "-ua_vec_view",ierr)) 359d8606c27SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, "-u_vec_view",ierr)) 360d8606c27SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, "-a_vec_view",ierr)) 361280aadf6SBlaise Bourdin ! restrict to UA2 362d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr)) 363d8606c27SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, "-ua2_vec_view",ierr)) 364280aadf6SBlaise Bourdin 365e2739ba6SAlexis Marboeuf ! Getting Natural Vec 366d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr)) 367d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr)) 368280aadf6SBlaise Bourdin 369d8606c27SBarry Smith PetscCallA(VecView(U, viewer,ierr)) 370d8606c27SBarry Smith PetscCallA(VecView(A, viewer,ierr)) 371280aadf6SBlaise Bourdin 372280aadf6SBlaise Bourdin !Saving U and Alpha in one shot. 373280aadf6SBlaise Bourdin !For this, we need to cheat and change the Vec's name 374280aadf6SBlaise Bourdin !Note that in the end we write variables one component at a time, 375280aadf6SBlaise Bourdin !so that there is no real value in doing this 376d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr)) 377d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr)) 378d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec,ierr)) 379d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "U",ierr)) 380d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 381280aadf6SBlaise Bourdin 382280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 383d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 384d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer,ierr)) 385d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr)) 386d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr)) 387280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 388280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 389280aadf6SBlaise Bourdin end if 390d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr)) 391280aadf6SBlaise Bourdin 392280aadf6SBlaise Bourdin ! same thing with the UA2 Vec obtained from the superDM 393d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr)) 394d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec,ierr)) 395d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "U",ierr)) 396d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr)) 397d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 398280aadf6SBlaise Bourdin 399280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 400d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 401d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 402d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr)) 403d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr)) 404280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 405280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 406280aadf6SBlaise Bourdin end if 407d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr)) 408280aadf6SBlaise Bourdin 409280aadf6SBlaise Bourdin ! Building and saving Sigma 410280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 411280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 412280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 413d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr)) 414d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr)) 415d8606c27SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, "Cell Sets", csIS,ierr)) 416d8606c27SBarry Smith PetscCallA(DMGetLabelSize(dmS, "Cell Sets",numCS,ierr)) 417d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID,ierr)) 418280aadf6SBlaise Bourdin 419280aadf6SBlaise Bourdin do set = 1, numCS 420d8606c27SBarry Smith PetscCallA(DMGetStratumIS(dmS, "Cell Sets", csID(set), cellIS,ierr)) 421d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 422d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells,ierr)) 423280aadf6SBlaise Bourdin do cell = 1,numCells 424d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 425d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 426280aadf6SBlaise Bourdin cval(1) = rank 427280aadf6SBlaise Bourdin cval(2) = csID(set) 428280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 429280aadf6SBlaise Bourdin do c = 1, size(xyz),sdim 430280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 431280aadf6SBlaise Bourdin end do 432280aadf6SBlaise Bourdin cval(3) = cval(3) * sdim / size(xyz) 433d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr)) 434d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 435d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 436280aadf6SBlaise Bourdin end do 437d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 438d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 439280aadf6SBlaise Bourdin end do 440d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 441d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 442d8606c27SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, "-s_vec_view",ierr)) 443280aadf6SBlaise Bourdin 444280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 445d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr)) 446d8606c27SBarry Smith PetscCallA(VecView(S,viewer,ierr)) 447280aadf6SBlaise Bourdin 448280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 449d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr)) 450d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 451d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "Sigma",ierr)) 452d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 453d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr)) 454d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr)) 455280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 456280aadf6SBlaise Bourdin write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm 457280aadf6SBlaise Bourdin end if 458d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr)) 459280aadf6SBlaise Bourdin 460d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr)) 461d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr)) 462d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S,ierr)) 463d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A,ierr)) 464d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U,ierr)) 465e2739ba6SAlexis Marboeuf PetscCallA(DMRestoreGlobalVector(pdm, X,ierr)) 466d8606c27SBarry Smith PetscCallA(DMDestroy(dmU,ierr)) 467d8606c27SBarry Smith PetscCallA(ISDestroy(isU,ierr)) 468d8606c27SBarry Smith PetscCallA(DMDestroy(dmA,ierr)) 469d8606c27SBarry Smith PetscCallA(ISDestroy(isA,ierr)) 470d8606c27SBarry Smith PetscCallA(DMDestroy(dmS,ierr)) 471d8606c27SBarry Smith PetscCallA(ISDestroy(isS,ierr)) 472d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA,ierr)) 473d8606c27SBarry Smith PetscCallA(ISDestroy(isUA,ierr)) 474d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2,ierr)) 475e2739ba6SAlexis Marboeuf PetscCallA(DMDestroy(pdm,ierr)) 476e2739ba6SAlexis Marboeuf if (numProc > 1) then 477d8606c27SBarry Smith PetscCallA(DMDestroy(dm,ierr)) 478e2739ba6SAlexis Marboeuf end if 479280aadf6SBlaise Bourdin 480280aadf6SBlaise Bourdin deallocate(pStartDepth) 481280aadf6SBlaise Bourdin deallocate(pEndDepth) 482280aadf6SBlaise Bourdin 483d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer,ierr)) 484d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 4852e8d78feSBlaise Bourdinend program ex62f90 486280aadf6SBlaise Bourdin 487280aadf6SBlaise Bourdin! /*TEST 488280aadf6SBlaise Bourdin! 489280aadf6SBlaise Bourdin! build: 490280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 491280aadf6SBlaise Bourdin! # 2D seq 492e2739ba6SAlexis Marboeuf! test: 493e2739ba6SAlexis Marboeuf! suffix: 0 494e2739ba6SAlexis 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 495e2739ba6SAlexis 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 496e2739ba6SAlexis Marboeuf! test: 497e2739ba6SAlexis Marboeuf! suffix: 1 498e2739ba6SAlexis 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 499e2739ba6SAlexis Marboeuf! 500e2739ba6SAlexis Marboeuf! test: 501e2739ba6SAlexis Marboeuf! suffix: 2 502e2739ba6SAlexis 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 503e2739ba6SAlexis 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 504e2739ba6SAlexis Marboeuf! test: 505e2739ba6SAlexis Marboeuf! suffix: 3 506e2739ba6SAlexis 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 507e2739ba6SAlexis Marboeuf! test: 508e2739ba6SAlexis Marboeuf! suffix: 4 509e2739ba6SAlexis 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 510e2739ba6SAlexis Marboeuf! test: 511e2739ba6SAlexis Marboeuf! suffix: 5 512e2739ba6SAlexis 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 513e2739ba6SAlexis Marboeuf! # 2D par 514e2739ba6SAlexis Marboeuf! test: 515e2739ba6SAlexis Marboeuf! suffix: 6 516e2739ba6SAlexis Marboeuf! nsize: 2 517e2739ba6SAlexis 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 518e2739ba6SAlexis 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 519e2739ba6SAlexis Marboeuf! test: 520e2739ba6SAlexis Marboeuf! suffix: 7 521e2739ba6SAlexis Marboeuf! nsize: 2 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! test: 524e2739ba6SAlexis Marboeuf! suffix: 8 525e2739ba6SAlexis Marboeuf! nsize: 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: invalid dimension ID or name 528e2739ba6SAlexis Marboeuf! test: 529e2739ba6SAlexis Marboeuf! suffix: 9 530e2739ba6SAlexis Marboeuf! nsize: 2 531e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 532e2739ba6SAlexis Marboeuf! test: 533e2739ba6SAlexis Marboeuf! suffix: 10 534e2739ba6SAlexis Marboeuf! nsize: 2 535e2739ba6SAlexis 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 536e2739ba6SAlexis Marboeuf! test: 537e2739ba6SAlexis Marboeuf! # Something is now broken with parallel read/write for wahtever shape H is 538e2739ba6SAlexis Marboeuf! TODO: broken 539e2739ba6SAlexis Marboeuf! suffix: 11 540e2739ba6SAlexis Marboeuf! nsize: 2 541e2739ba6SAlexis 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 542e2739ba6SAlexis Marboeuf 543e2739ba6SAlexis Marboeuf! #3d seq 544e2739ba6SAlexis Marboeuf! test: 545e2739ba6SAlexis Marboeuf! suffix: 12 546e2739ba6SAlexis 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 547e2739ba6SAlexis Marboeuf! test: 548e2739ba6SAlexis Marboeuf! suffix: 13 549e2739ba6SAlexis 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 550e2739ba6SAlexis 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 551e2739ba6SAlexis Marboeuf! test: 552e2739ba6SAlexis Marboeuf! suffix: 14 553e2739ba6SAlexis 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 554e2739ba6SAlexis Marboeuf! test: 555e2739ba6SAlexis Marboeuf! suffix: 15 556e2739ba6SAlexis 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 557e2739ba6SAlexis 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 558e2739ba6SAlexis Marboeuf! #3d par 559e2739ba6SAlexis Marboeuf! test: 560e2739ba6SAlexis Marboeuf! suffix: 16 561e2739ba6SAlexis Marboeuf! nsize: 2 562e2739ba6SAlexis 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 563e2739ba6SAlexis Marboeuf! test: 564e2739ba6SAlexis Marboeuf! suffix: 17 565e2739ba6SAlexis Marboeuf! nsize: 2 566e2739ba6SAlexis 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 567e2739ba6SAlexis 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 568e2739ba6SAlexis Marboeuf! test: 569e2739ba6SAlexis Marboeuf! suffix: 18 570e2739ba6SAlexis Marboeuf! nsize: 2 571e2739ba6SAlexis Marboeuf! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2 572e2739ba6SAlexis Marboeuf! test: 573e2739ba6SAlexis Marboeuf! suffix: 19 574e2739ba6SAlexis Marboeuf! nsize: 2 575e2739ba6SAlexis 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 576e2739ba6SAlexis 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 577e2739ba6SAlexis Marboeuf 578280aadf6SBlaise Bourdin! TEST*/ 579