xref: /petsc/src/dm/impls/plex/tests/ex62f90.F90 (revision dcb3e68992f1c4897946af7e8406e2b4165e50f2)
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))
88*dcb3e689SBarry Smith    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-i',ifilename,flg,ierr))
89*dcb3e689SBarry Smith    PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing input file name -i <input file name>')
90*dcb3e689SBarry Smith    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-o',ofilename,flg,ierr))
91*dcb3e689SBarry Smith    PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing output file name -o <output file name>')
92*dcb3e689SBarry Smith    PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-order',order,flg,ierr))
93280aadf6SBlaise Bourdin    if ((order > 2) .or. (order < 1)) then
94280aadf6SBlaise Bourdin        write(IOBuffer,'("Unsupported polynomial order ", I2, " not in [1,2]")') order
95280aadf6SBlaise Bourdin        SETERRQ(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))
103*dcb3e689SBarry Smith    PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,'-dm_view',ierr))
104280aadf6SBlaise Bourdin
105280aadf6SBlaise Bourdin    ! Create the exodus result file
106280aadf6SBlaise Bourdin
107a5b23f4aSJose E. Roman    ! enable exodus debugging information
108d8606c27SBarry Smith    PetscCallA(exopts(EXVRBS+EXDEBG,ierr))
109280aadf6SBlaise Bourdin    ! Create the exodus file
110d8606c27SBarry Smith    PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr))
111280aadf6SBlaise Bourdin    ! The long way would be
112280aadf6SBlaise Bourdin    !
113d8606c27SBarry Smith    ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr))
114d8606c27SBarry Smith    ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr))
115d8606c27SBarry Smith    ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr))
116d8606c27SBarry Smith    ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr))
117280aadf6SBlaise Bourdin
118280aadf6SBlaise Bourdin    ! set the mesh order
119d8606c27SBarry Smith    PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr))
120d8606c27SBarry Smith    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
121280aadf6SBlaise Bourdin    !
122280aadf6SBlaise Bourdin    !    Notice how the exodus file is actually NOT open at this point (exoid is -1)
123d5b43468SJose E. Roman    !    Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to
124280aadf6SBlaise Bourdin    !    write the geometry (the DM), which can only be done on a brand new file.
125280aadf6SBlaise Bourdin    !
126280aadf6SBlaise Bourdin
127280aadf6SBlaise Bourdin    ! Save the geometry to the file, erasing all previous content
128d8606c27SBarry Smith    PetscCallA(DMView(dm,viewer,ierr))
129d8606c27SBarry Smith    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
130280aadf6SBlaise Bourdin    !
131280aadf6SBlaise Bourdin    !    Note how the exodus file is now open
132280aadf6SBlaise Bourdin    !
133*dcb3e689SBarry Smith    ! 'Format' the exodus result file, i.e. allocate space for nodal and zonal variables
134280aadf6SBlaise Bourdin    select case(sdim)
135280aadf6SBlaise Bourdin    case(2)
136280aadf6SBlaise Bourdin        numNodalVar = 3
137*dcb3e689SBarry Smith        nodalVarName(1:numNodalVar) = ['U_x  ','U_y  ','Alpha']
138280aadf6SBlaise Bourdin        numZonalVar = 3
139*dcb3e689SBarry Smith        zonalVarName(1:numZonalVar) = ['Sigma_11','Sigma_22','Sigma_12']
140280aadf6SBlaise Bourdin    case(3)
141280aadf6SBlaise Bourdin        numNodalVar = 4
142*dcb3e689SBarry Smith        nodalVarName(1:numNodalVar) = ['U_x  ','U_y  ','U_z  ','Alpha']
143280aadf6SBlaise Bourdin        numZonalVar = 6
144*dcb3e689SBarry Smith        zonalVarName(1:numZonalVar) = ['Sigma_11','Sigma_22','Sigma_33','Sigma_23','Sigma_13','Sigma_12']
145280aadf6SBlaise Bourdin    case default
146280aadf6SBlaise Bourdin        write(IOBuffer,'("No layout for dimension ",I2)') sdim
147280aadf6SBlaise Bourdin    end select
148d8606c27SBarry Smith    PetscCallA(PetscViewerExodusIIGetId(viewer,exoid,ierr))
149*dcb3e689SBarry Smith    PetscCallA(expvp(exoid, 'E', numZonalVar,ierr))
150*dcb3e689SBarry Smith    PetscCallA(expvan(exoid, 'E', numZonalVar, zonalVarName,ierr))
151*dcb3e689SBarry Smith    PetscCallA(expvp(exoid, 'N', numNodalVar,ierr))
152*dcb3e689SBarry Smith    PetscCallA(expvan(exoid, 'N', numNodalVar, nodalVarName,ierr))
153d8606c27SBarry Smith    PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK,numCS,PETSC_NULL_REAL,sjunk,ierr))
154280aadf6SBlaise Bourdin
155280aadf6SBlaise Bourdin    !    An exodusII truth table specifies which fields are saved at which time step
156280aadf6SBlaise Bourdin    !    It speeds up I/O but reserving space for fields in the file ahead of time.
157280aadf6SBlaise Bourdin    allocate(truthtable(numCS,numZonalVar))
158280aadf6SBlaise Bourdin    truthtable = .true.
159d8606c27SBarry Smith    PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr))
160280aadf6SBlaise Bourdin    deallocate(truthtable)
161280aadf6SBlaise Bourdin
162280aadf6SBlaise Bourdin    !   Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
163280aadf6SBlaise Bourdin    do step = 1,numstep
164d8606c27SBarry Smith        PetscCallA(exptim(exoid,step,Real(step,kind=kPR),ierr))
165280aadf6SBlaise Bourdin    end do
166280aadf6SBlaise Bourdin
167d8606c27SBarry Smith    PetscCallA(DMSetUseNatural(dm, PETSC_TRUE, ierr))
168d8606c27SBarry Smith    PetscCallA(DMPlexGetPartitioner(dm,part,ierr))
169d8606c27SBarry Smith    PetscCallA(PetscPartitionerSetFromOptions(part,ierr))
170d8606c27SBarry Smith    PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr))
171280aadf6SBlaise Bourdin
172280aadf6SBlaise Bourdin    if (numProc > 1) then
173d8606c27SBarry Smith        PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr))
174d8606c27SBarry Smith        PetscCallA(PetscSFDestroy(migrationSF,ierr))
175e2739ba6SAlexis Marboeuf    else
176e2739ba6SAlexis Marboeuf        pdm = dm
177280aadf6SBlaise Bourdin    end if
178*dcb3e689SBarry Smith    PetscCallA(DMViewFromOptions(pdm,PETSC_NULL_OPTIONS,'-dm_view',ierr))
179280aadf6SBlaise Bourdin
180e2739ba6SAlexis Marboeuf    PetscCallA(PetscObjectGetComm(pdm,comm,ierr))
181d8606c27SBarry Smith    PetscCallA(PetscSectionCreate(comm, section,ierr))
182d8606c27SBarry Smith    PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr))
183*dcb3e689SBarry Smith    PetscCallA(PetscSectionSetFieldName(section, fieldU, 'U',ierr))
184*dcb3e689SBarry Smith    PetscCallA(PetscSectionSetFieldName(section, fieldA, 'Alpha',ierr))
185*dcb3e689SBarry Smith    PetscCallA(PetscSectionSetFieldName(section, fieldS, 'Sigma',ierr))
186e2739ba6SAlexis Marboeuf    PetscCallA(DMPlexGetChart(pdm, pStart, pEnd,ierr))
187d8606c27SBarry Smith    PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr))
188280aadf6SBlaise Bourdin
189280aadf6SBlaise Bourdin    allocate(pStartDepth(sdim+1))
190280aadf6SBlaise Bourdin    allocate(pEndDepth(sdim+1))
191280aadf6SBlaise Bourdin    do d = 1, sdim+1
192e2739ba6SAlexis Marboeuf        PetscCallA(DMPlexGetDepthStratum(pdm, d-1, pStartDepth(d), pEndDepth(d),ierr))
193280aadf6SBlaise Bourdin    end do
194280aadf6SBlaise Bourdin
195280aadf6SBlaise Bourdin    ! Vector field U, Scalar field Alpha, Tensor field Sigma
196e2739ba6SAlexis Marboeuf    PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr));
197e2739ba6SAlexis Marboeuf    PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr));
198e2739ba6SAlexis Marboeuf    PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr));
199280aadf6SBlaise Bourdin
200280aadf6SBlaise Bourdin    ! Going through cell sets then cells, and setting up storage for the sections
201*dcb3e689SBarry Smith    PetscCallA(DMGetLabelSize(pdm, 'Cell Sets', numCS, ierr))
202*dcb3e689SBarry Smith    PetscCallA(DMGetLabelIdIS(pdm, 'Cell Sets', csIS, ierr))
203d8606c27SBarry Smith    PetscCallA(ISGetIndicesF90(csIS, csID, ierr))
204280aadf6SBlaise Bourdin    do set = 1,numCS
205*dcb3e689SBarry Smith        PetscCallA(DMGetStratumSize(pdm, 'Cell Sets', csID(set), numCells,ierr))
206*dcb3e689SBarry Smith        PetscCallA(DMGetStratumIS(pdm, 'Cell Sets', csID(set), cellIS,ierr))
207280aadf6SBlaise Bourdin        if (numCells > 0) then
208280aadf6SBlaise Bourdin            select case(sdim)
209280aadf6SBlaise Bourdin            case(2)
210280aadf6SBlaise Bourdin                dofs => dofS2D
211280aadf6SBlaise Bourdin            case(3)
212280aadf6SBlaise Bourdin                dofs => dofS3D
213280aadf6SBlaise Bourdin            case default
214280aadf6SBlaise Bourdin                write(IOBuffer,'("No layout for dimension ",I2)') sdim
215280aadf6SBlaise Bourdin                SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer)
216280aadf6SBlaise Bourdin            end select ! sdim
217280aadf6SBlaise Bourdin
218280aadf6SBlaise Bourdin            ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
219280aadf6SBlaise Bourdin            ! It will not be enough to identify more exotic elements like pyramid or prisms...  */
220d8606c27SBarry Smith            PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr))
221280aadf6SBlaise Bourdin            nullify(closureA)
222e2739ba6SAlexis Marboeuf            PetscCallA(DMPlexGetTransitiveClosure(pdm,cellID(1), PETSC_TRUE, closureA,ierr))
223280aadf6SBlaise Bourdin            select case(size(closureA)/2)
224280aadf6SBlaise Bourdin            case(7) ! Tri
225280aadf6SBlaise Bourdin                if (order == 1) then
226280aadf6SBlaise Bourdin                    dofU => dofUP1Tri
227280aadf6SBlaise Bourdin                    dofA => dofAP1Tri
228280aadf6SBlaise Bourdin                else
229280aadf6SBlaise Bourdin                    dofU => dofUP2Tri
230280aadf6SBlaise Bourdin                    dofA => dofAP2Tri
231280aadf6SBlaise Bourdin                end if
232280aadf6SBlaise Bourdin            case(9) ! Quad
233280aadf6SBlaise Bourdin                if (order == 1) then
234280aadf6SBlaise Bourdin                    dofU => dofUP1Quad
235280aadf6SBlaise Bourdin                    dofA => dofAP1Quad
236280aadf6SBlaise Bourdin                else
237280aadf6SBlaise Bourdin                    dofU => dofUP2Quad
238280aadf6SBlaise Bourdin                    dofA => dofAP2Quad
239280aadf6SBlaise Bourdin                end if
240280aadf6SBlaise Bourdin            case(15) ! Tet
241280aadf6SBlaise Bourdin                if (order == 1) then
242280aadf6SBlaise Bourdin                    dofU => dofUP1Tet
243280aadf6SBlaise Bourdin                    dofA => dofAP1Tet
244280aadf6SBlaise Bourdin                else
245280aadf6SBlaise Bourdin                    dofU => dofUP2Tet
246280aadf6SBlaise Bourdin                    dofA => dofAP2Tet
247280aadf6SBlaise Bourdin                end if
248280aadf6SBlaise Bourdin            case(27) ! Hex
249280aadf6SBlaise Bourdin                if (order == 1) then
250280aadf6SBlaise Bourdin                    dofU => dofUP1Hex
251280aadf6SBlaise Bourdin                    dofA => dofAP1Hex
252280aadf6SBlaise Bourdin                else
253280aadf6SBlaise Bourdin                    dofU => dofUP2Hex
254280aadf6SBlaise Bourdin                    dofA => dofAP2Hex
255280aadf6SBlaise Bourdin                end if
256280aadf6SBlaise Bourdin            case default
257280aadf6SBlaise Bourdin                write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2
258280aadf6SBlaise Bourdin                SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer)
259280aadf6SBlaise Bourdin            end select
260e2739ba6SAlexis Marboeuf            PetscCallA(DMPlexRestoreTransitiveClosure(pdm, cellID(1), PETSC_TRUE,closureA,ierr))
261280aadf6SBlaise Bourdin            do cell = 1,numCells!
262280aadf6SBlaise Bourdin                nullify(closure)
263e2739ba6SAlexis Marboeuf                PetscCallA(DMPlexGetTransitiveClosure(pdm, cellID(cell), PETSC_TRUE, closure,ierr))
264280aadf6SBlaise Bourdin                do p = 1,size(closure),2
265280aadf6SBlaise Bourdin                    ! find the depth of p
266280aadf6SBlaise Bourdin                    do d = 1,sdim+1
267280aadf6SBlaise Bourdin                        if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then
268d8606c27SBarry Smith                            PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr))
269d8606c27SBarry Smith                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr))
270d8606c27SBarry Smith                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr))
271d8606c27SBarry Smith                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr))
272280aadf6SBlaise Bourdin                        end if ! closure(p)
273280aadf6SBlaise Bourdin                    end do ! d
274280aadf6SBlaise Bourdin                end do ! p
275e2739ba6SAlexis Marboeuf                PetscCallA(DMPlexRestoreTransitiveClosure(pdm, cellID(cell), PETSC_TRUE, closure,ierr))
276280aadf6SBlaise Bourdin            end do ! cell
277d8606c27SBarry Smith            PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr))
278d8606c27SBarry Smith            PetscCallA(ISDestroy(cellIS,ierr))
279280aadf6SBlaise Bourdin        end if ! numCells
280280aadf6SBlaise Bourdin    end do ! set
281d8606c27SBarry Smith    PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr))
282d8606c27SBarry Smith    PetscCallA(ISDestroy(csIS,ierr))
283d8606c27SBarry Smith    PetscCallA(PetscSectionSetUp(section,ierr))
284e2739ba6SAlexis Marboeuf    PetscCallA(DMSetLocalSection(pdm, section,ierr))
285*dcb3e689SBarry Smith    PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, '-dm_section_view',ierr))
286d8606c27SBarry Smith    PetscCallA(PetscSectionDestroy(section,ierr))
287280aadf6SBlaise Bourdin
288e2739ba6SAlexis Marboeuf    ! Creating section on the sequential DM + creating the GlobalToNatural SF
289e2739ba6SAlexis Marboeuf    if (numProc > 1) then
290e2739ba6SAlexis Marboeuf        PetscCallA(DMPlexGetMigrationSF(pdm, natPointSF, ierr))
291e2739ba6SAlexis Marboeuf        PetscCallA(DMGetLocalSection(pdm, rootSection, ierr))
292e2739ba6SAlexis Marboeuf        PetscCallA(PetscSFCreateInverseSF(natPointSF, natPointSFInv, ierr))
293e2739ba6SAlexis Marboeuf        PetscCallA(PetscSectionCreate(PETSC_COMM_WORLD, leafSection, ierr))
294e2739ba6SAlexis Marboeuf        PetscCallA(PetscSFDistributeSection(natPointSFInv, rootSection, PETSC_NULL_INTEGER, leafSection, ierr))
295e2739ba6SAlexis Marboeuf        PetscCallA(DMSetLocalSection(dm, leafSection, ierr))
296e2739ba6SAlexis Marboeuf        PetscCallA(DMPlexCreateGlobalToNaturalSF(pdm, leafSection, natPointSF, natSF, ierr))
297e2739ba6SAlexis Marboeuf        PetscCallA(PetscSFDestroy(natPointSFInv, ierr))
298e2739ba6SAlexis Marboeuf        PetscCallA(PetscSectionDestroy(leafSection,ierr))
299e2739ba6SAlexis Marboeuf        PetscCallA(DMSetNaturalSF(pdm, natSF, ierr))
300e2739ba6SAlexis Marboeuf        PetscCallA(PetscObjectDereference(natSF, ierr))
301e2739ba6SAlexis Marboeuf    end if
302e2739ba6SAlexis Marboeuf
303280aadf6SBlaise Bourdin    ! Get DM and IS for each field of dm
304e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldU,  isU,  dmU,ierr))
305e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldA,  isA,  dmA,ierr))
306e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldS,  isS,  dmS,ierr))
307e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA,ierr))
308280aadf6SBlaise Bourdin
309280aadf6SBlaise Bourdin    !Create the exodus result file
310280aadf6SBlaise Bourdin    allocate(dmList(2))
311280aadf6SBlaise Bourdin    dmList(1) = dmU;
312280aadf6SBlaise Bourdin    dmList(2) = dmA;
313d8606c27SBarry Smith    PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr))
314280aadf6SBlaise Bourdin    deallocate(dmList)
315280aadf6SBlaise Bourdin
316e2739ba6SAlexis Marboeuf    PetscCallA(DMGetGlobalVector(pdm,   X,ierr))
317d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmU,  U,ierr))
318d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmA,  A,ierr))
319d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmS,  S,ierr))
320d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA, UA,ierr))
321d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr))
322280aadf6SBlaise Bourdin
323*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(U,  'U',ierr))
324*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(A,  'Alpha',ierr))
325*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(S,  'Sigma',ierr))
326*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(UA, 'UAlpha',ierr))
327*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(UA2, 'UAlpha2',ierr))
328d8606c27SBarry Smith    PetscCallA(VecSet(X, -111.0_kPR,ierr))
329280aadf6SBlaise Bourdin
330280aadf6SBlaise 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 */
331d8606c27SBarry Smith    PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr))
332d8606c27SBarry Smith    PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr))
333d8606c27SBarry Smith    PetscCallA(VecGetArrayF90(UALoc, cval,ierr))
334d8606c27SBarry Smith    PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr))
335d8606c27SBarry Smith    PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr))
336d8606c27SBarry Smith    PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr))
337280aadf6SBlaise Bourdin
338280aadf6SBlaise Bourdin    do p = pStart,pEnd-1
339d8606c27SBarry Smith        PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr))
340280aadf6SBlaise Bourdin        if (dofUA > 0) then
341d8606c27SBarry Smith            PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr))
342d8606c27SBarry Smith            PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr))
343280aadf6SBlaise Bourdin            closureSize = size(xyz)
344280aadf6SBlaise Bourdin            do i = 1,sdim
345280aadf6SBlaise Bourdin                do j = 0, closureSize-1,sdim
346280aadf6SBlaise Bourdin                    cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i)
347280aadf6SBlaise Bourdin                end do
348280aadf6SBlaise Bourdin                cval(offUA+i) = cval(offUA+i) * sdim / closureSize;
349280aadf6SBlaise Bourdin                cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2
350280aadf6SBlaise Bourdin            end do
351d8606c27SBarry Smith            PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr))
352280aadf6SBlaise Bourdin        end if
353280aadf6SBlaise Bourdin    end do
354280aadf6SBlaise Bourdin
355d8606c27SBarry Smith    PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr))
356d8606c27SBarry Smith    PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr))
357d8606c27SBarry Smith    PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr))
358d8606c27SBarry Smith    PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr))
359280aadf6SBlaise Bourdin
360280aadf6SBlaise Bourdin    !Update X
361d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr))
362280aadf6SBlaise Bourdin    ! Restrict to U and Alpha
363d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr))
364d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr))
365*dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, '-ua_vec_view',ierr))
366*dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, '-u_vec_view',ierr))
367*dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, '-a_vec_view',ierr))
368280aadf6SBlaise Bourdin    ! restrict to UA2
369d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr))
370*dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, '-ua2_vec_view',ierr))
371280aadf6SBlaise Bourdin
372e2739ba6SAlexis Marboeuf    ! Getting Natural Vec
373d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr))
374d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr))
375280aadf6SBlaise Bourdin
376d8606c27SBarry Smith    PetscCallA(VecView(U, viewer,ierr))
377d8606c27SBarry Smith    PetscCallA(VecView(A, viewer,ierr))
378280aadf6SBlaise Bourdin
379280aadf6SBlaise Bourdin    ! Saving U and Alpha in one shot.
380280aadf6SBlaise Bourdin    ! For this, we need to cheat and change the Vec's name
381280aadf6SBlaise Bourdin    ! Note that in the end we write variables one component at a time,
382280aadf6SBlaise Bourdin    ! so that there is no real value in doing this
383d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr))
384d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr))
385d8606c27SBarry Smith    PetscCallA(VecCopy(UA, tmpVec,ierr))
386*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr))
387d8606c27SBarry Smith    PetscCallA(VecView(tmpVec, viewer,ierr))
388e2739ba6SAlexis Marboeuf
389280aadf6SBlaise Bourdin    ! Reading nodal variables in Exodus file
390d8606c27SBarry Smith    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
391d8606c27SBarry Smith    PetscCallA(VecLoad(tmpVec, viewer,ierr))
392d8606c27SBarry Smith    PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr))
393d8606c27SBarry Smith    PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr))
394280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
395280aadf6SBlaise Bourdin        write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm
396280aadf6SBlaise Bourdin    end if
397d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr))
398280aadf6SBlaise Bourdin
399e2739ba6SAlexis Marboeuf    ! same thing with the UA2 Vec obtained from the superDM
400d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr))
401d8606c27SBarry Smith    PetscCallA(VecCopy(UA2, tmpVec,ierr))
402*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr))
403d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr))
404d8606c27SBarry Smith    PetscCallA(VecView(tmpVec, viewer,ierr))
405e2739ba6SAlexis Marboeuf
406280aadf6SBlaise Bourdin    ! Reading nodal variables in Exodus file
407d8606c27SBarry Smith    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
408d8606c27SBarry Smith    PetscCallA(VecLoad(tmpVec,viewer,ierr))
409d8606c27SBarry Smith    PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr))
410d8606c27SBarry Smith    PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr))
411280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
412280aadf6SBlaise Bourdin        write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm
413280aadf6SBlaise Bourdin    end if
414d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr))
415280aadf6SBlaise Bourdin
416280aadf6SBlaise Bourdin    ! Building and saving Sigma
417280aadf6SBlaise Bourdin    !   We set sigma_0 = rank (to see partitioning)
418280aadf6SBlaise Bourdin    !          sigma_1 = cell set ID
419280aadf6SBlaise Bourdin    !          sigma_2 = x_coordinate of the cell center of mass
420d8606c27SBarry Smith    PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr))
421d8606c27SBarry Smith    PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr))
422*dcb3e689SBarry Smith    PetscCallA(DMGetLabelIdIS(dmS, 'Cell Sets', csIS,ierr))
423*dcb3e689SBarry Smith    PetscCallA(DMGetLabelSize(dmS, 'Cell Sets',numCS,ierr))
424d8606c27SBarry Smith    PetscCallA(ISGetIndicesF90(csIS, csID,ierr))
425280aadf6SBlaise Bourdin
426280aadf6SBlaise Bourdin    do set = 1, numCS
427*dcb3e689SBarry Smith        PetscCallA(DMGetStratumIS(dmS, 'Cell Sets', csID(set), cellIS,ierr))
428d8606c27SBarry Smith        PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr))
429d8606c27SBarry Smith        PetscCallA(ISGetSize(cellIS, numCells,ierr))
430280aadf6SBlaise Bourdin        do cell = 1,numCells
431d8606c27SBarry Smith            PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr))
432d8606c27SBarry Smith            PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr))
433280aadf6SBlaise Bourdin            cval(1) = rank
434280aadf6SBlaise Bourdin            cval(2) = csID(set)
435280aadf6SBlaise Bourdin            cval(3) = 0.0_kPR
436280aadf6SBlaise Bourdin            do c = 1, size(xyz),sdim
437280aadf6SBlaise Bourdin                cval(3) = cval(3) + xyz(c)
438280aadf6SBlaise Bourdin            end do
439280aadf6SBlaise Bourdin            cval(3) = cval(3) * sdim / size(xyz)
440d8606c27SBarry Smith            PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr))
441d8606c27SBarry Smith            PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr))
442d8606c27SBarry Smith            PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr))
443280aadf6SBlaise Bourdin        end do
444d8606c27SBarry Smith        PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr))
445d8606c27SBarry Smith        PetscCallA(ISDestroy(cellIS,ierr))
446280aadf6SBlaise Bourdin    end do
447d8606c27SBarry Smith    PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr))
448d8606c27SBarry Smith    PetscCallA(ISDestroy(csIS,ierr))
449*dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, '-s_vec_view',ierr))
450280aadf6SBlaise Bourdin
451280aadf6SBlaise Bourdin    ! Writing zonal variables in Exodus file
452d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr))
453d8606c27SBarry Smith    PetscCallA(VecView(S,viewer,ierr))
454280aadf6SBlaise Bourdin
455280aadf6SBlaise Bourdin    ! Reading zonal variables in Exodus file */
456d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr))
457d8606c27SBarry Smith    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
458*dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(tmpVec, 'Sigma',ierr))
459d8606c27SBarry Smith    PetscCallA(VecLoad(tmpVec,viewer,ierr))
460d8606c27SBarry Smith    PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr))
461d8606c27SBarry Smith    PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr))
462280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
463280aadf6SBlaise Bourdin       write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm
464280aadf6SBlaise Bourdin    end if
465d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr))
466280aadf6SBlaise Bourdin
467d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr))
468d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr))
469d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmS,  S,ierr))
470d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmA,  A,ierr))
471d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmU,  U,ierr))
472e2739ba6SAlexis Marboeuf    PetscCallA(DMRestoreGlobalVector(pdm,   X,ierr))
473d8606c27SBarry Smith    PetscCallA(DMDestroy(dmU,ierr))
474d8606c27SBarry Smith    PetscCallA(ISDestroy(isU,ierr))
475d8606c27SBarry Smith    PetscCallA(DMDestroy(dmA,ierr))
476d8606c27SBarry Smith    PetscCallA(ISDestroy(isA,ierr))
477d8606c27SBarry Smith    PetscCallA(DMDestroy(dmS,ierr))
478d8606c27SBarry Smith    PetscCallA(ISDestroy(isS,ierr))
479d8606c27SBarry Smith    PetscCallA(DMDestroy(dmUA,ierr))
480d8606c27SBarry Smith    PetscCallA(ISDestroy(isUA,ierr))
481d8606c27SBarry Smith    PetscCallA(DMDestroy(dmUA2,ierr))
482e2739ba6SAlexis Marboeuf    PetscCallA(DMDestroy(pdm,ierr))
483e2739ba6SAlexis Marboeuf    if (numProc > 1) then
484d8606c27SBarry Smith        PetscCallA(DMDestroy(dm,ierr))
485e2739ba6SAlexis Marboeuf    end if
486280aadf6SBlaise Bourdin
487280aadf6SBlaise Bourdin    deallocate(pStartDepth)
488280aadf6SBlaise Bourdin    deallocate(pEndDepth)
489280aadf6SBlaise Bourdin
490d8606c27SBarry Smith    PetscCallA(PetscViewerDestroy(viewer,ierr))
491d8606c27SBarry Smith    PetscCallA(PetscFinalize(ierr))
4922e8d78feSBlaise Bourdinend program ex62f90
493280aadf6SBlaise Bourdin
494280aadf6SBlaise Bourdin! /*TEST
495280aadf6SBlaise Bourdin!
496280aadf6SBlaise Bourdin! build:
497280aadf6SBlaise Bourdin!   requires: exodusii pnetcdf !complex
498280aadf6SBlaise Bourdin!   # 2D seq
499e2739ba6SAlexis Marboeuf! test:
500e2739ba6SAlexis Marboeuf!   suffix: 0
5012e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1
502e2739ba6SAlexis 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
503e2739ba6SAlexis Marboeuf! test:
504e2739ba6SAlexis Marboeuf!   suffix: 1
5052e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1
506e2739ba6SAlexis Marboeuf!
507e2739ba6SAlexis Marboeuf! test:
508e2739ba6SAlexis Marboeuf!   suffix: 2
5092e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1
510e2739ba6SAlexis 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
511e2739ba6SAlexis Marboeuf! test:
512e2739ba6SAlexis Marboeuf!   suffix: 3
5132e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2
514e2739ba6SAlexis Marboeuf! test:
515e2739ba6SAlexis Marboeuf!   suffix: 4
5162e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2
517e2739ba6SAlexis Marboeuf! test:
518e2739ba6SAlexis Marboeuf!   suffix: 5
5192e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2
520e2739ba6SAlexis Marboeuf!   # 2D par
521e2739ba6SAlexis Marboeuf! test:
522e2739ba6SAlexis Marboeuf!   suffix: 6
523e2739ba6SAlexis Marboeuf!   nsize: 2
5242e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1
525e2739ba6SAlexis 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
526e2739ba6SAlexis Marboeuf! test:
527e2739ba6SAlexis Marboeuf!   suffix: 7
528e2739ba6SAlexis Marboeuf!   nsize: 2
5292e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1
530e2739ba6SAlexis Marboeuf! test:
531e2739ba6SAlexis Marboeuf!   suffix: 8
532e2739ba6SAlexis Marboeuf!   nsize: 2
5332e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1
534e2739ba6SAlexis Marboeuf!   #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
535e2739ba6SAlexis Marboeuf! test:
536e2739ba6SAlexis Marboeuf!   suffix: 9
537e2739ba6SAlexis Marboeuf!   nsize: 2
5382e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2
539e2739ba6SAlexis Marboeuf! test:
540e2739ba6SAlexis Marboeuf!   suffix: 10
541e2739ba6SAlexis Marboeuf!   nsize: 2
5422e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2
543e2739ba6SAlexis Marboeuf! test:
544e2739ba6SAlexis Marboeuf!   # Something is now broken with parallel read/write for wahtever shape H is
545e2739ba6SAlexis Marboeuf!   TODO: broken
546e2739ba6SAlexis Marboeuf!   suffix: 11
547e2739ba6SAlexis Marboeuf!   nsize: 2
5482e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2
549e2739ba6SAlexis Marboeuf
550e2739ba6SAlexis Marboeuf!   #3d seq
551e2739ba6SAlexis Marboeuf! test:
552e2739ba6SAlexis Marboeuf!   suffix: 12
5532e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1
554e2739ba6SAlexis Marboeuf! test:
555e2739ba6SAlexis Marboeuf!   suffix: 13
5562e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1
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! test:
559e2739ba6SAlexis Marboeuf!   suffix: 14
5602e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2
561e2739ba6SAlexis Marboeuf! test:
562e2739ba6SAlexis Marboeuf!   suffix: 15
5632e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2
564e2739ba6SAlexis 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
565e2739ba6SAlexis Marboeuf!   #3d par
566e2739ba6SAlexis Marboeuf! test:
567e2739ba6SAlexis Marboeuf!   suffix: 16
568e2739ba6SAlexis Marboeuf!   nsize: 2
5692e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1
570e2739ba6SAlexis Marboeuf! test:
571e2739ba6SAlexis Marboeuf!   suffix: 17
572e2739ba6SAlexis Marboeuf!   nsize: 2
5732e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1
574e2739ba6SAlexis Marboeuf!   #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
575e2739ba6SAlexis Marboeuf! test:
576e2739ba6SAlexis Marboeuf!   suffix: 18
577e2739ba6SAlexis Marboeuf!   nsize: 2
5782e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2
579e2739ba6SAlexis Marboeuf! test:
580e2739ba6SAlexis Marboeuf!   suffix: 19
581e2739ba6SAlexis Marboeuf!   nsize: 2
5822e8d78feSBlaise Bourdin!   args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2
583e2739ba6SAlexis 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
584e2739ba6SAlexis Marboeuf
585280aadf6SBlaise Bourdin! TEST*/
586