xref: /petsc/src/dm/impls/plex/tests/ex95f90.F90 (revision 49c89c767329b4815515fe9047927c4fa12c2a65)
1*49c89c76SBlaise Bourdinprogram ex95f90
2*49c89c76SBlaise Bourdin#include "petsc/finclude/petsc.h"
3*49c89c76SBlaise Bourdin    use petsc
4*49c89c76SBlaise Bourdin    implicit none
5*49c89c76SBlaise Bourdin#include "exodusII.inc"
6*49c89c76SBlaise Bourdin
7*49c89c76SBlaise Bourdin    ! Get the Fortran kind associated with PetscInt and PetscReal so that we can use literal constants.
8*49c89c76SBlaise Bourdin    PetscInt                             :: dummyPetscInt
9*49c89c76SBlaise Bourdin    PetscReal                            :: dummyPetscreal
10*49c89c76SBlaise Bourdin    PetscBool                            :: flg
11*49c89c76SBlaise Bourdin    integer,parameter                    :: kPI = kind(dummyPetscInt)
12*49c89c76SBlaise Bourdin    integer,parameter                    :: kPR = kind(dummyPetscReal)
13*49c89c76SBlaise Bourdin    integer                              :: nNodalVar = 4
14*49c89c76SBlaise Bourdin    integer                              :: nZonalVar = 3
15*49c89c76SBlaise Bourdin    integer                              :: i
16*49c89c76SBlaise Bourdin
17*49c89c76SBlaise Bourdin    PetscErrorCode                       :: ierr
18*49c89c76SBlaise Bourdin    type(tDM)                            :: dm, pdm
19*49c89c76SBlaise Bourdin    character(len=PETSC_MAX_PATH_LEN)    :: ifilename,ofilename,IOBuffer
20*49c89c76SBlaise Bourdin    integer                              :: order = 1
21*49c89c76SBlaise Bourdin    type(tPetscViewer)                   :: viewer
22*49c89c76SBlaise Bourdin    character(len = MXNAME),dimension(4) :: nodalVarName = ["U_x  ", &
23*49c89c76SBlaise Bourdin                                                            "U_y  ", &
24*49c89c76SBlaise Bourdin                                                            "Alpha", &
25*49c89c76SBlaise Bourdin                                                            "Beta "]
26*49c89c76SBlaise Bourdin    character(len = MXNAME),dimension(3) :: zonalVarName = ["Sigma_11", &
27*49c89c76SBlaise Bourdin                                                            "Sigma_12", &
28*49c89c76SBlaise Bourdin                                                            "Sigma_22"]
29*49c89c76SBlaise Bourdin    character(len = MXNAME)              :: varName
30*49c89c76SBlaise Bourdin
31*49c89c76SBlaise Bourdin    PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER,ierr))
32*49c89c76SBlaise Bourdin    if (ierr /= 0) then
33*49c89c76SBlaise Bourdin      print*,'Unable to initialize PETSc'
34*49c89c76SBlaise Bourdin      stop
35*49c89c76SBlaise Bourdin    endif
36*49c89c76SBlaise Bourdin
37*49c89c76SBlaise Bourdin    PetscCallA(PetscOptionsBegin(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,'PetscViewer_ExodusII test','ex95f90',ierr))
38*49c89c76SBlaise Bourdin    PetscCallA(PetscOptionsString("-i", "Filename to read", "ex95f90", ifilename, ifilename, flg, ierr));
39*49c89c76SBlaise Bourdin    PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing input file name -i <input file name>')
40*49c89c76SBlaise Bourdin    PetscCallA(PetscOptionsString("-o", "Filename to write", "ex95f90", ofilename, ofilename, flg, ierr));
41*49c89c76SBlaise Bourdin    PetscCheckA(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,'missing output file name -o <output file name>')
42*49c89c76SBlaise Bourdin    PetscCallA(PetscOptionsEnd(ierr))
43*49c89c76SBlaise Bourdin
44*49c89c76SBlaise Bourdin    ! Read the mesh in any supported format
45*49c89c76SBlaise Bourdin    PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr))
46*49c89c76SBlaise Bourdin    PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE,ierr));
47*49c89c76SBlaise Bourdin    PetscCallA(DMSetFromOptions(dm,ierr));
48*49c89c76SBlaise Bourdin    PetscCallA(PetscObjectSetName(dm, "ex95f90", ierr));
49*49c89c76SBlaise Bourdin    PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OPTIONS,'-dm_view',ierr));
50*49c89c76SBlaise Bourdin
51*49c89c76SBlaise Bourdin    ! enable exodus debugging information
52*49c89c76SBlaise Bourdin#ifdef PETSC_USE_DEBUG
53*49c89c76SBlaise Bourdin    PetscCallA(exopts(EXVRBS+EXDEBG,ierr))
54*49c89c76SBlaise Bourdin#endif
55*49c89c76SBlaise Bourdin
56*49c89c76SBlaise Bourdin    ! Create the exodus file
57*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,viewer,ierr))
58*49c89c76SBlaise Bourdin
59*49c89c76SBlaise Bourdin    ! Save the geometry to the file, erasing all previous content
60*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIISetOrder(viewer,order,ierr))
61*49c89c76SBlaise Bourdin    PetscCallA(DMView(dm,viewer,ierr))
62*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD,ierr))
63*49c89c76SBlaise Bourdin    PetscCall(PetscViewerFlush(viewer,ierr))
64*49c89c76SBlaise Bourdin
65*49c89c76SBlaise Bourdin    PetscCallA(DMPlexDistribute(dm,0_kPI,PETSC_NULL_SF,pdm,ierr))
66*49c89c76SBlaise Bourdin    if (pdm /= PETSC_NULL_DM) Then
67*49c89c76SBlaise Bourdin      pdm = dm
68*49c89c76SBlaise Bourdin    end if
69*49c89c76SBlaise Bourdin
70*49c89c76SBlaise Bourdin    ! Testing Variable Number
71*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIISetZonalVariable(viewer, nZonalVar, ierr))
72*49c89c76SBlaise Bourdin    nZonalVar = -1
73*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIIGetZonalVariable(viewer, nZonalVar, ierr))
74*49c89c76SBlaise Bourdin    Write(IOBuffer,'("Number of zonal variables:", I2,"\n")') nZonalVar
75*49c89c76SBlaise Bourdin    PetscCallA(PetscPrintf(PETSC_COMM_WORLD,IOBuffer,ierr))
76*49c89c76SBlaise Bourdin
77*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIISetNodalVariable(viewer, nNodalVar, ierr))
78*49c89c76SBlaise Bourdin    nNodalVar = -1
79*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIIGetNodalVariable(viewer, nNodalVar, ierr))
80*49c89c76SBlaise Bourdin    Write(IOBuffer,'("Number of nodal variables:", I2,"\n")') nNodalVar
81*49c89c76SBlaise Bourdin    PetscCallA(PetscPrintf(PETSC_COMM_WORLD,IOBuffer,ierr))
82*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr))
83*49c89c76SBlaise Bourdin
84*49c89c76SBlaise Bourdin    ! Test of PetscViewerExodusIISet[Nodal/Zonal]VariableName
85*49c89c76SBlaise Bourdin    PetscCallA(PetscPrintf(PETSC_COMM_WORLD, "Testing PetscViewerExodusIISet[Nodal/Zonal]VariableName\n", ierr))
86*49c89c76SBlaise Bourdin    do i = 1, nZonalVar
87*49c89c76SBlaise Bourdin        PetscCallA(PetscViewerExodusIISetZonalVariableName(viewer, i-1, zonalVarName(i), ierr))
88*49c89c76SBlaise Bourdin    end do
89*49c89c76SBlaise Bourdin    do i = 1, nNodalVar
90*49c89c76SBlaise Bourdin        PetscCallA(PetscViewerExodusIISetNodalVariableName(viewer, i-1, nodalVarName(i), ierr))
91*49c89c76SBlaise Bourdin    end do
92*49c89c76SBlaise Bourdin    PetscCall(PetscViewerFlush(viewer,ierr))
93*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD, ierr))
94*49c89c76SBlaise Bourdin
95*49c89c76SBlaise Bourdin    do i = 1, nZonalVar
96*49c89c76SBlaise Bourdin        PetscCallA(PetscViewerExodusIIGetZonalVariableName(viewer, i-1, varName, ierr))
97*49c89c76SBlaise Bourdin        Write(IOBuffer,'("   zonal variable:", I2,": ",A,"\n")') i,varName
98*49c89c76SBlaise Bourdin        PetscCallA(PetscPrintf(PETSC_COMM_WORLD,IOBuffer,ierr))
99*49c89c76SBlaise Bourdin    end do
100*49c89c76SBlaise Bourdin    do i = 1, nNodalVar
101*49c89c76SBlaise Bourdin        PetscCallA(PetscViewerExodusIIGetNodalVariableName(viewer, i-1, varName, ierr))
102*49c89c76SBlaise Bourdin        Write(IOBuffer,'("   nodal variable:", I2,": ",A,"\n")') i,varName
103*49c89c76SBlaise Bourdin        PetscCallA(PetscPrintf(PETSC_COMM_WORLD,IOBuffer,ierr))
104*49c89c76SBlaise Bourdin    end do
105*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerDestroy(viewer, ierr))
106*49c89c76SBlaise Bourdin
107*49c89c76SBlaise Bourdin    ! Test of PetscViewerExodusIIGet[Nodal/Zonal]VariableName
108*49c89c76SBlaise Bourdin    nZonalVar = -1
109*49c89c76SBlaise Bourdin    nNodalVar = -1
110*49c89c76SBlaise Bourdin    PetscCallA(PetscPrintf(PETSC_COMM_WORLD, "\n\nReopenning the output file in Read-only mode\n", ierr))
111*49c89c76SBlaise Bourdin    PetscCallA(PetscPrintf(PETSC_COMM_WORLD, "Testing PetscViewerExodusIIGet[Nodal/Zonal]VariableName\n", ierr))
112*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_APPEND, viewer, ierr))
113*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIISetOrder(viewer, order, ierr))
114*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIIGetZonalVariable(viewer, nZonalVar, ierr))
115*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerExodusIIGetNodalVariable(viewer, nNodalVar, ierr))
116*49c89c76SBlaise Bourdin
117*49c89c76SBlaise Bourdin    do i = 1, nZonalVar
118*49c89c76SBlaise Bourdin        PetscCallA(PetscViewerExodusIIGetZonalVariableName(viewer, i-1, varName, ierr))
119*49c89c76SBlaise Bourdin        Write(IOBuffer,'("   zonal variable:", I2,": ",A,"\n")') i,varName
120*49c89c76SBlaise Bourdin        PetscCallA(PetscPrintf(PETSC_COMM_WORLD,IOBuffer,ierr))
121*49c89c76SBlaise Bourdin    end do
122*49c89c76SBlaise Bourdin    do i = 1, nNodalVar
123*49c89c76SBlaise Bourdin        PetscCallA(PetscViewerExodusIIGetNodalVariableName(viewer, i-1, varName, ierr))
124*49c89c76SBlaise Bourdin        Write(IOBuffer,'("   nodal variable:", I2,": ",A,"\n")') i,varName
125*49c89c76SBlaise Bourdin        PetscCallA(PetscPrintf(PETSC_COMM_WORLD,IOBuffer,ierr))
126*49c89c76SBlaise Bourdin    end do
127*49c89c76SBlaise Bourdin
128*49c89c76SBlaise Bourdin    PetscCallA(DMDestroy(dm,ierr))
129*49c89c76SBlaise Bourdin    PetscCallA(PetscViewerDestroy(viewer, ierr))
130*49c89c76SBlaise Bourdin    PetscCallA(PetscFinalize(ierr))
131*49c89c76SBlaise Bourdinend program ex95f90
132*49c89c76SBlaise Bourdin
133*49c89c76SBlaise Bourdin/*TEST
134*49c89c76SBlaise Bourdin
135*49c89c76SBlaise Bourdin  build:
136*49c89c76SBlaise Bourdin    requires: exodusii pnetcdf !complex
137*49c89c76SBlaise Bourdin  test:
138*49c89c76SBlaise Bourdin    suffix: 0
139*49c89c76SBlaise Bourdin    nsize: 1
140*49c89c76SBlaise Bourdin    args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo
141*49c89c76SBlaise Bourdin
142*49c89c76SBlaise BourdinTEST*/
143