xref: /petsc/src/dm/impls/plex/tests/ex62f90.F90 (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
1program ex62f90
2#include "petsc/finclude/petsc.h"
3    use petsc
4    implicit none
5#include "exodusII.inc"
6
7    ! Get the fortran kind associated with PetscInt and PetscReal so that we can use literal constants.
8    PetscInt                           :: dummyPetscInt
9    PetscReal                          :: dummyPetscreal
10    integer,parameter                  :: kPI = kind(dummyPetscInt)
11    integer,parameter                  :: kPR = kind(dummyPetscReal)
12
13    type(tDM)                          :: dm,dmU,dmA,dmS,dmUA,dmUA2,pDM
14    type(tDM),dimension(:),pointer     :: dmList
15    type(tVec)                         :: X,U,A,S,UA,UA2
16    type(tIS)                          :: isU,isA,isS,isUA
17    type(tPetscSection)                :: section
18    PetscInt,dimension(1)              :: fieldU = [0]
19    PetscInt,dimension(1)              :: fieldA = [2]
20    PetscInt,dimension(1)              :: fieldS = [1]
21    PetscInt,dimension(2)              :: fieldUA = [0,2]
22    character(len=PETSC_MAX_PATH_LEN)  :: ifilename,ofilename,IOBuffer
23    integer                            :: exoid = -1
24    type(tIS)                          :: csIS
25    PetscInt,dimension(:),pointer      :: csID
26    PetscInt,dimension(:),pointer      :: pStartDepth,pEndDepth
27    PetscInt                           :: order = 1
28    PetscInt                           :: sdim,d,pStart,pEnd,p,numCS,set,i,j
29    PetscMPIInt                        :: rank,numProc
30    PetscBool                          :: flg
31    PetscErrorCode                     :: ierr
32    MPI_Comm                           :: comm
33    type(tPetscViewer)                 :: viewer
34
35    Character(len=MXSTLN)              :: sJunk
36    PetscInt                           :: numstep = 3, step
37    PetscInt                           :: numNodalVar,numZonalVar
38    character(len=MXSTLN)              :: nodalVarName(4)
39    character(len=MXSTLN)              :: zonalVarName(6)
40    logical,dimension(:,:),pointer     :: truthtable
41
42    type(tIS)                          :: cellIS
43    PetscInt,dimension(:),pointer      :: cellID
44    PetscInt                           :: numCells, cell, closureSize
45    PetscInt,dimension(:),pointer      :: closureA,closure
46
47    type(tPetscSection)                :: sectionUA,coordSection
48    type(tVec)                         :: UALoc,coord
49    PetscScalar,dimension(:),pointer   :: cval,xyz
50    PetscInt                           :: dofUA,offUA,c
51
52    ! dof layout ordered by increasing height in the DAG: cell, face, edge, vertex
53    PetscInt,dimension(3),target        :: dofS2D     = [0, 0, 3]
54    PetscInt,dimension(3),target        :: dofUP1Tri  = [2, 0, 0]
55    PetscInt,dimension(3),target        :: dofAP1Tri  = [1, 0, 0]
56    PetscInt,dimension(3),target        :: dofUP2Tri  = [2, 2, 0]
57    PetscInt,dimension(3),target        :: dofAP2Tri  = [1, 1, 0]
58    PetscInt,dimension(3),target        :: dofUP1Quad = [2, 0, 0]
59    PetscInt,dimension(3),target        :: dofAP1Quad = [1, 0, 0]
60    PetscInt,dimension(3),target        :: dofUP2Quad = [2, 2, 2]
61    PetscInt,dimension(3),target        :: dofAP2Quad = [1, 1, 1]
62    PetscInt,dimension(4),target        :: dofS3D     = [0, 0, 0, 6]
63    PetscInt,dimension(4),target        :: dofUP1Tet  = [3, 0, 0, 0]
64    PetscInt,dimension(4),target        :: dofAP1Tet  = [1, 0, 0, 0]
65    PetscInt,dimension(4),target        :: dofUP2Tet  = [3, 3, 0, 0]
66    PetscInt,dimension(4),target        :: dofAP2Tet  = [1, 1, 0, 0]
67    PetscInt,dimension(4),target        :: dofUP1Hex  = [3, 0, 0, 0]
68    PetscInt,dimension(4),target        :: dofAP1Hex  = [1, 0, 0, 0]
69    PetscInt,dimension(4),target        :: dofUP2Hex  = [3, 3, 3, 3]
70    PetscInt,dimension(4),target        :: dofAP2Hex  = [1, 1, 1, 1]
71    PetscInt,dimension(:),pointer       :: dofU,dofA,dofS
72
73    type(tPetscSF)                      :: migrationSF
74    PetscPartitioner                    :: part
75
76    type(tVec)                          :: tmpVec
77    PetscReal                           :: norm
78    PetscReal                           :: time = 1.234_kPR
79
80    PetscCallA(PetscInitialize(ierr))
81
82    PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
83    PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,numProc,ierr))
84    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-i",ifilename,flg,ierr))
85    if (.not. flg) then
86        SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing input file name -i <input file name>")
87    end if
88    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-o",ofilename,flg,ierr))
89    if (.not. flg) then
90        SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing output file name -o <output file name>")
91    end if
92    PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-order",order,flg,ierr))
93    if ((order > 2) .or. (order < 1)) then
94        write(IOBuffer,'("Unsupported polynomial order ", I2, " not in [1,2]")') order
95        SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,IOBuffer)
96    end if
97
98    ! Read the mesh in any supported format
99    PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr))
100    PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr))
101    PetscCallA(DMSetFromOptions(dm,ierr))
102    PetscCallA(DMGetDimension(dm, sdim,ierr))
103    PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,"-dm_view",ierr))
104
105    ! Create the exodus result file
106
107    ! enable exodus debugging information
108    PetscCallA(exopts(EXVRBS+EXDEBG,ierr))
109    ! Create the exodus file
110    PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr))
111    ! The long way would be
112    !
113    ! PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr))
114    ! PetscCallA(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII,ierr))
115    ! PetscCallA(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE,ierr))
116    ! PetscCallA(PetscViewerFileSetName(viewer,ofilename,ierr))
117
118    ! set the mesh order
119    PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr))
120    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
121    !
122    !    Notice how the exodus file is actually NOT open at this point (exoid is -1)
123    !    Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to
124    !    write the geometry (the DM), which can only be done on a brand new file.
125    !
126
127    ! Save the geometry to the file, erasing all previous content
128    PetscCallA(DMView(dm,viewer,ierr))
129    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
130    !
131    !    Note how the exodus file is now open
132    !
133    ! "Format" the exodus result file, i.e. allocate space for nodal and zonal variables
134    select case(sdim)
135    case(2)
136        numNodalVar = 3
137        nodalVarName(1:numNodalVar) = ["U_x  ","U_y  ","Alpha"]
138        numZonalVar = 3
139        zonalVarName(1:numZonalVar) = ["Sigma_11","Sigma_22","Sigma_12"]
140    case(3)
141        numNodalVar = 4
142        nodalVarName(1:numNodalVar) = ["U_x  ","U_y  ","U_z  ","Alpha"]
143        numZonalVar = 6
144        zonalVarName(1:numZonalVar) = ["Sigma_11","Sigma_22","Sigma_33","Sigma_23","Sigma_13","Sigma_12"]
145    case default
146        write(IOBuffer,'("No layout for dimension ",I2)') sdim
147    end select
148    PetscCallA(PetscViewerExodusIIGetId(viewer,exoid,ierr))
149    PetscCallA(expvp(exoid, "E", numZonalVar,ierr))
150    PetscCallA(expvan(exoid, "E", numZonalVar, zonalVarName,ierr))
151    PetscCallA(expvp(exoid, "N", numNodalVar,ierr))
152    PetscCallA(expvan(exoid, "N", numNodalVar, nodalVarName,ierr))
153    PetscCallA(exinq(exoid, EX_INQ_ELEM_BLK,numCS,PETSC_NULL_REAL,sjunk,ierr))
154
155    !    An exodusII truth table specifies which fields are saved at which time step
156    !    It speeds up I/O but reserving space for fields in the file ahead of time.
157    allocate(truthtable(numCS,numZonalVar))
158    truthtable = .true.
159    PetscCallA(expvtt(exoid, numCS, numZonalVar, truthtable, ierr))
160    deallocate(truthtable)
161
162    !   Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
163    do step = 1,numstep
164        PetscCallA(exptim(exoid,step,Real(step,kind=kPR),ierr))
165    end do
166
167    PetscCallA(DMSetUseNatural(dm,PETSC_TRUE,ierr))
168    PetscCallA(DMPlexGetPartitioner(dm,part,ierr))
169    PetscCallA(PetscPartitionerSetFromOptions(part,ierr))
170    PetscCallA(DMPlexDistribute(dm,0_kPI,migrationSF,pdm,ierr))
171
172    if (numProc > 1) then
173        PetscCallA(DMPlexSetMigrationSF(pdm,migrationSF,ierr))
174        PetscCallA(PetscSFDestroy(migrationSF,ierr))
175        PetscCallA(DMDestroy(dm,ierr))
176        dm = pdm
177    end if
178    PetscCallA(DMViewFromOptions(dm,PETSC_NULL_OPTIONS,"-dm_view",ierr))
179
180    PetscCallA(PetscObjectGetComm(dm,comm,ierr))
181    PetscCallA(PetscSectionCreate(comm, section,ierr))
182    PetscCallA(PetscSectionSetNumFields(section, 3_kPI,ierr))
183    PetscCallA(PetscSectionSetFieldName(section, fieldU, "U",ierr))
184    PetscCallA(PetscSectionSetFieldName(section, fieldA, "Alpha",ierr))
185    PetscCallA(PetscSectionSetFieldName(section, fieldS, "Sigma",ierr))
186    PetscCallA(DMPlexGetChart(dm, pStart, pEnd,ierr))
187    PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr))
188
189    allocate(pStartDepth(sdim+1))
190    allocate(pEndDepth(sdim+1))
191    do d = 1, sdim+1
192        PetscCallA(DMPlexGetDepthStratum(dm, d-1, pStartDepth(d), pEndDepth(d),ierr))
193    end do
194
195    ! Vector field U, Scalar field Alpha, Tensor field Sigma
196    PetscCallA(PetscSectionSetFieldComponents(section, fieldU, sdim,ierr))
197    PetscCallA(PetscSectionSetFieldComponents(section, fieldA, 1_kPI,ierr))
198    PetscCallA(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2,ierr))
199
200    ! Going through cell sets then cells, and setting up storage for the sections
201    PetscCallA(DMGetLabelSize(dm, "Cell Sets", numCS, ierr))
202    PetscCallA(DMGetLabelIdIS(dm, "Cell Sets", csIS, ierr))
203    PetscCallA(ISGetIndicesF90(csIS, csID, ierr))
204    do set = 1,numCS
205        PetscCallA(DMGetStratumSize(dm, "Cell Sets", csID(set), numCells,ierr))
206        PetscCallA(DMGetStratumIS(dm, "Cell Sets", csID(set), cellIS,ierr))
207        if (numCells > 0) then
208            select case(sdim)
209            case(2)
210                dofs => dofS2D
211            case(3)
212                dofs => dofS3D
213            case default
214                write(IOBuffer,'("No layout for dimension ",I2)') sdim
215                SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,IOBuffer)
216            end select ! sdim
217
218            ! Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
219            ! It will not be enough to identify more exotic elements like pyramid or prisms...  */
220            PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr))
221            nullify(closureA)
222            PetscCallA(DMPlexGetTransitiveClosure(dm,cellID(1), PETSC_TRUE, closureA,ierr))
223            select case(size(closureA)/2)
224            case(7) ! Tri
225                if (order == 1) then
226                    dofU => dofUP1Tri
227                    dofA => dofAP1Tri
228                else
229                    dofU => dofUP2Tri
230                    dofA => dofAP2Tri
231                end if
232            case(9) ! Quad
233                if (order == 1) then
234                    dofU => dofUP1Quad
235                    dofA => dofAP1Quad
236                else
237                    dofU => dofUP2Quad
238                    dofA => dofAP2Quad
239                end if
240            case(15) ! Tet
241                if (order == 1) then
242                    dofU => dofUP1Tet
243                    dofA => dofAP1Tet
244                else
245                    dofU => dofUP2Tet
246                    dofA => dofAP2Tet
247                end if
248            case(27) ! Hex
249                if (order == 1) then
250                    dofU => dofUP1Hex
251                    dofA => dofAP1Hex
252                else
253                    dofU => dofUP2Hex
254                    dofA => dofAP2Hex
255                end if
256            case default
257                write(IOBuffer,'("Unknown element with closure size ",I2)') size(closureA)/2
258                SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,IOBuffer)
259            end select
260            PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(1), PETSC_TRUE,closureA,ierr))
261            do cell = 1,numCells!
262                nullify(closure)
263                PetscCallA(DMPlexGetTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr))
264                do p = 1,size(closure),2
265                    ! find the depth of p
266                    do d = 1,sdim+1
267                        if ((closure(p) >= pStartDepth(d)) .and. (closure(p) < pEndDepth(d))) then
268                            PetscCallA(PetscSectionSetDof(section, closure(p), dofU(d)+dofA(d)+dofS(d),ierr))
269                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldU, dofU(d),ierr))
270                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldA, dofA(d),ierr))
271                            PetscCallA(PetscSectionSetFieldDof(section, closure(p), fieldS, dofS(d),ierr))
272                        end if ! closure(p)
273                    end do ! d
274                end do ! p
275                PetscCallA(DMPlexRestoreTransitiveClosure(dm, cellID(cell), PETSC_TRUE, closure,ierr))
276            end do ! cell
277            PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr))
278            PetscCallA(ISDestroy(cellIS,ierr))
279        end if ! numCells
280    end do ! set
281    PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr))
282    PetscCallA(ISDestroy(csIS,ierr))
283    PetscCallA(PetscSectionSetUp(section,ierr))
284    PetscCallA(DMSetLocalSection(dm, section,ierr))
285    PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_SECTION, "-dm_section_view",ierr))
286    PetscCallA(PetscSectionDestroy(section,ierr))
287
288    ! Get DM and IS for each field of dm
289    PetscCallA(DMCreateSubDM(dm, 1_kPI, fieldU,  isU,  dmU,ierr))
290    PetscCallA(DMCreateSubDM(dm, 1_kPI, fieldA,  isA,  dmA,ierr))
291    PetscCallA(DMCreateSubDM(dm, 1_kPI, fieldS,  isS,  dmS,ierr))
292    PetscCallA(DMCreateSubDM(dm, 2_kPI, fieldUA, isUA, dmUA,ierr))
293
294    !Create the exodus result file
295    allocate(dmList(2))
296    dmList(1) = dmU;
297    dmList(2) = dmA;
298    PetscCallA(DMCreateSuperDM(dmList,2_kPI,PETSC_NULL_IS,dmUA2,ierr))
299    deallocate(dmList)
300
301    PetscCallA(DMGetGlobalVector(dm,   X,ierr))
302    PetscCallA(DMGetGlobalVector(dmU,  U,ierr))
303    PetscCallA(DMGetGlobalVector(dmA,  A,ierr))
304    PetscCallA(DMGetGlobalVector(dmS,  S,ierr))
305    PetscCallA(DMGetGlobalVector(dmUA, UA,ierr))
306    PetscCallA(DMGetGlobalVector(dmUA2, UA2,ierr))
307
308    PetscCallA(PetscObjectSetName(U,  "U",ierr))
309    PetscCallA(PetscObjectSetName(A,  "Alpha",ierr))
310    PetscCallA(PetscObjectSetName(S,  "Sigma",ierr))
311    PetscCallA(PetscObjectSetName(UA, "UAlpha",ierr))
312    PetscCallA(PetscObjectSetName(UA2, "UAlpha2",ierr))
313    PetscCallA(VecSet(X, -111.0_kPR,ierr))
314
315    ! 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 */
316    PetscCallA(DMGetLocalSection(dmUA, sectionUA,ierr))
317    PetscCallA(DMGetLocalVector(dmUA, UALoc,ierr))
318    PetscCallA(VecGetArrayF90(UALoc, cval,ierr))
319    PetscCallA(DMGetCoordinateSection(dmUA, coordSection,ierr))
320    PetscCallA(DMGetCoordinatesLocal(dmUA, coord,ierr))
321    PetscCallA(DMPlexGetChart(dmUA, pStart, pEnd,ierr))
322
323    do p = pStart,pEnd-1
324        PetscCallA(PetscSectionGetDof(sectionUA, p, dofUA,ierr))
325        if (dofUA > 0) then
326            PetscCallA(PetscSectionGetOffset(sectionUA, p, offUA,ierr))
327            PetscCallA(DMPlexVecGetClosure(dmUA, coordSection, coord, p, xyz,ierr))
328            closureSize = size(xyz)
329            do i = 1,sdim
330                do j = 0, closureSize-1,sdim
331                    cval(offUA+i) = cval(offUA+i) + xyz(j/sdim+i)
332                end do
333                cval(offUA+i) = cval(offUA+i) * sdim / closureSize;
334                cval(offUA+sdim+1) = cval(offUA+sdim+1) + cval(offUA+i)**2
335            end do
336            PetscCallA(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, xyz,ierr))
337        end if
338    end do
339
340    PetscCallA(VecRestoreArrayF90(UALoc, cval,ierr))
341    PetscCallA(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA,ierr))
342    PetscCallA(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA,ierr))
343    PetscCallA(DMRestoreLocalVector(dmUA, UALoc,ierr))
344
345    !Update X
346    PetscCallA(VecISCopy(X, isUA, SCATTER_FORWARD, UA,ierr))
347    ! Restrict to U and Alpha
348    PetscCallA(VecISCopy(X, isU, SCATTER_REVERSE, U,ierr))
349    PetscCallA(VecISCopy(X, isA, SCATTER_REVERSE, A,ierr))
350    PetscCallA(VecViewFromOptions(UA, PETSC_NULL_OPTIONS, "-ua_vec_view",ierr))
351    PetscCallA(VecViewFromOptions(U, PETSC_NULL_OPTIONS, "-u_vec_view",ierr))
352    PetscCallA(VecViewFromOptions(A, PETSC_NULL_OPTIONS, "-a_vec_view",ierr))
353    ! restrict to UA2
354    PetscCallA(VecISCopy(X, isUA, SCATTER_REVERSE, UA2,ierr))
355    PetscCallA(VecViewFromOptions(UA2, PETSC_NULL_OPTIONS, "-ua2_vec_view",ierr))
356
357    ! Writing nodal variables to ExodusII file
358    PetscCallA(DMSetOutputSequenceNumber(dmU,0_kPI,time,ierr))
359    PetscCallA(DMSetOutputSequenceNumber(dmA,0_kPI,time,ierr))
360
361    PetscCallA(VecView(U, viewer,ierr))
362    PetscCallA(VecView(A, viewer,ierr))
363
364    ! Saving U and Alpha in one shot.
365    ! For this, we need to cheat and change the Vec's name
366    ! Note that in the end we write variables one component at a time,
367    ! so that there is no real value in doing this
368
369    PetscCallA(DMSetOutputSequenceNumber(dmUA,1_kPI,time,ierr))
370    PetscCallA(DMGetGlobalVector(dmUA, tmpVec,ierr))
371    PetscCallA(VecCopy(UA, tmpVec,ierr))
372    PetscCallA(PetscObjectSetName(tmpVec, "U",ierr))
373    PetscCallA(VecView(tmpVec, viewer,ierr))
374    ! Reading nodal variables in Exodus file
375    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
376    PetscCallA(VecLoad(tmpVec, viewer,ierr))
377    PetscCallA(VecAXPY(UA, -1.0_kPR, tmpVec,ierr))
378    PetscCallA(VecNorm(UA, NORM_INFINITY, norm,ierr))
379    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
380        write(IOBuffer,'("UAlpha ||Vin - Vout|| = ",ES12.5)') norm
381    end if
382    PetscCallA(DMRestoreGlobalVector(dmUA, tmpVec,ierr))
383
384    ! ! same thing with the UA2 Vec obtained from the superDM
385    PetscCallA(DMGetGlobalVector(dmUA2, tmpVec,ierr))
386    PetscCallA(VecCopy(UA2, tmpVec,ierr))
387    PetscCallA(PetscObjectSetName(tmpVec, "U",ierr))
388    PetscCallA(DMSetOutputSequenceNumber(dmUA2,2_kPI,time,ierr))
389    PetscCallA(VecView(tmpVec, viewer,ierr))
390    ! Reading nodal variables in Exodus file
391    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
392    PetscCallA(VecLoad(tmpVec,viewer,ierr))
393    PetscCallA(VecAXPY(UA2, -1.0_kPR, tmpVec,ierr))
394    PetscCallA(VecNorm(UA2, NORM_INFINITY, norm,ierr))
395    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
396        write(IOBuffer,'("UAlpha2 ||Vin - Vout|| = ",ES12.5)') norm
397    end if
398    PetscCallA(DMRestoreGlobalVector(dmUA2, tmpVec,ierr))
399
400    ! Building and saving Sigma
401    !   We set sigma_0 = rank (to see partitioning)
402    !          sigma_1 = cell set ID
403    !          sigma_2 = x_coordinate of the cell center of mass
404    PetscCallA(DMGetCoordinateSection(dmS, coordSection,ierr))
405    PetscCallA(DMGetCoordinatesLocal(dmS, coord,ierr))
406    PetscCallA(DMGetLabelIdIS(dmS, "Cell Sets", csIS,ierr))
407    PetscCallA(DMGetLabelSize(dmS, "Cell Sets",numCS,ierr))
408    PetscCallA(ISGetIndicesF90(csIS, csID,ierr))
409
410    do set = 1, numCS
411        PetscCallA(DMGetStratumIS(dmS, "Cell Sets", csID(set), cellIS,ierr))
412        PetscCallA(ISGetIndicesF90(cellIS, cellID,ierr))
413        PetscCallA(ISGetSize(cellIS, numCells,ierr))
414        do cell = 1,numCells
415            PetscCallA(DMPlexVecGetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr))
416            PetscCallA(DMPlexVecGetClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr))
417            cval(1) = rank
418            cval(2) = csID(set)
419            cval(3) = 0.0_kPR
420            do c = 1, size(xyz),sdim
421                cval(3) = cval(3) + xyz(c)
422            end do
423            cval(3) = cval(3) * sdim / size(xyz)
424            PetscCallA(DMPlexVecSetClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval, INSERT_ALL_VALUES,ierr))
425            PetscCallA(DMPlexVecRestoreClosure(dmS, PETSC_NULL_SECTION, S, cellID(cell), cval,ierr))
426            PetscCallA(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID(cell), xyz,ierr))
427        end do
428        PetscCallA(ISRestoreIndicesF90(cellIS, cellID,ierr))
429        PetscCallA(ISDestroy(cellIS,ierr))
430    end do
431    PetscCallA(ISRestoreIndicesF90(csIS, csID,ierr))
432    PetscCallA(ISDestroy(csIS,ierr))
433    PetscCallA(VecViewFromOptions(S, PETSC_NULL_OPTIONS, "-s_vec_view",ierr))
434
435    ! Writing zonal variables in Exodus file
436    PetscCallA(DMSetOutputSequenceNumber(dmS,0_kPI,time,ierr))
437    PetscCallA(VecView(S,viewer,ierr))
438
439    ! Reading zonal variables in Exodus file */
440    PetscCallA(DMGetGlobalVector(dmS, tmpVec,ierr))
441    PetscCallA(VecSet(tmpVec, -1000.0_kPR,ierr))
442    PetscCallA(PetscObjectSetName(tmpVec, "Sigma",ierr))
443    PetscCallA(VecLoad(tmpVec,viewer,ierr))
444    PetscCallA(VecAXPY(S, -1.0_kPR, tmpVec,ierr))
445    PetscCallA(VecNorm(S, NORM_INFINITY,norm,ierr))
446    if (norm > PETSC_SQRT_MACHINE_EPSILON) then
447       write(IOBuffer,'("Sigma ||Vin - Vout|| = ",ES12.5)') norm
448    end if
449    PetscCallA(DMRestoreGlobalVector(dmS, tmpVec,ierr))
450
451    PetscCallA(DMRestoreGlobalVector(dmUA2, UA2,ierr))
452    PetscCallA(DMRestoreGlobalVector(dmUA, UA,ierr))
453    PetscCallA(DMRestoreGlobalVector(dmS,  S,ierr))
454    PetscCallA(DMRestoreGlobalVector(dmA,  A,ierr))
455    PetscCallA(DMRestoreGlobalVector(dmU,  U,ierr))
456    PetscCallA(DMRestoreGlobalVector(dm,   X,ierr))
457    PetscCallA(DMDestroy(dmU,ierr))
458    PetscCallA(ISDestroy(isU,ierr))
459    PetscCallA(DMDestroy(dmA,ierr))
460    PetscCallA(ISDestroy(isA,ierr))
461    PetscCallA(DMDestroy(dmS,ierr))
462    PetscCallA(ISDestroy(isS,ierr))
463    PetscCallA(DMDestroy(dmUA,ierr))
464    PetscCallA(ISDestroy(isUA,ierr))
465    PetscCallA(DMDestroy(dmUA2,ierr))
466    PetscCallA(DMDestroy(dm,ierr))
467
468    deallocate(pStartDepth)
469    deallocate(pEndDepth)
470
471    PetscCallA(PetscViewerDestroy(viewer,ierr))
472    PetscCallA(PetscFinalize(ierr))
473end program ex62f90
474
475! /*TEST
476!
477! build:
478!   requires: exodusii pnetcdf !complex
479!   # 2D seq
480! test:
481!   suffix: 0
482!   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
483!   #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
484! test:
485!   suffix: 1
486!   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
487!
488! test:
489!   suffix: 2
490!   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
491!   #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
492! test:
493!   suffix: 3
494!   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
495! test:
496!   suffix: 4
497!   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
498! test:
499!   suffix: 5
500!   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
501!   # 2D par
502! test:
503!   suffix: 6
504!   nsize: 2
505!   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
506!   #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
507! test:
508!   suffix: 7
509!   nsize: 2
510!   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
511! test:
512!   suffix: 8
513!   nsize: 2
514!   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
515!   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
516! test:
517!   suffix: 9
518!   nsize: 2
519!   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
520! test:
521!   suffix: 10
522!   nsize: 2
523!   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
524! test:
525!   # Something is now broken with parallel read/write for wahtever shape H is
526!   TODO: broken
527!   suffix: 11
528!   nsize: 2
529!   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
530
531!   #3d seq
532! test:
533!   suffix: 12
534!   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
535! test:
536!   suffix: 13
537!   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
538!   #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
539! test:
540!   suffix: 14
541!   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
542! test:
543!   suffix: 15
544!   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
545!   #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
546!   #3d par
547! test:
548!   suffix: 16
549!   nsize: 2
550!   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
551! test:
552!   suffix: 17
553!   nsize: 2
554!   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
555!   #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
556! test:
557!   suffix: 18
558!   nsize: 2
559!   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
560! test:
561!   suffix: 19
562!   nsize: 2
563!   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
564!   #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
565! TEST*/
566