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 17e2739ba6SAlexis Marboeuf type(tPetscSection) :: section, rootSection, leafSection 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 73e2739ba6SAlexis Marboeuf type(tPetscSF) :: migrationSF, natSF, natPointSF, natPointSFInv 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 86e2739ba6SAlexis Marboeuf PetscCallA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 87e2739ba6SAlexis Marboeuf PetscCallA(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)) 104d8606c27SBarry Smith PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr)) 105d8606c27SBarry Smith PetscCallA(DMSetFromOptions(dm,ierr)) 106d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, sdim,ierr)) 107d8606c27SBarry Smith 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) 127d5b43468SJose E. Roman ! Since we are overwriting 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(DMSetUseNatural(dm, PETSC_TRUE, ierr)) 172d8606c27SBarry Smith PetscCallA(DMPlexGetPartitioner(dm,part,ierr)) 173d8606c27SBarry Smith PetscCallA(PetscPartitionerSetFromOptions(part,ierr)) 174d8606c27SBarry Smith PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr)) 175280aadf6SBlaise Bourdin 176280aadf6SBlaise Bourdin if (numProc > 1) then 177d8606c27SBarry Smith PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr)) 178d8606c27SBarry Smith PetscCallA(PetscSFDestroy(migrationSF,ierr)) 179e2739ba6SAlexis Marboeuf else 180e2739ba6SAlexis Marboeuf pdm = dm 181280aadf6SBlaise Bourdin end if 182e2739ba6SAlexis Marboeuf PetscCallA(DMViewFromOptions(pdm,PETSC_NULL_OPTIONS,"-dm_view",ierr)) 183280aadf6SBlaise Bourdin 184e2739ba6SAlexis Marboeuf PetscCallA(PetscObjectGetComm(pdm,comm,ierr)) 185d8606c27SBarry Smith PetscCallA(PetscSectionCreate(comm, section,ierr)) 186d8606c27SBarry Smith PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr)) 187d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldU, "U",ierr)) 188d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldA, "Alpha",ierr)) 189d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldName(section, fieldS, "Sigma",ierr)) 190e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetChart(pdm, pStart, pEnd,ierr)) 191d8606c27SBarry Smith PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr)) 192280aadf6SBlaise Bourdin 193280aadf6SBlaise Bourdin allocate(pStartDepth(sdim+1)) 194280aadf6SBlaise Bourdin allocate(pEndDepth(sdim+1)) 195280aadf6SBlaise Bourdin do d = 1, sdim+1 196e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetDepthStratum(pdm, d-1, pStartDepth(d), pEndDepth(d),ierr)) 197280aadf6SBlaise Bourdin end do 198280aadf6SBlaise Bourdin 199280aadf6SBlaise Bourdin ! Vector field U, Scalar field Alpha, Tensor field Sigma 200e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr)); 201e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr)); 202e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr)); 203280aadf6SBlaise Bourdin 204280aadf6SBlaise Bourdin ! Going through cell sets then cells, and setting up storage for the sections 205e2739ba6SAlexis Marboeuf PetscCallA(DMGetLabelSize(pdm, "Cell Sets", numCS, ierr)) 206e2739ba6SAlexis Marboeuf PetscCallA(DMGetLabelIdIS(pdm, "Cell Sets", csIS, ierr)) 207d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID, ierr)) 208280aadf6SBlaise Bourdin do set = 1,numCS 209e2739ba6SAlexis Marboeuf PetscCallA(DMGetStratumSize(pdm, "Cell Sets", csID(set), numCells,ierr)) 210e2739ba6SAlexis Marboeuf PetscCallA(DMGetStratumIS(pdm, "Cell Sets", csID(set), cellIS,ierr)) 211280aadf6SBlaise Bourdin if (numCells > 0) then 212280aadf6SBlaise Bourdin select case(sdim) 213280aadf6SBlaise Bourdin case(2) 214280aadf6SBlaise Bourdin dofs => dofS2D 215280aadf6SBlaise Bourdin case(3) 216280aadf6SBlaise Bourdin dofs => dofS3D 217280aadf6SBlaise Bourdin case default 218280aadf6SBlaise Bourdin write(IOBuffer,'("No layout for dimension ",I2)') sdim 219280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer) 220280aadf6SBlaise Bourdin end select ! sdim 221280aadf6SBlaise Bourdin 222280aadf6SBlaise Bourdin ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 223280aadf6SBlaise Bourdin ! It will not be enough to identify more exotic elements like pyramid or prisms... */ 224d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 225280aadf6SBlaise Bourdin nullify(closureA) 226e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetTransitiveClosure(pdm,cellID(1), PETSC_TRUE, closureA,ierr)) 227280aadf6SBlaise Bourdin select case(size(closureA)/2) 228280aadf6SBlaise Bourdin case(7) ! Tri 229280aadf6SBlaise Bourdin if (order == 1) then 230280aadf6SBlaise Bourdin dofU => dofUP1Tri 231280aadf6SBlaise Bourdin dofA => dofAP1Tri 232280aadf6SBlaise Bourdin else 233280aadf6SBlaise Bourdin dofU => dofUP2Tri 234280aadf6SBlaise Bourdin dofA => dofAP2Tri 235280aadf6SBlaise Bourdin end if 236280aadf6SBlaise Bourdin case(9) ! Quad 237280aadf6SBlaise Bourdin if (order == 1) then 238280aadf6SBlaise Bourdin dofU => dofUP1Quad 239280aadf6SBlaise Bourdin dofA => dofAP1Quad 240280aadf6SBlaise Bourdin else 241280aadf6SBlaise Bourdin dofU => dofUP2Quad 242280aadf6SBlaise Bourdin dofA => dofAP2Quad 243280aadf6SBlaise Bourdin end if 244280aadf6SBlaise Bourdin case(15) ! Tet 245280aadf6SBlaise Bourdin if (order == 1) then 246280aadf6SBlaise Bourdin dofU => dofUP1Tet 247280aadf6SBlaise Bourdin dofA => dofAP1Tet 248280aadf6SBlaise Bourdin else 249280aadf6SBlaise Bourdin dofU => dofUP2Tet 250280aadf6SBlaise Bourdin dofA => dofAP2Tet 251280aadf6SBlaise Bourdin end if 252280aadf6SBlaise Bourdin case(27) ! Hex 253280aadf6SBlaise Bourdin if (order == 1) then 254280aadf6SBlaise Bourdin dofU => dofUP1Hex 255280aadf6SBlaise Bourdin dofA => dofAP1Hex 256280aadf6SBlaise Bourdin else 257280aadf6SBlaise Bourdin dofU => dofUP2Hex 258280aadf6SBlaise Bourdin dofA => dofAP2Hex 259280aadf6SBlaise Bourdin end if 260280aadf6SBlaise Bourdin case default 261280aadf6SBlaise Bourdin write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2 262280aadf6SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer) 263280aadf6SBlaise Bourdin end select 264e2739ba6SAlexis Marboeuf PetscCallA(DMPlexRestoreTransitiveClosure(pdm, cellID(1), PETSC_TRUE,closureA,ierr)) 265280aadf6SBlaise Bourdin do cell = 1,numCells! 266280aadf6SBlaise Bourdin nullify(closure) 267e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetTransitiveClosure(pdm, cellID(cell), PETSC_TRUE, closure,ierr)) 268280aadf6SBlaise Bourdin do p = 1,size(closure),2 269280aadf6SBlaise Bourdin ! find the depth of p 270280aadf6SBlaise Bourdin do d = 1,sdim+1 271280aadf6SBlaise Bourdin if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then 272d8606c27SBarry Smith PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr)) 273d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr)) 274d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr)) 275d8606c27SBarry Smith PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr)) 276280aadf6SBlaise Bourdin end if ! closure(p) 277280aadf6SBlaise Bourdin end do ! d 278280aadf6SBlaise Bourdin end do ! p 279e2739ba6SAlexis Marboeuf PetscCallA(DMPlexRestoreTransitiveClosure(pdm, cellID(cell), PETSC_TRUE, closure,ierr)) 280280aadf6SBlaise Bourdin end do ! cell 281d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 282d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 283280aadf6SBlaise Bourdin end if ! numCells 284280aadf6SBlaise Bourdin end do ! set 285d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 286d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 287d8606c27SBarry Smith PetscCallA(PetscSectionSetUp(section,ierr)) 288e2739ba6SAlexis Marboeuf PetscCallA(DMSetLocalSection(pdm, section,ierr)) 289d8606c27SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, "-dm_section_view",ierr)) 290d8606c27SBarry Smith PetscCallA(PetscSectionDestroy(section,ierr)) 291280aadf6SBlaise Bourdin 292e2739ba6SAlexis Marboeuf ! Creating section on the sequential DM + creating the GlobalToNatural SF 293e2739ba6SAlexis Marboeuf if (numProc > 1) then 294e2739ba6SAlexis Marboeuf PetscCallA(DMPlexGetMigrationSF(pdm, natPointSF, ierr)) 295e2739ba6SAlexis Marboeuf PetscCallA(DMGetLocalSection(pdm, rootSection, ierr)) 296e2739ba6SAlexis Marboeuf PetscCallA(PetscSFCreateInverseSF(natPointSF, natPointSFInv, ierr)) 297e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionCreate(PETSC_COMM_WORLD, leafSection, ierr)) 298e2739ba6SAlexis Marboeuf PetscCallA(PetscSFDistributeSection(natPointSFInv, rootSection, PETSC_NULL_INTEGER, leafSection, ierr)) 299e2739ba6SAlexis Marboeuf PetscCallA(DMSetLocalSection(dm, leafSection, ierr)) 300e2739ba6SAlexis Marboeuf PetscCallA(DMPlexCreateGlobalToNaturalSF(pdm, leafSection, natPointSF, natSF, ierr)) 301e2739ba6SAlexis Marboeuf PetscCallA(PetscSFDestroy(natPointSFInv, ierr)) 302e2739ba6SAlexis Marboeuf PetscCallA(PetscSectionDestroy(leafSection,ierr)) 303e2739ba6SAlexis Marboeuf PetscCallA(DMSetNaturalSF(pdm, natSF, ierr)) 304e2739ba6SAlexis Marboeuf PetscCallA(PetscObjectDereference(natSF, ierr)) 305e2739ba6SAlexis Marboeuf end if 306e2739ba6SAlexis Marboeuf 307280aadf6SBlaise Bourdin ! Get DM and IS for each field of dm 308e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldU, isU, dmU,ierr)) 309e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldA, isA, dmA,ierr)) 310e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldS, isS, dmS,ierr)) 311e2739ba6SAlexis Marboeuf PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA,ierr)) 312280aadf6SBlaise Bourdin 313280aadf6SBlaise Bourdin !Create the exodus result file 314280aadf6SBlaise Bourdin allocate(dmList(2)) 315280aadf6SBlaise Bourdin dmList(1) = dmU; 316280aadf6SBlaise Bourdin dmList(2) = dmA; 317d8606c27SBarry Smith PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr)) 318280aadf6SBlaise Bourdin deallocate(dmList) 319280aadf6SBlaise Bourdin 320e2739ba6SAlexis Marboeuf PetscCallA(DMGetGlobalVector(pdm, X,ierr)) 321d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmU, U,ierr)) 322d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmA, A,ierr)) 323d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, S,ierr)) 324d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, UA,ierr)) 325d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr)) 326280aadf6SBlaise Bourdin 327d8606c27SBarry Smith PetscCallA(PetscObjectSetName(U, "U",ierr)) 328d8606c27SBarry Smith PetscCallA(PetscObjectSetName(A, "Alpha",ierr)) 329d8606c27SBarry Smith PetscCallA(PetscObjectSetName(S, "Sigma",ierr)) 330d8606c27SBarry Smith PetscCallA(PetscObjectSetName(UA, "UAlpha",ierr)) 331d8606c27SBarry Smith PetscCallA(PetscObjectSetName(UA2, "UAlpha2",ierr)) 332d8606c27SBarry Smith PetscCallA(VecSet(X, -111.0_kPR,ierr)) 333280aadf6SBlaise Bourdin 334280aadf6SBlaise 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 */ 335d8606c27SBarry Smith PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr)) 336d8606c27SBarry Smith PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr)) 337d8606c27SBarry Smith PetscCallA(VecGetArrayF90(UALoc, cval,ierr)) 338d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr)) 339d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr)) 340d8606c27SBarry Smith PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr)) 341280aadf6SBlaise Bourdin 342280aadf6SBlaise Bourdin do p = pStart,pEnd-1 343d8606c27SBarry Smith PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr)) 344280aadf6SBlaise Bourdin if (dofUA > 0) then 345d8606c27SBarry Smith PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr)) 346d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr)) 347280aadf6SBlaise Bourdin closureSize = size(xyz) 348280aadf6SBlaise Bourdin do i = 1,sdim 349280aadf6SBlaise Bourdin do j = 0, closureSize-1,sdim 350280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i) 351280aadf6SBlaise Bourdin end do 352280aadf6SBlaise Bourdin cval(offUA+i) = cval(offUA+i) * sdim / closureSize; 353280aadf6SBlaise Bourdin cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2 354280aadf6SBlaise Bourdin end do 355d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr)) 356280aadf6SBlaise Bourdin end if 357280aadf6SBlaise Bourdin end do 358280aadf6SBlaise Bourdin 359d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr)) 360d8606c27SBarry Smith PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 361d8606c27SBarry Smith PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr)) 362d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr)) 363280aadf6SBlaise Bourdin 364280aadf6SBlaise Bourdin !Update X 365d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr)) 366280aadf6SBlaise Bourdin ! Restrict to U and Alpha 367d8606c27SBarry Smith PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr)) 368d8606c27SBarry Smith PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr)) 369d8606c27SBarry Smith PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, "-ua_vec_view",ierr)) 370d8606c27SBarry Smith PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, "-u_vec_view",ierr)) 371d8606c27SBarry Smith PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, "-a_vec_view",ierr)) 372280aadf6SBlaise Bourdin ! restrict to UA2 373d8606c27SBarry Smith PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr)) 374d8606c27SBarry Smith PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, "-ua2_vec_view",ierr)) 375280aadf6SBlaise Bourdin 376e2739ba6SAlexis Marboeuf ! Getting Natural Vec 377d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr)) 378d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr)) 379280aadf6SBlaise Bourdin 380d8606c27SBarry Smith PetscCallA(VecView(U, viewer,ierr)) 381d8606c27SBarry Smith PetscCallA(VecView(A, viewer,ierr)) 382280aadf6SBlaise Bourdin 383280aadf6SBlaise Bourdin ! Saving U and Alpha in one shot. 384280aadf6SBlaise Bourdin ! For this, we need to cheat and change the Vec's name 385280aadf6SBlaise Bourdin ! Note that in the end we write variables one component at a time, 386280aadf6SBlaise Bourdin ! so that there is no real value in doing this 387d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr)) 388d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr)) 389d8606c27SBarry Smith PetscCallA(VecCopy(UA, tmpVec,ierr)) 390d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "U",ierr)) 391d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 392e2739ba6SAlexis Marboeuf 393280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 394d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 395d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec, viewer,ierr)) 396d8606c27SBarry Smith PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr)) 397d8606c27SBarry Smith PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr)) 398280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 399280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm 400280aadf6SBlaise Bourdin end if 401d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr)) 402280aadf6SBlaise Bourdin 403e2739ba6SAlexis Marboeuf ! same thing with the UA2 Vec obtained from the superDM 404d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr)) 405d8606c27SBarry Smith PetscCallA(VecCopy(UA2, tmpVec,ierr)) 406d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "U",ierr)) 407d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr)) 408d8606c27SBarry Smith PetscCallA(VecView(tmpVec, viewer,ierr)) 409e2739ba6SAlexis Marboeuf 410280aadf6SBlaise Bourdin ! Reading nodal variables in Exodus file 411d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 412d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 413d8606c27SBarry Smith PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr)) 414d8606c27SBarry Smith PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr)) 415280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 416280aadf6SBlaise Bourdin write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm 417280aadf6SBlaise Bourdin end if 418d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr)) 419280aadf6SBlaise Bourdin 420280aadf6SBlaise Bourdin ! Building and saving Sigma 421280aadf6SBlaise Bourdin ! We set sigma_0 = rank (to see partitioning) 422280aadf6SBlaise Bourdin ! sigma_1 = cell set ID 423280aadf6SBlaise Bourdin ! sigma_2 = x_coordinate of the cell center of mass 424d8606c27SBarry Smith PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr)) 425d8606c27SBarry Smith PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr)) 426d8606c27SBarry Smith PetscCallA(DMGetLabelIdIS(dmS, "Cell Sets", csIS,ierr)) 427d8606c27SBarry Smith PetscCallA(DMGetLabelSize(dmS, "Cell Sets",numCS,ierr)) 428d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(csIS, csID,ierr)) 429280aadf6SBlaise Bourdin 430280aadf6SBlaise Bourdin do set = 1, numCS 431d8606c27SBarry Smith PetscCallA(DMGetStratumIS(dmS, "Cell Sets", csID(set), cellIS,ierr)) 432d8606c27SBarry Smith PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr)) 433d8606c27SBarry Smith PetscCallA(ISGetSize(cellIS, numCells,ierr)) 434280aadf6SBlaise Bourdin do cell = 1,numCells 435d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 436d8606c27SBarry Smith PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 437280aadf6SBlaise Bourdin cval(1) = rank 438280aadf6SBlaise Bourdin cval(2) = csID(set) 439280aadf6SBlaise Bourdin cval(3) = 0.0_kPR 440280aadf6SBlaise Bourdin do c = 1, size(xyz),sdim 441280aadf6SBlaise Bourdin cval(3) = cval(3) + xyz(c) 442280aadf6SBlaise Bourdin end do 443280aadf6SBlaise Bourdin cval(3) = cval(3) * sdim / size(xyz) 444d8606c27SBarry Smith PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr)) 445d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr)) 446d8606c27SBarry Smith PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr)) 447280aadf6SBlaise Bourdin end do 448d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr)) 449d8606c27SBarry Smith PetscCallA(ISDestroy(cellIS,ierr)) 450280aadf6SBlaise Bourdin end do 451d8606c27SBarry Smith PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr)) 452d8606c27SBarry Smith PetscCallA(ISDestroy(csIS,ierr)) 453d8606c27SBarry Smith PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, "-s_vec_view",ierr)) 454280aadf6SBlaise Bourdin 455280aadf6SBlaise Bourdin ! Writing zonal variables in Exodus file 456d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr)) 457d8606c27SBarry Smith PetscCallA(VecView(S,viewer,ierr)) 458280aadf6SBlaise Bourdin 459280aadf6SBlaise Bourdin ! Reading zonal variables in Exodus file */ 460d8606c27SBarry Smith PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr)) 461d8606c27SBarry Smith PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr)) 462d8606c27SBarry Smith PetscCallA(PetscObjectSetName(tmpVec, "Sigma",ierr)) 463d8606c27SBarry Smith PetscCallA(VecLoad(tmpVec,viewer,ierr)) 464d8606c27SBarry Smith PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr)) 465d8606c27SBarry Smith PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr)) 466280aadf6SBlaise Bourdin if (norm > PETSC_SQRT_MACHINE_EPSILON) then 467280aadf6SBlaise Bourdin write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm 468280aadf6SBlaise Bourdin end if 469d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr)) 470280aadf6SBlaise Bourdin 471d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr)) 472d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr)) 473d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmS, S,ierr)) 474d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmA, A,ierr)) 475d8606c27SBarry Smith PetscCallA(DMRestoreGlobalVector(dmU, U,ierr)) 476e2739ba6SAlexis Marboeuf PetscCallA(DMRestoreGlobalVector(pdm, X,ierr)) 477d8606c27SBarry Smith PetscCallA(DMDestroy(dmU,ierr)) 478d8606c27SBarry Smith PetscCallA(ISDestroy(isU,ierr)) 479d8606c27SBarry Smith PetscCallA(DMDestroy(dmA,ierr)) 480d8606c27SBarry Smith PetscCallA(ISDestroy(isA,ierr)) 481d8606c27SBarry Smith PetscCallA(DMDestroy(dmS,ierr)) 482d8606c27SBarry Smith PetscCallA(ISDestroy(isS,ierr)) 483d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA,ierr)) 484d8606c27SBarry Smith PetscCallA(ISDestroy(isUA,ierr)) 485d8606c27SBarry Smith PetscCallA(DMDestroy(dmUA2,ierr)) 486e2739ba6SAlexis Marboeuf PetscCallA(DMDestroy(pdm,ierr)) 487e2739ba6SAlexis Marboeuf if (numProc > 1) then 488d8606c27SBarry Smith PetscCallA(DMDestroy(dm,ierr)) 489e2739ba6SAlexis Marboeuf end if 490280aadf6SBlaise Bourdin 491280aadf6SBlaise Bourdin deallocate(pStartDepth) 492280aadf6SBlaise Bourdin deallocate(pEndDepth) 493280aadf6SBlaise Bourdin 494d8606c27SBarry Smith PetscCallA(PetscViewerDestroy(viewer,ierr)) 495d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 4962e8d78feSBlaise Bourdinend program ex62f90 497280aadf6SBlaise Bourdin 498280aadf6SBlaise Bourdin! /*TEST 499280aadf6SBlaise Bourdin! 500280aadf6SBlaise Bourdin! build: 501280aadf6SBlaise Bourdin! requires: exodusii pnetcdf !complex 502280aadf6SBlaise Bourdin! # 2D seq 503e2739ba6SAlexis Marboeuf! test: 504e2739ba6SAlexis Marboeuf! suffix: 0 5052e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1 506e2739ba6SAlexis 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 507e2739ba6SAlexis Marboeuf! test: 508e2739ba6SAlexis Marboeuf! suffix: 1 5092e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1 510e2739ba6SAlexis Marboeuf! 511e2739ba6SAlexis Marboeuf! test: 512e2739ba6SAlexis Marboeuf! suffix: 2 5132e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_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: 3 5172e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2 518e2739ba6SAlexis Marboeuf! test: 519e2739ba6SAlexis Marboeuf! suffix: 4 5202e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2 521e2739ba6SAlexis Marboeuf! test: 522e2739ba6SAlexis Marboeuf! suffix: 5 5232e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2 524e2739ba6SAlexis Marboeuf! # 2D par 525e2739ba6SAlexis Marboeuf! test: 526e2739ba6SAlexis Marboeuf! suffix: 6 527e2739ba6SAlexis Marboeuf! nsize: 2 5282e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1 529e2739ba6SAlexis 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 530e2739ba6SAlexis Marboeuf! test: 531e2739ba6SAlexis Marboeuf! suffix: 7 532e2739ba6SAlexis Marboeuf! nsize: 2 5332e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1 534e2739ba6SAlexis Marboeuf! test: 535e2739ba6SAlexis Marboeuf! suffix: 8 536e2739ba6SAlexis Marboeuf! nsize: 2 5372e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1 538e2739ba6SAlexis Marboeuf! #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 539e2739ba6SAlexis Marboeuf! test: 540e2739ba6SAlexis Marboeuf! suffix: 9 541e2739ba6SAlexis Marboeuf! nsize: 2 5422e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2 543e2739ba6SAlexis Marboeuf! test: 544e2739ba6SAlexis Marboeuf! suffix: 10 545e2739ba6SAlexis Marboeuf! nsize: 2 5462e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2 547e2739ba6SAlexis Marboeuf! test: 548*9d2ffcd0SPierre Jolivet! # Something is now broken with parallel read/write for whatever shape H is 549e2739ba6SAlexis Marboeuf! TODO: broken 550e2739ba6SAlexis Marboeuf! suffix: 11 551e2739ba6SAlexis Marboeuf! nsize: 2 5522e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2 553e2739ba6SAlexis Marboeuf 554e2739ba6SAlexis Marboeuf! #3d seq 555e2739ba6SAlexis Marboeuf! test: 556e2739ba6SAlexis Marboeuf! suffix: 12 5572e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1 558e2739ba6SAlexis Marboeuf! test: 559e2739ba6SAlexis Marboeuf! suffix: 13 5602e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1 561e2739ba6SAlexis 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 562e2739ba6SAlexis Marboeuf! test: 563e2739ba6SAlexis Marboeuf! suffix: 14 5642e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2 565e2739ba6SAlexis Marboeuf! test: 566e2739ba6SAlexis Marboeuf! suffix: 15 5672e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2 568e2739ba6SAlexis 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 569e2739ba6SAlexis Marboeuf! #3d par 570e2739ba6SAlexis Marboeuf! test: 571e2739ba6SAlexis Marboeuf! suffix: 16 572e2739ba6SAlexis Marboeuf! nsize: 2 5732e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1 574e2739ba6SAlexis Marboeuf! test: 575e2739ba6SAlexis Marboeuf! suffix: 17 576e2739ba6SAlexis Marboeuf! nsize: 2 5772e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1 578e2739ba6SAlexis 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 579e2739ba6SAlexis Marboeuf! test: 580e2739ba6SAlexis Marboeuf! suffix: 18 581e2739ba6SAlexis Marboeuf! nsize: 2 5822e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2 583e2739ba6SAlexis Marboeuf! test: 584e2739ba6SAlexis Marboeuf! suffix: 19 585e2739ba6SAlexis Marboeuf! nsize: 2 5862e8d78feSBlaise Bourdin! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2 587e2739ba6SAlexis 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 588e2739ba6SAlexis Marboeuf 589280aadf6SBlaise Bourdin! TEST*/ 590