xref: /petsc/src/dm/impls/plex/tests/ex62f90.F90 (revision 6aad120caa16b1027d343a5f30f73d01448e4dc0)
1280aadf6SBlaise Bourdinprogram ex62f90
2280aadf6SBlaise Bourdin#include "petsc/finclude/petsc.h"
3280aadf6SBlaise Bourdin    use petsc
4280aadf6SBlaise Bourdin    implicit none
5280aadf6SBlaise Bourdin#include "exodusII.inc"
6280aadf6SBlaise Bourdin
7280aadf6SBlaise Bourdin    ! Get the fortran kind associated with PetscInt and PetscReal so that we can use literal constants.
8280aadf6SBlaise Bourdin    PetscInt                           :: dummyPetscInt
9280aadf6SBlaise Bourdin    PetscReal                          :: dummyPetscreal
10280aadf6SBlaise Bourdin    integer,parameter                  :: kPI = kind(dummyPetscInt)
11280aadf6SBlaise Bourdin    integer,parameter                  :: kPR = kind(dummyPetscReal)
12280aadf6SBlaise Bourdin
13280aadf6SBlaise Bourdin    type(tDM)                          :: dm,dmU,dmA,dmS,dmUA,dmUA2,pDM
14280aadf6SBlaise Bourdin    type(tDM),dimension(:),pointer     :: dmList
15280aadf6SBlaise Bourdin    type(tVec)                         :: X,U,A,S,UA,UA2
16280aadf6SBlaise Bourdin    type(tIS)                          :: isU,isA,isS,isUA
17280aadf6SBlaise Bourdin    type(tPetscSection)                :: section
18280aadf6SBlaise Bourdin    PetscInt,dimension(1)              :: fieldU = [0]
19280aadf6SBlaise Bourdin    PetscInt,dimension(1)              :: fieldA = [2]
20280aadf6SBlaise Bourdin    PetscInt,dimension(1)              :: fieldS = [1]
21280aadf6SBlaise Bourdin    PetscInt,dimension(2)              :: fieldUA = [0,2]
22280aadf6SBlaise Bourdin    character(len=PETSC_MAX_PATH_LEN)  :: ifilename,ofilename,IOBuffer
23280aadf6SBlaise Bourdin    integer                            :: exoid = -1
24280aadf6SBlaise Bourdin    type(tIS)                          :: csIS
25280aadf6SBlaise Bourdin    PetscInt,dimension(:),pointer      :: csID
26280aadf6SBlaise Bourdin    PetscInt,dimension(:),pointer      :: pStartDepth,pEndDepth
27280aadf6SBlaise Bourdin    PetscInt                           :: order = 1
28280aadf6SBlaise Bourdin    PetscInt                           :: sdim,d,pStart,pEnd,p,numCS,set,i,j
29280aadf6SBlaise Bourdin    PetscMPIInt                        :: rank,numProc
30280aadf6SBlaise Bourdin    PetscBool                          :: flg
31280aadf6SBlaise Bourdin    PetscErrorCode                     :: ierr
32280aadf6SBlaise Bourdin    MPI_Comm                           :: comm
33280aadf6SBlaise Bourdin    type(tPetscViewer)                 :: viewer
34280aadf6SBlaise Bourdin
35280aadf6SBlaise Bourdin    Character(len=MXSTLN)              :: sJunk
36280aadf6SBlaise Bourdin    PetscInt                           :: numstep = 3, step
37280aadf6SBlaise Bourdin    PetscInt                           :: numNodalVar,numZonalVar
38280aadf6SBlaise Bourdin    character(len=MXSTLN)              :: nodalVarName(4)
39280aadf6SBlaise Bourdin    character(len=MXSTLN)              :: zonalVarName(6)
40280aadf6SBlaise Bourdin    logical,dimension(:,:),pointer     :: truthtable
41280aadf6SBlaise Bourdin
42280aadf6SBlaise Bourdin    type(tIS)                          :: cellIS
43280aadf6SBlaise Bourdin    PetscInt,dimension(:),pointer      :: cellID
44280aadf6SBlaise Bourdin    PetscInt                           :: numCells, cell, closureSize
45280aadf6SBlaise Bourdin    PetscInt,dimension(:),pointer      :: closureA,closure
46280aadf6SBlaise Bourdin
47280aadf6SBlaise Bourdin    type(tPetscSection)                :: sectionUA,coordSection
48280aadf6SBlaise Bourdin    type(tVec)                         :: UALoc,coord
49280aadf6SBlaise Bourdin    PetscScalar,dimension(:),pointer   :: cval,xyz
50280aadf6SBlaise Bourdin    PetscInt                           :: dofUA,offUA,c
51280aadf6SBlaise Bourdin
52280aadf6SBlaise Bourdin    ! dof layout ordered by increasing height in the DAG: cell, face, edge, vertex
53280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofS2D     = [0, 0, 3]
54280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofUP1Tri  = [2, 0, 0]
55280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofAP1Tri  = [1, 0, 0]
56280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofUP2Tri  = [2, 2, 0]
57280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofAP2Tri  = [1, 1, 0]
58280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofUP1Quad = [2, 0, 0]
59280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofAP1Quad = [1, 0, 0]
60280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofUP2Quad = [2, 2, 2]
61280aadf6SBlaise Bourdin    PetscInt,dimension(3),target        :: dofAP2Quad = [1, 1, 1]
62280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofS3D     = [0, 0, 0, 6]
63280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofUP1Tet  = [3, 0, 0, 0]
64280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofAP1Tet  = [1, 0, 0, 0]
65280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofUP2Tet  = [3, 3, 0, 0]
66280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofAP2Tet  = [1, 1, 0, 0]
67280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofUP1Hex  = [3, 0, 0, 0]
68280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofAP1Hex  = [1, 0, 0, 0]
69280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofUP2Hex  = [3, 3, 3, 3]
70280aadf6SBlaise Bourdin    PetscInt,dimension(4),target        :: dofAP2Hex  = [1, 1, 1, 1]
71280aadf6SBlaise Bourdin    PetscInt,dimension(:),pointer       :: dofU,dofA,dofS
72280aadf6SBlaise Bourdin
73280aadf6SBlaise Bourdin    type(tPetscSF)                      :: migrationSF
74280aadf6SBlaise Bourdin    PetscPartitioner                    :: part
75280aadf6SBlaise Bourdin
76280aadf6SBlaise Bourdin    type(tVec)                          :: tmpVec
77280aadf6SBlaise Bourdin    PetscReal                           :: norm
78280aadf6SBlaise Bourdin    PetscReal                           :: time = 1.234_kPR
79280aadf6SBlaise Bourdin
80280aadf6SBlaise Bourdin    call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
81280aadf6SBlaise Bourdin    if (ierr /= 0) then
82280aadf6SBlaise Bourdin      print*,'Unable to initialize PETSc'
83280aadf6SBlaise Bourdin      stop
84280aadf6SBlaise Bourdin    endif
85280aadf6SBlaise Bourdin
86280aadf6SBlaise Bourdin    call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
87280aadf6SBlaise Bourdin    call MPI_Comm_size(PETSC_COMM_WORLD,numProc,ierr)
88280aadf6SBlaise Bourdin    call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-i",ifilename,flg,ierr);CHKERRA(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
92280aadf6SBlaise Bourdin    call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-o",ofilename,flg,ierr);CHKERRA(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
96280aadf6SBlaise Bourdin    call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-order",order,flg,ierr);CHKERRA(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
103ed5e4e85SVaclav Hapla    call DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr);CHKERRA(ierr)
104e600fa54SMatthew G. Knepley    call DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr);CHKERRA(ierr);
105280aadf6SBlaise Bourdin    call DMSetFromOptions(dm,ierr);CHKERRA(ierr);
106280aadf6SBlaise Bourdin    call DMGetDimension(dm, sdim,ierr);CHKERRA(ierr)
107280aadf6SBlaise Bourdin    call DMViewFromOptions(dm, PETSC_NULL_OPTIONS,"-dm_view",ierr);CHKERRA(ierr);
108280aadf6SBlaise Bourdin
109280aadf6SBlaise Bourdin    ! Create the exodus result file
110280aadf6SBlaise Bourdin
111a5b23f4aSJose E. Roman    ! enable exodus debugging information
112280aadf6SBlaise Bourdin    call exopts(EXVRBS+EXDEBG,ierr)
113280aadf6SBlaise Bourdin    ! Create the exodus file
114280aadf6SBlaise Bourdin    call PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr);CHKERRA(ierr)
115280aadf6SBlaise Bourdin    ! The long way would be
116280aadf6SBlaise Bourdin    !
117280aadf6SBlaise Bourdin    ! call PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr);CHKERRA(ierr)
118280aadf6SBlaise Bourdin    ! call PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr);CHKERRA(ierr)
119280aadf6SBlaise Bourdin    ! call PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr);CHKERRA(ierr)
120280aadf6SBlaise Bourdin    ! call PetscViewerFileSetName(viewer,ofilename,ierr);CHKERRA(ierr)
121280aadf6SBlaise Bourdin
122280aadf6SBlaise Bourdin    ! set the mesh order
123280aadf6SBlaise Bourdin    call PetscViewerExodusIISetOrder(viewer,order,ierr);CHKERRA(ierr)
124280aadf6SBlaise Bourdin    call PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr)
125280aadf6SBlaise Bourdin    !
126280aadf6SBlaise Bourdin    !    Notice how the exodus file is actually NOT open at this point (exoid is -1)
127*6aad120cSJose 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
132280aadf6SBlaise Bourdin    call DMView(dm,viewer,ierr);CHKERRA(ierr)
133280aadf6SBlaise Bourdin    call PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(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
152280aadf6SBlaise Bourdin    call PetscViewerExodusIIGetId(viewer,exoid,ierr);CHKERRA(ierr)
153280aadf6SBlaise Bourdin    call expvp(exoid, "E", numZonalVar,ierr)
154280aadf6SBlaise Bourdin    call expvan(exoid, "E", numZonalVar, zonalVarName,ierr)
155280aadf6SBlaise Bourdin    call expvp(exoid, "N", numNodalVar,ierr)
156280aadf6SBlaise Bourdin    call expvan(exoid, "N", numNodalVar, nodalVarName,ierr)
157280aadf6SBlaise Bourdin    call 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.
163280aadf6SBlaise Bourdin    call 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
168280aadf6SBlaise Bourdin        call exptim(exoid,step,Real(step,kind=kPR),ierr)
169280aadf6SBlaise Bourdin    end do
170280aadf6SBlaise Bourdin
171280aadf6SBlaise Bourdin    call DMSetUseNatural(dm,PETSC_TRUE,ierr);CHKERRA(ierr)
172280aadf6SBlaise Bourdin    call DMPlexGetPartitioner(dm,part,ierr);CHKERRA(ierr)
173280aadf6SBlaise Bourdin    call PetscPartitionerSetFromOptions(part,ierr);CHKERRA(ierr)
174280aadf6SBlaise Bourdin    call DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr);CHKERRA(ierr)
175280aadf6SBlaise Bourdin
176280aadf6SBlaise Bourdin    if (numProc > 1) then
177280aadf6SBlaise Bourdin        call DMPlexSetMigrationSF(pdm,migrationSF,ierr);CHKERRA(ierr)
178280aadf6SBlaise Bourdin        call PetscSFDestroy(migrationSF,ierr);CHKERRA(ierr)
179280aadf6SBlaise Bourdin        call DMDestroy(dm,ierr);CHKERRA(ierr)
180280aadf6SBlaise Bourdin        dm = pdm
181280aadf6SBlaise Bourdin    end if
182280aadf6SBlaise Bourdin    call DMViewFromOptions(dm,PETSC_NULL_OPTIONS,"-dm_view",ierr);CHKERRA(ierr)
183280aadf6SBlaise Bourdin
184280aadf6SBlaise Bourdin    call PetscObjectGetComm(dm,comm,ierr);CHKERRA(ierr)
185280aadf6SBlaise Bourdin    call PetscSectionCreate(comm, section,ierr);CHKERRA(ierr)
186280aadf6SBlaise Bourdin    call PetscSectionSetNumFields(section, 3_kPI,ierr);CHKERRA(ierr)
187280aadf6SBlaise Bourdin    call PetscSectionSetFieldName(section, fieldU, "U",ierr);CHKERRA(ierr)
188280aadf6SBlaise Bourdin    call PetscSectionSetFieldName(section, fieldA, "Alpha",ierr);CHKERRA(ierr)
189280aadf6SBlaise Bourdin    call PetscSectionSetFieldName(section, fieldS, "Sigma",ierr);CHKERRA(ierr)
190280aadf6SBlaise Bourdin    call DMPlexGetChart(dm, pStart, pEnd,ierr);CHKERRA(ierr)
191280aadf6SBlaise Bourdin    call PetscSectionSetChart(section, pStart, pEnd,ierr);CHKERRA(ierr)
192280aadf6SBlaise Bourdin
193280aadf6SBlaise Bourdin    allocate(pStartDepth(sdim+1))
194280aadf6SBlaise Bourdin    allocate(pEndDepth(sdim+1))
195280aadf6SBlaise Bourdin    do d = 1, sdim+1
196280aadf6SBlaise Bourdin        call DMPlexGetDepthStratum(dm, d-1, pStartDepth(d), pEndDepth(d),ierr);CHKERRA(ierr)
197280aadf6SBlaise Bourdin    end do
198280aadf6SBlaise Bourdin
199280aadf6SBlaise Bourdin    ! Vector field U, Scalar field Alpha, Tensor field Sigma
200280aadf6SBlaise Bourdin    call PetscSectionSetFieldComponents(section, fieldU, sdim,ierr);CHKERRA(ierr);
201280aadf6SBlaise Bourdin    call PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr);CHKERRA(ierr);
202280aadf6SBlaise Bourdin    call PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr);CHKERRA(ierr);
203280aadf6SBlaise Bourdin
204280aadf6SBlaise Bourdin    ! Going through cell sets then cells, and setting up storage for the sections
205280aadf6SBlaise Bourdin    call DMGetLabelSize(dm, "Cell Sets", numCS, ierr);CHKERRA(ierr)
206280aadf6SBlaise Bourdin    call DMGetLabelIdIS(dm, "Cell Sets", csIS, ierr);CHKERRA(ierr)
207280aadf6SBlaise Bourdin    call ISGetIndicesF90(csIS, csID, ierr);CHKERRA(ierr)
208280aadf6SBlaise Bourdin    do set = 1,numCS
209280aadf6SBlaise Bourdin        call DMGetStratumSize(dm, "Cell Sets", csID(set), numCells,ierr);CHKERRA(ierr)
210280aadf6SBlaise Bourdin        call DMGetStratumIS(dm, "Cell Sets", csID(set), cellIS,ierr);CHKERRA(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...  */
224280aadf6SBlaise Bourdin            call ISGetIndicesF90(cellIS, cellID,ierr);CHKERRA(ierr)
225280aadf6SBlaise Bourdin            nullify(closureA)
226280aadf6SBlaise Bourdin            call DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr);CHKERRA(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
264280aadf6SBlaise Bourdin            call DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr);CHKERRA(ierr)
265280aadf6SBlaise Bourdin            do cell = 1,numCells!
266280aadf6SBlaise Bourdin                nullify(closure)
267280aadf6SBlaise Bourdin                call DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr);CHKERRA(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
272280aadf6SBlaise Bourdin                            call PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr);CHKERRA(ierr)
273280aadf6SBlaise Bourdin                            call PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr);CHKERRA(ierr)
274280aadf6SBlaise Bourdin                            call PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr);CHKERRA(ierr)
275280aadf6SBlaise Bourdin                            call PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr);CHKERRA(ierr)
276280aadf6SBlaise Bourdin                        end if ! closure(p)
277280aadf6SBlaise Bourdin                    end do ! d
278280aadf6SBlaise Bourdin                end do ! p
279280aadf6SBlaise Bourdin                call DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr);CHKERRA(ierr)
280280aadf6SBlaise Bourdin            end do ! cell
281280aadf6SBlaise Bourdin            call ISRestoreIndicesF90(cellIS, cellID,ierr);CHKERRA(ierr)
282280aadf6SBlaise Bourdin            call ISDestroy(cellIS,ierr);CHKERRA(ierr)
283280aadf6SBlaise Bourdin        end if ! numCells
284280aadf6SBlaise Bourdin    end do ! set
285280aadf6SBlaise Bourdin    call ISRestoreIndicesF90(csIS, csID,ierr);CHKERRA(ierr)
286280aadf6SBlaise Bourdin    call ISDestroy(csIS,ierr);CHKERRA(ierr)
287280aadf6SBlaise Bourdin    call PetscSectionSetUp(section,ierr);CHKERRA(ierr)
288280aadf6SBlaise Bourdin    call DMSetLocalSection(dm, section,ierr);CHKERRA(ierr)
289ba213955SJose E. Roman    call PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, "-dm_section_view",ierr);CHKERRQ(ierr)
290280aadf6SBlaise Bourdin    call PetscSectionDestroy(section,ierr);CHKERRA(ierr)
291280aadf6SBlaise Bourdin
292280aadf6SBlaise Bourdin    ! Get DM and IS for each field of dm
293280aadf6SBlaise Bourdin    call DMCreateSubDM(dm, 1_kPI, fieldU,  isU,  dmU,ierr);CHKERRA(ierr)
294280aadf6SBlaise Bourdin    call DMCreateSubDM(dm, 1_kPI, fieldA,  isA,  dmA,ierr);CHKERRA(ierr)
295280aadf6SBlaise Bourdin    call DMCreateSubDM(dm, 1_kPI, fieldS,  isS,  dmS,ierr);CHKERRA(ierr)
296280aadf6SBlaise Bourdin    call DMCreateSubDM(dm, 2_kPI, fieldUA, isUA, dmUA,ierr);CHKERRA(ierr)
297280aadf6SBlaise Bourdin
298280aadf6SBlaise Bourdin    !Create the exodus result file
299280aadf6SBlaise Bourdin    allocate(dmList(2))
300280aadf6SBlaise Bourdin    dmList(1) = dmU;
301280aadf6SBlaise Bourdin    dmList(2) = dmA;
302280aadf6SBlaise Bourdin    call DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr);CHKERRA(ierr)
303280aadf6SBlaise Bourdin    deallocate(dmList)
304280aadf6SBlaise Bourdin
305280aadf6SBlaise Bourdin    call DMGetGlobalVector(dm,   X,ierr);CHKERRA(ierr)
306280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmU,  U,ierr);CHKERRA(ierr)
307280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmA,  A,ierr);CHKERRA(ierr)
308280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmS,  S,ierr);CHKERRA(ierr)
309280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmUA, UA,ierr);CHKERRA(ierr)
310280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmUA2, UA2,ierr);CHKERRA(ierr)
311280aadf6SBlaise Bourdin
312280aadf6SBlaise Bourdin    call PetscObjectSetName(U,  "U",ierr);CHKERRA(ierr)
313280aadf6SBlaise Bourdin    call PetscObjectSetName(A,  "Alpha",ierr);CHKERRA(ierr)
314280aadf6SBlaise Bourdin    call PetscObjectSetName(S,  "Sigma",ierr);CHKERRA(ierr)
315280aadf6SBlaise Bourdin    call PetscObjectSetName(UA, "UAlpha",ierr);CHKERRA(ierr)
316280aadf6SBlaise Bourdin    call PetscObjectSetName(UA2, "UAlpha2",ierr);CHKERRA(ierr)
317280aadf6SBlaise Bourdin    call VecSet(X, -111.0_kPR,ierr);CHKERRA(ierr)
318280aadf6SBlaise Bourdin
319280aadf6SBlaise 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 */
320280aadf6SBlaise Bourdin    call DMGetLocalSection(dmUA, sectionUA,ierr);CHKERRA(ierr)
321280aadf6SBlaise Bourdin    call DMGetLocalVector(dmUA, UALoc,ierr);CHKERRA(ierr)
322280aadf6SBlaise Bourdin    call VecGetArrayF90(UALoc, cval,ierr);CHKERRA(ierr)
323280aadf6SBlaise Bourdin    call DMGetCoordinateSection(dmUA, coordSection,ierr);CHKERRA(ierr)
324280aadf6SBlaise Bourdin    call DMGetCoordinatesLocal(dmUA, coord,ierr);CHKERRA(ierr)
325280aadf6SBlaise Bourdin    call DMPlexGetChart(dmUA, pStart, pEnd,ierr);CHKERRA(ierr)
326280aadf6SBlaise Bourdin
327280aadf6SBlaise Bourdin    do p = pStart,pEnd-1
328280aadf6SBlaise Bourdin        call PetscSectionGetDof(sectionUA, p, dofUA,ierr);CHKERRA(ierr)
329280aadf6SBlaise Bourdin        if (dofUA > 0) then
330280aadf6SBlaise Bourdin            call PetscSectionGetOffset(sectionUA, p, offUA,ierr);CHKERRA(ierr)
331280aadf6SBlaise Bourdin            call DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr);CHKERRA(ierr)
332280aadf6SBlaise Bourdin            closureSize = size(xyz)
333280aadf6SBlaise Bourdin            do i = 1,sdim
334280aadf6SBlaise Bourdin                do j = 0, closureSize-1,sdim
335280aadf6SBlaise Bourdin                    cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i)
336280aadf6SBlaise Bourdin                end do
337280aadf6SBlaise Bourdin                cval(offUA+i) = cval(offUA+i) * sdim / closureSize;
338280aadf6SBlaise Bourdin                cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2
339280aadf6SBlaise Bourdin            end do
340280aadf6SBlaise Bourdin            call DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr);CHKERRA(ierr)
341280aadf6SBlaise Bourdin        end if
342280aadf6SBlaise Bourdin    end do
343280aadf6SBlaise Bourdin
344280aadf6SBlaise Bourdin    call VecRestoreArrayF90(UALoc, cval,ierr);CHKERRA(ierr)
345280aadf6SBlaise Bourdin    call DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr);CHKERRA(ierr)
346280aadf6SBlaise Bourdin    call DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr);CHKERRA(ierr)
347280aadf6SBlaise Bourdin    call DMRestoreLocalVector(dmUA, UALoc,ierr);CHKERRA(ierr)
348280aadf6SBlaise Bourdin
349280aadf6SBlaise Bourdin    !Update X
350280aadf6SBlaise Bourdin    call VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr);CHKERRA(ierr)
351280aadf6SBlaise Bourdin    ! Restrict to U and Alpha
352280aadf6SBlaise Bourdin    call VecISCopy(X, isU, SCATTER_REVERSE, U,ierr);CHKERRA(ierr)
353280aadf6SBlaise Bourdin    call VecISCopy(X, isA, SCATTER_REVERSE, A,ierr);CHKERRA(ierr)
354280aadf6SBlaise Bourdin    call VecViewFromOptions(UA, PETSC_NULL_OPTIONS, "-ua_vec_view",ierr);CHKERRA(ierr)
355280aadf6SBlaise Bourdin    call VecViewFromOptions(U, PETSC_NULL_OPTIONS, "-u_vec_view",ierr);CHKERRA(ierr)
356280aadf6SBlaise Bourdin    call VecViewFromOptions(A, PETSC_NULL_OPTIONS, "-a_vec_view",ierr);CHKERRA(ierr)
357280aadf6SBlaise Bourdin    ! restrict to UA2
358280aadf6SBlaise Bourdin    call VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr);CHKERRA(ierr)
359280aadf6SBlaise Bourdin    call VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, "-ua2_vec_view",ierr);CHKERRA(ierr)
360280aadf6SBlaise Bourdin
361280aadf6SBlaise Bourdin    ! Writing nodal variables to ExodusII file
362280aadf6SBlaise Bourdin    call DMSetOutputSequenceNumber(dmU,0_kPI,time,ierr);CHKERRA(ierr)
363280aadf6SBlaise Bourdin    call DMSetOutputSequenceNumber(dmA,0_kPI,time,ierr);CHKERRA(ierr)
364280aadf6SBlaise Bourdin
365280aadf6SBlaise Bourdin    call VecView(U, viewer,ierr);CHKERRA(ierr)
366280aadf6SBlaise Bourdin    call VecView(A, viewer,ierr);CHKERRA(ierr)
367280aadf6SBlaise Bourdin
368280aadf6SBlaise Bourdin    ! Saving U and Alpha in one shot.
369280aadf6SBlaise Bourdin    ! For this, we need to cheat and change the Vec's name
370280aadf6SBlaise Bourdin    ! Note that in the end we write variables one component at a time,
371280aadf6SBlaise Bourdin    ! so that there is no real value in doing this
372280aadf6SBlaise Bourdin
373280aadf6SBlaise Bourdin    call DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr);CHKERRA(ierr)
374280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmUA, tmpVec,ierr);CHKERRA(ierr)
375280aadf6SBlaise Bourdin    call VecCopy(UA, tmpVec,ierr);CHKERRA(ierr)
376280aadf6SBlaise Bourdin    call PetscObjectSetName(tmpVec, "U",ierr);CHKERRA(ierr)
377280aadf6SBlaise Bourdin    call VecView(tmpVec, viewer,ierr);CHKERRA(ierr)
378280aadf6SBlaise Bourdin    ! Reading nodal variables in Exodus file
379280aadf6SBlaise Bourdin    call VecSet(tmpVec, -1000.0_kPR,ierr);CHKERRA(ierr)
380280aadf6SBlaise Bourdin    call VecLoad(tmpVec, viewer,ierr);CHKERRA(ierr)
381280aadf6SBlaise Bourdin    call VecAXPY(UA, -1.0_kPR, tmpVec,ierr);CHKERRA(ierr)
382280aadf6SBlaise Bourdin    call VecNorm(UA, NORM_INFINITY, norm,ierr);CHKERRA(ierr)
383280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
384280aadf6SBlaise Bourdin        write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm
385280aadf6SBlaise Bourdin    end if
386280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmUA, tmpVec,ierr);CHKERRA(ierr)
387280aadf6SBlaise Bourdin
388280aadf6SBlaise Bourdin    ! ! same thing with the UA2 Vec obtained from the superDM
389280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmUA2, tmpVec,ierr);CHKERRA(ierr)
390280aadf6SBlaise Bourdin    call VecCopy(UA2, tmpVec,ierr);CHKERRA(ierr)
391280aadf6SBlaise Bourdin    call PetscObjectSetName(tmpVec, "U",ierr);CHKERRA(ierr)
392280aadf6SBlaise Bourdin    call DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr);CHKERRA(ierr)
393280aadf6SBlaise Bourdin    call VecView(tmpVec, viewer,ierr);CHKERRA(ierr)
394280aadf6SBlaise Bourdin    ! Reading nodal variables in Exodus file
395280aadf6SBlaise Bourdin    call VecSet(tmpVec, -1000.0_kPR,ierr);CHKERRA(ierr)
396280aadf6SBlaise Bourdin    call VecLoad(tmpVec,viewer,ierr);CHKERRA(ierr)
397280aadf6SBlaise Bourdin    call VecAXPY(UA2, -1.0_kPR, tmpVec,ierr);CHKERRA(ierr)
398280aadf6SBlaise Bourdin    call VecNorm(UA2, NORM_INFINITY, norm,ierr);CHKERRA(ierr)
399280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
400280aadf6SBlaise Bourdin        write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm
401280aadf6SBlaise Bourdin    end if
402280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmUA2, tmpVec,ierr);CHKERRA(ierr)
403280aadf6SBlaise Bourdin
404280aadf6SBlaise Bourdin    ! Building and saving Sigma
405280aadf6SBlaise Bourdin    !   We set sigma_0 = rank (to see partitioning)
406280aadf6SBlaise Bourdin    !          sigma_1 = cell set ID
407280aadf6SBlaise Bourdin    !          sigma_2 = x_coordinate of the cell center of mass
408280aadf6SBlaise Bourdin    call DMGetCoordinateSection(dmS, coordSection,ierr);CHKERRA(ierr)
409280aadf6SBlaise Bourdin    call DMGetCoordinatesLocal(dmS, coord,ierr);CHKERRA(ierr)
410280aadf6SBlaise Bourdin    call DMGetLabelIdIS(dmS, "Cell Sets", csIS,ierr);CHKERRA(ierr)
411280aadf6SBlaise Bourdin    call DMGetLabelSize(dmS, "Cell Sets",numCS,ierr);CHKERRA(ierr)
412280aadf6SBlaise Bourdin    call ISGetIndicesF90(csIS, csID,ierr);CHKERRA(ierr)
413280aadf6SBlaise Bourdin
414280aadf6SBlaise Bourdin    do set = 1, numCS
415280aadf6SBlaise Bourdin        call DMGetStratumIS(dmS, "Cell Sets", csID(set), cellIS,ierr);CHKERRA(ierr)
416280aadf6SBlaise Bourdin        call ISGetIndicesF90(cellIS, cellID,ierr);CHKERRA(ierr)
417280aadf6SBlaise Bourdin        call ISGetSize(cellIS, numCells,ierr);CHKERRA(ierr)
418280aadf6SBlaise Bourdin        do cell = 1,numCells
419280aadf6SBlaise Bourdin            call DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr);CHKERRA(ierr)
420280aadf6SBlaise Bourdin            call DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr);CHKERRA(ierr)
421280aadf6SBlaise Bourdin            cval(1) = rank
422280aadf6SBlaise Bourdin            cval(2) = csID(set)
423280aadf6SBlaise Bourdin            cval(3) = 0.0_kPR
424280aadf6SBlaise Bourdin            do c = 1, size(xyz),sdim
425280aadf6SBlaise Bourdin                cval(3) = cval(3) + xyz(c)
426280aadf6SBlaise Bourdin            end do
427280aadf6SBlaise Bourdin            cval(3) = cval(3) * sdim / size(xyz)
428280aadf6SBlaise Bourdin            call DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr);CHKERRA(ierr)
429280aadf6SBlaise Bourdin            call DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr);CHKERRA(ierr)
430280aadf6SBlaise Bourdin            call DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr);CHKERRA(ierr)
431280aadf6SBlaise Bourdin        end do
432280aadf6SBlaise Bourdin        call ISRestoreIndicesF90(cellIS, cellID,ierr);CHKERRA(ierr)
433280aadf6SBlaise Bourdin        call ISDestroy(cellIS,ierr);CHKERRA(ierr)
434280aadf6SBlaise Bourdin    end do
435280aadf6SBlaise Bourdin    call ISRestoreIndicesF90(csIS, csID,ierr);CHKERRA(ierr)
436280aadf6SBlaise Bourdin    call ISDestroy(csIS,ierr);CHKERRA(ierr)
437280aadf6SBlaise Bourdin    call VecViewFromOptions(S, PETSC_NULL_OPTIONS, "-s_vec_view",ierr);CHKERRA(ierr)
438280aadf6SBlaise Bourdin
439280aadf6SBlaise Bourdin    ! Writing zonal variables in Exodus file
440280aadf6SBlaise Bourdin    call DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr);CHKERRA(ierr)
441280aadf6SBlaise Bourdin    call VecView(S,viewer,ierr);CHKERRA(ierr)
442280aadf6SBlaise Bourdin
443280aadf6SBlaise Bourdin    ! Reading zonal variables in Exodus file */
444280aadf6SBlaise Bourdin    call DMGetGlobalVector(dmS, tmpVec,ierr);CHKERRA(ierr)
445280aadf6SBlaise Bourdin    call VecSet(tmpVec, -1000.0_kPR,ierr);CHKERRA(ierr)
446280aadf6SBlaise Bourdin    call PetscObjectSetName(tmpVec, "Sigma",ierr);CHKERRA(ierr)
447280aadf6SBlaise Bourdin    call VecLoad(tmpVec,viewer,ierr);CHKERRA(ierr)
448280aadf6SBlaise Bourdin    call VecAXPY(S, -1.0_kPR, tmpVec,ierr);CHKERRA(ierr)
449ba213955SJose E. Roman    call VecNorm(S, NORM_INFINITY,norm,ierr);CHKERRQ(ierr)
450280aadf6SBlaise Bourdin    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
451280aadf6SBlaise Bourdin       write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm
452280aadf6SBlaise Bourdin    end if
453280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmS, tmpVec,ierr);CHKERRA(ierr)
454280aadf6SBlaise Bourdin
455280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmUA2, UA2,ierr);CHKERRA(ierr)
456280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmUA, UA,ierr);CHKERRA(ierr)
457280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmS,  S,ierr);CHKERRA(ierr)
458280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmA,  A,ierr);CHKERRA(ierr)
459280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dmU,  U,ierr);CHKERRA(ierr)
460280aadf6SBlaise Bourdin    call DMRestoreGlobalVector(dm,   X,ierr);CHKERRA(ierr)
461280aadf6SBlaise Bourdin    call DMDestroy(dmU,ierr);CHKERRA(ierr);
462280aadf6SBlaise Bourdin    call ISDestroy(isU,ierr);CHKERRA(ierr)
463280aadf6SBlaise Bourdin    call DMDestroy(dmA,ierr);CHKERRA(ierr);
464280aadf6SBlaise Bourdin    call ISDestroy(isA,ierr);CHKERRA(ierr)
465280aadf6SBlaise Bourdin    call DMDestroy(dmS,ierr);CHKERRA(ierr);
466280aadf6SBlaise Bourdin    call ISDestroy(isS,ierr);CHKERRA(ierr)
467280aadf6SBlaise Bourdin    call DMDestroy(dmUA,ierr);CHKERRA(ierr)
468280aadf6SBlaise Bourdin    call ISDestroy(isUA,ierr);CHKERRA(ierr)
469280aadf6SBlaise Bourdin    call DMDestroy(dmUA2,ierr);CHKERRA(ierr)
470280aadf6SBlaise Bourdin    call DMDestroy(dm,ierr);CHKERRA(ierr)
471280aadf6SBlaise Bourdin
472280aadf6SBlaise Bourdin    deallocate(pStartDepth)
473280aadf6SBlaise Bourdin    deallocate(pEndDepth)
474280aadf6SBlaise Bourdin
475280aadf6SBlaise Bourdin    call PetscViewerDestroy(viewer,ierr);CHKERRA(ierr)
476280aadf6SBlaise Bourdin    call PetscFinalize(ierr)
477280aadf6SBlaise Bourdinend program ex62f90
478280aadf6SBlaise Bourdin
479280aadf6SBlaise Bourdin! /*TEST
480280aadf6SBlaise Bourdin!
481280aadf6SBlaise Bourdin! build:
482280aadf6SBlaise Bourdin!   requires: exodusii pnetcdf !complex
483280aadf6SBlaise Bourdin!   # 2D seq
484280aadf6SBlaise Bourdin! test:
485280aadf6SBlaise Bourdin!   suffix: 0
486280aadf6SBlaise Bourdin!   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
487280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
488280aadf6SBlaise Bourdin! test:
489280aadf6SBlaise Bourdin!   suffix: 1
490280aadf6SBlaise Bourdin!   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
491280aadf6SBlaise Bourdin!
492280aadf6SBlaise Bourdin! test:
493280aadf6SBlaise Bourdin!   suffix: 2
494280aadf6SBlaise Bourdin!   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
495280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
496280aadf6SBlaise Bourdin! test:
497280aadf6SBlaise Bourdin!   suffix: 3
498280aadf6SBlaise Bourdin!   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
499280aadf6SBlaise Bourdin! test:
500280aadf6SBlaise Bourdin!   suffix: 4
501280aadf6SBlaise Bourdin!   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
502280aadf6SBlaise Bourdin! test:
503280aadf6SBlaise Bourdin!   suffix: 5
504280aadf6SBlaise Bourdin!   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
505280aadf6SBlaise Bourdin!   # 2D par
506280aadf6SBlaise Bourdin! test:
507280aadf6SBlaise Bourdin!   suffix: 6
508280aadf6SBlaise Bourdin!   nsize: 2
509280aadf6SBlaise Bourdin!   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
510280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
511280aadf6SBlaise Bourdin! test:
512280aadf6SBlaise Bourdin!   suffix: 7
513280aadf6SBlaise Bourdin!   nsize: 2
514280aadf6SBlaise Bourdin!   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
515280aadf6SBlaise Bourdin! test:
516280aadf6SBlaise Bourdin!   suffix: 8
517280aadf6SBlaise Bourdin!   nsize: 2
518280aadf6SBlaise Bourdin!   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
519280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
520280aadf6SBlaise Bourdin! test:
521280aadf6SBlaise Bourdin!   suffix: 9
522280aadf6SBlaise Bourdin!   nsize: 2
523280aadf6SBlaise Bourdin!   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
524280aadf6SBlaise Bourdin! test:
525280aadf6SBlaise Bourdin!   suffix: 10
526280aadf6SBlaise Bourdin!   nsize: 2
527280aadf6SBlaise Bourdin!   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
528280aadf6SBlaise Bourdin! test:
529280aadf6SBlaise Bourdin!   # Something is now broken with parallel read/write for wahtever shape H is
530280aadf6SBlaise Bourdin!   TODO: broken
531280aadf6SBlaise Bourdin!   suffix: 11
532280aadf6SBlaise Bourdin!   nsize: 2
533280aadf6SBlaise Bourdin!   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
534280aadf6SBlaise Bourdin
535280aadf6SBlaise Bourdin!   #3d seq
536280aadf6SBlaise Bourdin! test:
537280aadf6SBlaise Bourdin!   suffix: 12
538280aadf6SBlaise Bourdin!   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
539280aadf6SBlaise Bourdin! test:
540280aadf6SBlaise Bourdin!   suffix: 13
541280aadf6SBlaise Bourdin!   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
542280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
543280aadf6SBlaise Bourdin! test:
544280aadf6SBlaise Bourdin!   suffix: 14
545280aadf6SBlaise Bourdin!   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
546280aadf6SBlaise Bourdin! test:
547280aadf6SBlaise Bourdin!   suffix: 15
548280aadf6SBlaise Bourdin!   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
549280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
550280aadf6SBlaise Bourdin!   #3d par
551280aadf6SBlaise Bourdin! test:
552280aadf6SBlaise Bourdin!   suffix: 16
553280aadf6SBlaise Bourdin!   nsize: 2
554280aadf6SBlaise Bourdin!   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
555280aadf6SBlaise Bourdin! test:
556280aadf6SBlaise Bourdin!   suffix: 17
557280aadf6SBlaise Bourdin!   nsize: 2
558280aadf6SBlaise Bourdin!   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
559280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
560280aadf6SBlaise Bourdin! test:
561280aadf6SBlaise Bourdin!   suffix: 18
562280aadf6SBlaise Bourdin!   nsize: 2
563280aadf6SBlaise Bourdin!   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
564280aadf6SBlaise Bourdin! test:
565280aadf6SBlaise Bourdin!   suffix: 19
566280aadf6SBlaise Bourdin!   nsize: 2
567280aadf6SBlaise Bourdin!   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
568280aadf6SBlaise Bourdin!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
569280aadf6SBlaise Bourdin! TEST*/
570