xref: /petsc/src/dm/impls/plex/tests/ex26f90.F90 (revision f51a5268c6b8e4f2d40024d7e7d67866ed2bae1c)
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
7dd01b7e5SBarry Smith    ! 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
75*f51a5268SBarry Smith    PetscLayout                         :: layout
76280aadf6SBlaise Bourdin
77280aadf6SBlaise Bourdin    type(tVec)                          :: tmpVec
78280aadf6SBlaise Bourdin    PetscReal                           :: norm
79280aadf6SBlaise Bourdin    PetscReal                           :: time = 1.234_kPR
80280aadf6SBlaise Bourdin
81e2739ba6SAlexis Marboeuf    PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER,ierr))
82e2739ba6SAlexis Marboeuf    if (ierr /= 0) then
83e2739ba6SAlexis Marboeuf      print*,'Unable to initialize PETSc'
84e2739ba6SAlexis Marboeuf      stop
85e2739ba6SAlexis Marboeuf    endif
86280aadf6SBlaise Bourdin
87d8606c27SBarry Smith    PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
88d8606c27SBarry Smith    PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,numProc,ierr))
89dcb3e689SBarry Smith    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-i',ifilename,flg,ierr))
90dcb3e689SBarry Smith    PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing input file name -i <input file name>')
91dcb3e689SBarry Smith    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-o',ofilename,flg,ierr))
92dcb3e689SBarry Smith    PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing output file name -o <output file name>')
93dcb3e689SBarry Smith    PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-order',order,flg,ierr))
94280aadf6SBlaise Bourdin    if ((order > 2) .or. (order < 1)) then
95280aadf6SBlaise Bourdin        write(IOBuffer,'("Unsupported polynomial order ", I2, " not in [1,2]")') order
96de401611SBlaise Bourdin        SETERRA(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,IOBuffer)
97280aadf6SBlaise Bourdin    end if
98280aadf6SBlaise Bourdin
99280aadf6SBlaise Bourdin    ! Read the mesh in any supported format
100d8606c27SBarry Smith    PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr))
101e2739ba6SAlexis Marboeuf    PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr));
102e2739ba6SAlexis Marboeuf    PetscCallA(DMSetFromOptions(dm,ierr));
103d8606c27SBarry Smith    PetscCallA(DMGetDimension(dm, sdim,ierr))
104dcb3e689SBarry Smith    PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,'-dm_view',ierr));
105280aadf6SBlaise Bourdin
106280aadf6SBlaise Bourdin    ! Create the exodus result file
107280aadf6SBlaise Bourdin
108a5b23f4aSJose E. Roman    ! enable exodus debugging information
109d8606c27SBarry Smith    PetscCallA(exopts(EXVRBS+EXDEBG,ierr))
110280aadf6SBlaise Bourdin    ! Create the exodus file
111d8606c27SBarry Smith    PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr))
112280aadf6SBlaise Bourdin    ! The long way would be
113280aadf6SBlaise Bourdin    !
114d8606c27SBarry Smith    ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr))
115d8606c27SBarry Smith    ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr))
116d8606c27SBarry Smith    ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr))
117d8606c27SBarry Smith    ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr))
118280aadf6SBlaise Bourdin
119280aadf6SBlaise Bourdin    ! set the mesh order
120d8606c27SBarry Smith    PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr))
121d8606c27SBarry Smith    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
122280aadf6SBlaise Bourdin    !
123280aadf6SBlaise Bourdin    !    Notice how the exodus file is actually NOT open at this point (exoid is -1)
124d5b43468SJose E. Roman    !    Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to
125280aadf6SBlaise Bourdin    !    write the geometry (the DM), which can only be done on a brand new file.
126280aadf6SBlaise Bourdin    !
127280aadf6SBlaise Bourdin
128280aadf6SBlaise Bourdin    ! Save the geometry to the file, erasing all previous content
129d8606c27SBarry Smith    PetscCallA(DMView(dm,viewer,ierr))
130d8606c27SBarry Smith    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
131280aadf6SBlaise Bourdin    !
132280aadf6SBlaise Bourdin    !    Note how the exodus file is now open
133280aadf6SBlaise Bourdin    !
134dcb3e689SBarry Smith    ! 'Format' the exodus result file, i.e. allocate space for nodal and zonal variables
135280aadf6SBlaise Bourdin    select case(sdim)
136280aadf6SBlaise Bourdin    case(2)
137280aadf6SBlaise Bourdin        numNodalVar = 3
138dcb3e689SBarry Smith        nodalVarName(1:numNodalVar) = ['U_x  ','U_y  ','Alpha']
139280aadf6SBlaise Bourdin        numZonalVar = 3
140dcb3e689SBarry Smith        zonalVarName(1:numZonalVar) = ['Sigma_11','Sigma_22','Sigma_12']
141280aadf6SBlaise Bourdin    case(3)
142280aadf6SBlaise Bourdin        numNodalVar = 4
143dcb3e689SBarry Smith        nodalVarName(1:numNodalVar) = ['U_x  ','U_y  ','U_z  ','Alpha']
144280aadf6SBlaise Bourdin        numZonalVar = 6
145dcb3e689SBarry Smith        zonalVarName(1:numZonalVar) = ['Sigma_11','Sigma_22','Sigma_33','Sigma_23','Sigma_13','Sigma_12']
146280aadf6SBlaise Bourdin    case default
147280aadf6SBlaise Bourdin        write(IOBuffer,'("No layout for dimension ",I2)') sdim
148280aadf6SBlaise Bourdin    end select
149d8606c27SBarry Smith    PetscCallA(PetscViewerExodusIIGetId(viewer,exoid,ierr))
150dcb3e689SBarry Smith    PetscCallA(expvp(exoid, 'E', numZonalVar,ierr))
151dcb3e689SBarry Smith    PetscCallA(expvan(exoid, 'E', numZonalVar, zonalVarName,ierr))
152dcb3e689SBarry Smith    PetscCallA(expvp(exoid, 'N', numNodalVar,ierr))
153dcb3e689SBarry Smith    PetscCallA(expvan(exoid, 'N', numNodalVar, nodalVarName,ierr))
154d8606c27SBarry Smith    PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK,numCS,PETSC_NULL_REAL,sjunk,ierr))
155280aadf6SBlaise Bourdin
156280aadf6SBlaise Bourdin    !    An exodusII truth table specifies which fields are saved at which time step
157280aadf6SBlaise Bourdin    !    It speeds up I/O but reserving space for fields in the file ahead of time.
158280aadf6SBlaise Bourdin    allocate(truthtable(numCS,numZonalVar))
159280aadf6SBlaise Bourdin    truthtable = .true.
160d8606c27SBarry Smith    PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr))
161280aadf6SBlaise Bourdin    deallocate(truthtable)
162280aadf6SBlaise Bourdin
163280aadf6SBlaise Bourdin    !   Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
164280aadf6SBlaise Bourdin    do step = 1,numstep
165d8606c27SBarry Smith        PetscCallA(exptim(exoid,step,Real(step,kind=kPR),ierr))
166280aadf6SBlaise Bourdin    end do
167280aadf6SBlaise Bourdin
168d8606c27SBarry Smith    PetscCallA(PetscObjectGetComm(dm,comm,ierr))
169d8606c27SBarry Smith    PetscCallA(PetscSectionCreate(comm, section,ierr))
170d8606c27SBarry Smith    PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr))
171dcb3e689SBarry Smith    PetscCallA(PetscSectionSetFieldName(section, fieldU, 'U',ierr))
172dcb3e689SBarry Smith    PetscCallA(PetscSectionSetFieldName(section, fieldA, 'Alpha',ierr))
173dcb3e689SBarry Smith    PetscCallA(PetscSectionSetFieldName(section, fieldS, 'Sigma',ierr))
174d8606c27SBarry Smith    PetscCallA(DMPlexGetChart(dm, pStart, pEnd,ierr))
175d8606c27SBarry Smith    PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr))
176280aadf6SBlaise Bourdin
177280aadf6SBlaise Bourdin    allocate(pStartDepth(sdim+1))
178280aadf6SBlaise Bourdin    allocate(pEndDepth(sdim+1))
179280aadf6SBlaise Bourdin    do d = 1, sdim+1
180d8606c27SBarry Smith        PetscCallA(DMPlexGetDepthStratum(dm, d-1, pStartDepth(d), pEndDepth(d),ierr))
181280aadf6SBlaise Bourdin    end do
182280aadf6SBlaise Bourdin
183280aadf6SBlaise Bourdin    ! Vector field U, Scalar field Alpha, Tensor field Sigma
184e2739ba6SAlexis Marboeuf    PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr));
185e2739ba6SAlexis Marboeuf    PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr));
186e2739ba6SAlexis Marboeuf    PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr));
187280aadf6SBlaise Bourdin
188280aadf6SBlaise Bourdin    ! Going through cell sets then cells, and setting up storage for the sections
189dcb3e689SBarry Smith    PetscCallA(DMGetLabelSize(dm, 'Cell Sets', numCS, ierr))
190dcb3e689SBarry Smith    PetscCallA(DMGetLabelIdIS(dm, 'Cell Sets', csIS, ierr))
191d8606c27SBarry Smith    PetscCallA(ISGetIndicesF90(csIS, csID, ierr))
192280aadf6SBlaise Bourdin    do set = 1,numCS
1938f2079feSBlaise Bourdin        !!! This is pointless but will avoid a warning with gfortran and -Werror=maybe-uninitialized
1948f2079feSBlaise Bourdin        nullify(dofA)
1958f2079feSBlaise Bourdin        nullify(dofU)
1968f2079feSBlaise Bourdin        nullify(dofS)
197dcb3e689SBarry Smith        PetscCallA(DMGetStratumSize(dm, 'Cell Sets', csID(set), numCells,ierr))
198dcb3e689SBarry Smith        PetscCallA(DMGetStratumIS(dm, 'Cell Sets', csID(set), cellIS,ierr))
199280aadf6SBlaise Bourdin        if (numCells > 0) then
200280aadf6SBlaise Bourdin            select case(sdim)
201280aadf6SBlaise Bourdin            case(2)
202280aadf6SBlaise Bourdin                dofs => dofS2D
203280aadf6SBlaise Bourdin            case(3)
204280aadf6SBlaise Bourdin                dofs => dofS3D
205280aadf6SBlaise Bourdin            case default
206280aadf6SBlaise Bourdin                write(IOBuffer,'("No layout for dimension ",I2)') sdim
207de401611SBlaise Bourdin                SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer)
208280aadf6SBlaise Bourdin            end select ! sdim
209280aadf6SBlaise Bourdin
210280aadf6SBlaise Bourdin            ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
211280aadf6SBlaise Bourdin            ! It will not be enough to identify more exotic elements like pyramid or prisms...  */
212d8606c27SBarry Smith            PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr))
213280aadf6SBlaise Bourdin            nullify(closureA)
214d8606c27SBarry Smith            PetscCallA(DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr))
215280aadf6SBlaise Bourdin            select case(size(closureA)/2)
216280aadf6SBlaise Bourdin            case(7) ! Tri
217280aadf6SBlaise Bourdin                if (order == 1) then
218280aadf6SBlaise Bourdin                    dofU => dofUP1Tri
219280aadf6SBlaise Bourdin                    dofA => dofAP1Tri
220280aadf6SBlaise Bourdin                else
221280aadf6SBlaise Bourdin                    dofU => dofUP2Tri
222280aadf6SBlaise Bourdin                    dofA => dofAP2Tri
223280aadf6SBlaise Bourdin                end if
224280aadf6SBlaise Bourdin            case(9) ! Quad
225280aadf6SBlaise Bourdin                if (order == 1) then
226280aadf6SBlaise Bourdin                    dofU => dofUP1Quad
227280aadf6SBlaise Bourdin                    dofA => dofAP1Quad
228280aadf6SBlaise Bourdin                else
229280aadf6SBlaise Bourdin                    dofU => dofUP2Quad
230280aadf6SBlaise Bourdin                    dofA => dofAP2Quad
231280aadf6SBlaise Bourdin                end if
232280aadf6SBlaise Bourdin            case(15) ! Tet
233280aadf6SBlaise Bourdin                if (order == 1) then
234280aadf6SBlaise Bourdin                    dofU => dofUP1Tet
235280aadf6SBlaise Bourdin                    dofA => dofAP1Tet
236280aadf6SBlaise Bourdin                else
237280aadf6SBlaise Bourdin                    dofU => dofUP2Tet
238280aadf6SBlaise Bourdin                    dofA => dofAP2Tet
239280aadf6SBlaise Bourdin                end if
240280aadf6SBlaise Bourdin            case(27) ! Hex
241280aadf6SBlaise Bourdin                if (order == 1) then
242280aadf6SBlaise Bourdin                    dofU => dofUP1Hex
243280aadf6SBlaise Bourdin                    dofA => dofAP1Hex
244280aadf6SBlaise Bourdin                else
245280aadf6SBlaise Bourdin                    dofU => dofUP2Hex
246280aadf6SBlaise Bourdin                    dofA => dofAP2Hex
247280aadf6SBlaise Bourdin                end if
248280aadf6SBlaise Bourdin            case default
249280aadf6SBlaise Bourdin                write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2
250de401611SBlaise Bourdin                SETERRA(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer)
251280aadf6SBlaise Bourdin            end select
252d8606c27SBarry Smith            PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr))
253280aadf6SBlaise Bourdin            do cell = 1,numCells!
254280aadf6SBlaise Bourdin                nullify(closure)
255d8606c27SBarry Smith                PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr))
256280aadf6SBlaise Bourdin                do p = 1,size(closure),2
257280aadf6SBlaise Bourdin                    ! find the depth of p
258280aadf6SBlaise Bourdin                    do d = 1,sdim+1
259280aadf6SBlaise Bourdin                        if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then
260d8606c27SBarry Smith                            PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr))
261d8606c27SBarry Smith                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr))
262d8606c27SBarry Smith                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr))
263d8606c27SBarry Smith                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr))
264280aadf6SBlaise Bourdin                        end if ! closure(p)
265280aadf6SBlaise Bourdin                    end do ! d
266280aadf6SBlaise Bourdin                end do ! p
267d8606c27SBarry Smith                PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr))
268280aadf6SBlaise Bourdin            end do ! cell
269d8606c27SBarry Smith            PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr))
270d8606c27SBarry Smith            PetscCallA(ISDestroy(cellIS,ierr))
271280aadf6SBlaise Bourdin        end if ! numCells
272280aadf6SBlaise Bourdin    end do ! set
273d8606c27SBarry Smith    PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr))
274d8606c27SBarry Smith    PetscCallA(ISDestroy(csIS,ierr))
275d8606c27SBarry Smith    PetscCallA(PetscSectionSetUp(section,ierr))
276d8606c27SBarry Smith    PetscCallA(DMSetLocalSection(dm, section,ierr))
277dcb3e689SBarry Smith    PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, '-dm_section_view',ierr))
278*f51a5268SBarry Smith    PetscCallA(PetscSectionGetValueLayout(PETSC_COMM_WORLD, section, layout, ierr))
279*f51a5268SBarry Smith    PetscCallA(PetscLayoutDestroy(layout,ierr))
280d8606c27SBarry Smith    PetscCallA(PetscSectionDestroy(section,ierr))
281280aadf6SBlaise Bourdin
282d8606c27SBarry Smith    PetscCallA(DMSetUseNatural(dm,PETSC_TRUE,ierr))
283d8606c27SBarry Smith    PetscCallA(DMPlexGetPartitioner(dm,part,ierr))
284d8606c27SBarry Smith    PetscCallA(PetscPartitionerSetFromOptions(part,ierr))
285d8606c27SBarry Smith    PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr))
286280aadf6SBlaise Bourdin
287280aadf6SBlaise Bourdin    if (numProc > 1) then
288d8606c27SBarry Smith        PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr))
289d8606c27SBarry Smith        PetscCallA(PetscSFDestroy(migrationSF,ierr))
290e2739ba6SAlexis Marboeuf    else
291e2739ba6SAlexis Marboeuf        pdm = dm
292280aadf6SBlaise Bourdin    end if
293dcb3e689SBarry Smith    PetscCallA(DMViewFromOptions(pdm,PETSC_NULL_OPTIONS,'-dm_view',ierr))
294280aadf6SBlaise Bourdin
295280aadf6SBlaise Bourdin    ! Get DM and IS for each field of dm
296e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldU,  isU,  dmU,ierr))
297e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldA,  isA,  dmA,ierr))
298e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 1_kPI, fieldS,  isS,  dmS,ierr))
299e2739ba6SAlexis Marboeuf    PetscCallA(DMCreateSubDM(pdm, 2_kPI, fieldUA, isUA, dmUA,ierr))
300280aadf6SBlaise Bourdin
301280aadf6SBlaise Bourdin    !Create the exodus result file
302280aadf6SBlaise Bourdin    allocate(dmList(2))
303280aadf6SBlaise Bourdin    dmList(1) = dmU;
304280aadf6SBlaise Bourdin    dmList(2) = dmA;
305d8606c27SBarry Smith    PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr))
306280aadf6SBlaise Bourdin    deallocate(dmList)
307280aadf6SBlaise Bourdin
308e2739ba6SAlexis Marboeuf    PetscCallA(DMGetGlobalVector(pdm,   X,ierr))
309d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmU,  U,ierr))
310d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmA,  A,ierr))
311d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmS,  S,ierr))
312d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA, UA,ierr))
313d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr))
314280aadf6SBlaise Bourdin
315dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(U,  'U',ierr))
316dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(A,  'Alpha',ierr))
317dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(S,  'Sigma',ierr))
318dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(UA, 'UAlpha',ierr))
319dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(UA2, 'UAlpha2',ierr))
320d8606c27SBarry Smith    PetscCallA(VecSet(X, -111.0_kPR,ierr))
321280aadf6SBlaise Bourdin
322280aadf6SBlaise 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 */
323d8606c27SBarry Smith    PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr))
324d8606c27SBarry Smith    PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr))
325d8606c27SBarry Smith    PetscCallA(VecGetArrayF90(UALoc, cval,ierr))
326d8606c27SBarry Smith    PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr))
327d8606c27SBarry Smith    PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr))
328d8606c27SBarry Smith    PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr))
329280aadf6SBlaise Bourdin
330280aadf6SBlaise Bourdin    do p = pStart,pEnd-1
331d8606c27SBarry Smith        PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr))
332280aadf6SBlaise Bourdin        if (dofUA > 0) then
333d8606c27SBarry Smith            PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr))
334d8606c27SBarry Smith            PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr))
335280aadf6SBlaise Bourdin            closureSize = size(xyz)
336280aadf6SBlaise Bourdin            do i = 1,sdim
337280aadf6SBlaise Bourdin                do j = 0, closureSize-1,sdim
338280aadf6SBlaise Bourdin                    cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i)
339280aadf6SBlaise Bourdin                end do
340280aadf6SBlaise Bourdin                cval(offUA+i) = cval(offUA+i) * sdim / closureSize;
341280aadf6SBlaise Bourdin                cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2
342280aadf6SBlaise Bourdin            end do
343d8606c27SBarry Smith            PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr))
344280aadf6SBlaise Bourdin        end if
345280aadf6SBlaise Bourdin    end do
346280aadf6SBlaise Bourdin
347d8606c27SBarry Smith    PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr))
348d8606c27SBarry Smith    PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr))
349d8606c27SBarry Smith    PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr))
350d8606c27SBarry Smith    PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr))
351280aadf6SBlaise Bourdin
352280aadf6SBlaise Bourdin    !Update X
353d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr))
354280aadf6SBlaise Bourdin    ! Restrict to U and Alpha
355d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr))
356d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr))
357dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, '-ua_vec_view',ierr))
358dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, '-u_vec_view',ierr))
359dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, '-a_vec_view',ierr))
360280aadf6SBlaise Bourdin    ! restrict to UA2
361d8606c27SBarry Smith    PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr))
362dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, '-ua2_vec_view',ierr))
363280aadf6SBlaise Bourdin
364e2739ba6SAlexis Marboeuf    ! Getting Natural Vec
365d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmU, 0_kPI, time, ierr))
366d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmA, 0_kPI, time, ierr))
367280aadf6SBlaise Bourdin
368d8606c27SBarry Smith    PetscCallA(VecView(U, viewer,ierr))
369d8606c27SBarry Smith    PetscCallA(VecView(A, viewer,ierr))
370280aadf6SBlaise Bourdin
371280aadf6SBlaise Bourdin    !Saving U and Alpha in one shot.
372280aadf6SBlaise Bourdin    !For this, we need to cheat and change the Vec's name
373280aadf6SBlaise Bourdin    !Note that in the end we write variables one component at a time,
374280aadf6SBlaise Bourdin    !so that there is no real value in doing this
375d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr))
376d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr))
377d8606c27SBarry Smith    PetscCallA(VecCopy(UA, tmpVec,ierr))
378dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr))
379d8606c27SBarry Smith    PetscCallA(VecView(tmpVec, viewer,ierr))
380280aadf6SBlaise Bourdin
381280aadf6SBlaise Bourdin    ! Reading nodal variables in Exodus file
382d8606c27SBarry Smith    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
383d8606c27SBarry Smith    PetscCallA(VecLoad(tmpVec, viewer,ierr))
384d8606c27SBarry Smith    PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr))
385d8606c27SBarry Smith    PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr))
386280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
387280aadf6SBlaise Bourdin        write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm
388280aadf6SBlaise Bourdin    end if
389d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr))
390280aadf6SBlaise Bourdin
391280aadf6SBlaise Bourdin    ! same thing with the UA2 Vec obtained from the superDM
392d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr))
393d8606c27SBarry Smith    PetscCallA(VecCopy(UA2, tmpVec,ierr))
394dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(tmpVec, 'U',ierr))
395d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr))
396d8606c27SBarry Smith    PetscCallA(VecView(tmpVec, viewer,ierr))
397280aadf6SBlaise Bourdin
398280aadf6SBlaise Bourdin    ! Reading nodal variables in Exodus file
399d8606c27SBarry Smith    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
400d8606c27SBarry Smith    PetscCallA(VecLoad(tmpVec,viewer,ierr))
401d8606c27SBarry Smith    PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr))
402d8606c27SBarry Smith    PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr))
403280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
404280aadf6SBlaise Bourdin        write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm
405280aadf6SBlaise Bourdin    end if
406d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr))
407280aadf6SBlaise Bourdin
408280aadf6SBlaise Bourdin    ! Building and saving Sigma
409280aadf6SBlaise Bourdin    !   We set sigma_0 = rank (to see partitioning)
410280aadf6SBlaise Bourdin    !          sigma_1 = cell set ID
411280aadf6SBlaise Bourdin    !          sigma_2 = x_coordinate of the cell center of mass
412d8606c27SBarry Smith    PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr))
413d8606c27SBarry Smith    PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr))
414dcb3e689SBarry Smith    PetscCallA(DMGetLabelIdIS(dmS, 'Cell Sets', csIS,ierr))
415dcb3e689SBarry Smith    PetscCallA(DMGetLabelSize(dmS, 'Cell Sets',numCS,ierr))
416d8606c27SBarry Smith    PetscCallA(ISGetIndicesF90(csIS, csID,ierr))
417280aadf6SBlaise Bourdin
418280aadf6SBlaise Bourdin    do set = 1, numCS
419dcb3e689SBarry Smith        PetscCallA(DMGetStratumIS(dmS, 'Cell Sets', csID(set), cellIS,ierr))
420d8606c27SBarry Smith        PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr))
421d8606c27SBarry Smith        PetscCallA(ISGetSize(cellIS, numCells,ierr))
422280aadf6SBlaise Bourdin        do cell = 1,numCells
423d8606c27SBarry Smith            PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr))
424d8606c27SBarry Smith            PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr))
425280aadf6SBlaise Bourdin            cval(1) = rank
426280aadf6SBlaise Bourdin            cval(2) = csID(set)
427280aadf6SBlaise Bourdin            cval(3) = 0.0_kPR
428280aadf6SBlaise Bourdin            do c = 1, size(xyz),sdim
429280aadf6SBlaise Bourdin                cval(3) = cval(3) + xyz(c)
430280aadf6SBlaise Bourdin            end do
431280aadf6SBlaise Bourdin            cval(3) = cval(3) * sdim / size(xyz)
432d8606c27SBarry Smith            PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr))
433d8606c27SBarry Smith            PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr))
434d8606c27SBarry Smith            PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr))
435280aadf6SBlaise Bourdin        end do
436d8606c27SBarry Smith        PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr))
437d8606c27SBarry Smith        PetscCallA(ISDestroy(cellIS,ierr))
438280aadf6SBlaise Bourdin    end do
439d8606c27SBarry Smith    PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr))
440d8606c27SBarry Smith    PetscCallA(ISDestroy(csIS,ierr))
441dcb3e689SBarry Smith    PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, '-s_vec_view',ierr))
442280aadf6SBlaise Bourdin
443280aadf6SBlaise Bourdin    ! Writing zonal variables in Exodus file
444d8606c27SBarry Smith    PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr))
445d8606c27SBarry Smith    PetscCallA(VecView(S,viewer,ierr))
446280aadf6SBlaise Bourdin
447280aadf6SBlaise Bourdin    ! Reading zonal variables in Exodus file */
448d8606c27SBarry Smith    PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr))
449d8606c27SBarry Smith    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
450dcb3e689SBarry Smith    PetscCallA(PetscObjectSetName(tmpVec, 'Sigma',ierr))
451d8606c27SBarry Smith    PetscCallA(VecLoad(tmpVec,viewer,ierr))
452d8606c27SBarry Smith    PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr))
453d8606c27SBarry Smith    PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr))
454280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
455280aadf6SBlaise Bourdin       write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm
456280aadf6SBlaise Bourdin    end if
457d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr))
458280aadf6SBlaise Bourdin
459d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr))
460d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr))
461d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmS,  S,ierr))
462d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmA,  A,ierr))
463d8606c27SBarry Smith    PetscCallA(DMRestoreGlobalVector(dmU,  U,ierr))
464e2739ba6SAlexis Marboeuf    PetscCallA(DMRestoreGlobalVector(pdm,   X,ierr))
465d8606c27SBarry Smith    PetscCallA(DMDestroy(dmU,ierr))
466d8606c27SBarry Smith    PetscCallA(ISDestroy(isU,ierr))
467d8606c27SBarry Smith    PetscCallA(DMDestroy(dmA,ierr))
468d8606c27SBarry Smith    PetscCallA(ISDestroy(isA,ierr))
469d8606c27SBarry Smith    PetscCallA(DMDestroy(dmS,ierr))
470d8606c27SBarry Smith    PetscCallA(ISDestroy(isS,ierr))
471d8606c27SBarry Smith    PetscCallA(DMDestroy(dmUA,ierr))
472d8606c27SBarry Smith    PetscCallA(ISDestroy(isUA,ierr))
473d8606c27SBarry Smith    PetscCallA(DMDestroy(dmUA2,ierr))
474e2739ba6SAlexis Marboeuf    PetscCallA(DMDestroy(pdm,ierr))
475e2739ba6SAlexis Marboeuf    if (numProc > 1) then
476d8606c27SBarry Smith        PetscCallA(DMDestroy(dm,ierr))
477e2739ba6SAlexis Marboeuf    end if
478280aadf6SBlaise Bourdin
479280aadf6SBlaise Bourdin    deallocate(pStartDepth)
480280aadf6SBlaise Bourdin    deallocate(pEndDepth)
481280aadf6SBlaise Bourdin
482d8606c27SBarry Smith    PetscCallA(PetscViewerDestroy(viewer,ierr))
483d8606c27SBarry Smith    PetscCallA(PetscFinalize(ierr))
4842e8d78feSBlaise Bourdinend program ex62f90
485280aadf6SBlaise Bourdin
486280aadf6SBlaise Bourdin! /*TEST
487280aadf6SBlaise Bourdin!
488280aadf6SBlaise Bourdin! build:
489280aadf6SBlaise Bourdin!   requires: exodusii pnetcdf !complex
490280aadf6SBlaise Bourdin!   # 2D seq
491e2739ba6SAlexis Marboeuf! test:
492e2739ba6SAlexis Marboeuf!   suffix: 0
493e2739ba6SAlexis 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
494e2739ba6SAlexis 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
495e2739ba6SAlexis Marboeuf! test:
496e2739ba6SAlexis Marboeuf!   suffix: 1
497e2739ba6SAlexis 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
498e2739ba6SAlexis Marboeuf!
499e2739ba6SAlexis Marboeuf! test:
500e2739ba6SAlexis Marboeuf!   suffix: 2
501e2739ba6SAlexis 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
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: 3
505e2739ba6SAlexis 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
506e2739ba6SAlexis Marboeuf! test:
507e2739ba6SAlexis Marboeuf!   suffix: 4
508e2739ba6SAlexis 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
509e2739ba6SAlexis Marboeuf! test:
510e2739ba6SAlexis Marboeuf!   suffix: 5
511e2739ba6SAlexis 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
512e2739ba6SAlexis Marboeuf!   # 2D par
513e2739ba6SAlexis Marboeuf! test:
514e2739ba6SAlexis Marboeuf!   suffix: 6
515e2739ba6SAlexis Marboeuf!   nsize: 2
516e2739ba6SAlexis 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
517e2739ba6SAlexis 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
518e2739ba6SAlexis Marboeuf! test:
519e2739ba6SAlexis Marboeuf!   suffix: 7
520e2739ba6SAlexis Marboeuf!   nsize: 2
521e2739ba6SAlexis 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
522e2739ba6SAlexis Marboeuf! test:
523e2739ba6SAlexis Marboeuf!   suffix: 8
524e2739ba6SAlexis Marboeuf!   nsize: 2
525e2739ba6SAlexis 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
526e2739ba6SAlexis Marboeuf!   #TODO: bug in PetscCallA(to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
527e2739ba6SAlexis Marboeuf! test:
528e2739ba6SAlexis Marboeuf!   suffix: 9
529e2739ba6SAlexis Marboeuf!   nsize: 2
530e2739ba6SAlexis 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
531e2739ba6SAlexis Marboeuf! test:
532e2739ba6SAlexis Marboeuf!   suffix: 10
533e2739ba6SAlexis Marboeuf!   nsize: 2
534e2739ba6SAlexis 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
535e2739ba6SAlexis Marboeuf! test:
5369c89aa79SPierre Jolivet!   # Something is now broken with parallel read/write for whatever shape H is
537e2739ba6SAlexis Marboeuf!   TODO: broken
538e2739ba6SAlexis Marboeuf!   suffix: 11
539e2739ba6SAlexis Marboeuf!   nsize: 2
540e2739ba6SAlexis 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
541e2739ba6SAlexis Marboeuf
542e2739ba6SAlexis Marboeuf!   #3d seq
543e2739ba6SAlexis Marboeuf! test:
544e2739ba6SAlexis Marboeuf!   suffix: 12
545e2739ba6SAlexis 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
546e2739ba6SAlexis Marboeuf! test:
547e2739ba6SAlexis Marboeuf!   suffix: 13
548e2739ba6SAlexis 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
549e2739ba6SAlexis 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
550e2739ba6SAlexis Marboeuf! test:
551e2739ba6SAlexis Marboeuf!   suffix: 14
552e2739ba6SAlexis 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
553e2739ba6SAlexis Marboeuf! test:
554e2739ba6SAlexis Marboeuf!   suffix: 15
555e2739ba6SAlexis 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
556e2739ba6SAlexis 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
557e2739ba6SAlexis Marboeuf!   #3d par
558e2739ba6SAlexis Marboeuf! test:
559e2739ba6SAlexis Marboeuf!   suffix: 16
560e2739ba6SAlexis Marboeuf!   nsize: 2
561e2739ba6SAlexis 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
562e2739ba6SAlexis Marboeuf! test:
563e2739ba6SAlexis Marboeuf!   suffix: 17
564e2739ba6SAlexis Marboeuf!   nsize: 2
565e2739ba6SAlexis 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
566e2739ba6SAlexis 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
567e2739ba6SAlexis Marboeuf! test:
568e2739ba6SAlexis Marboeuf!   suffix: 18
569e2739ba6SAlexis Marboeuf!   nsize: 2
570e2739ba6SAlexis 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
571e2739ba6SAlexis Marboeuf! test:
572e2739ba6SAlexis Marboeuf!   suffix: 19
573e2739ba6SAlexis Marboeuf!   nsize: 2
574e2739ba6SAlexis 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
575e2739ba6SAlexis 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
576e2739ba6SAlexis Marboeuf
577280aadf6SBlaise Bourdin! TEST*/
578