MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
prt-obs.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
5  use basedismodule, only: disbasetype
6  use observemodule, only: observetype
7  use obsmodule, only: obstype
8  use simmodule, only: count_errors, store_error, &
10  implicit none
11 
12  private
13  public :: prtobstype, prt_obs_cr
14 
15  type, extends(obstype) :: prtobstype
16  ! -- Private members
17  real(dp), dimension(:), pointer, contiguous, private :: x => null() !< concentration
18  real(dp), dimension(:), pointer, contiguous, private :: flowja => null() !< intercell flows
19  contains
20  ! -- Public procedures
21  procedure, public :: prt_obs_ar
22  procedure, public :: obs_bd => prt_obs_bd
23  procedure, public :: obs_df => prt_obs_df
24  procedure, public :: obs_rp => prt_obs_rp
25  procedure, public :: obs_da => prt_obs_da
26  ! -- Private procedures
27  procedure, private :: set_pointers
28  end type prtobstype
29 
30 contains
31 
32  !> @brief Create a new PrtObsType object
33  subroutine prt_obs_cr(obs, inobs)
34  ! -- dummy
35  type(prtobstype), pointer, intent(out) :: obs
36  integer(I4B), pointer, intent(in) :: inobs
37  !
38  allocate (obs)
39  call obs%allocate_scalars()
40  obs%active = .false.
41  obs%inputFilename = ''
42  obs%inUnitObs => inobs
43 
44  end subroutine prt_obs_cr
45 
46  !> @brief Allocate and read
47  subroutine prt_obs_ar(this, x, flowja)
48  ! -- dummy
49  class(prtobstype), intent(inout) :: this
50  real(DP), dimension(:), pointer, contiguous, intent(in) :: x
51  real(DP), dimension(:), pointer, contiguous, intent(in) :: flowja
52  !
53  ! Call ar method of parent class
54  call this%obs_ar()
55  !
56  ! set pointers
57  call this%set_pointers(x, flowja)
58 
59  end subroutine prt_obs_ar
60 
61  !> @brief Define package
62  subroutine prt_obs_df(this, iout, pkgname, filtyp, dis)
63  ! -- dummy
64  class(prtobstype), intent(inout) :: this
65  integer(I4B), intent(in) :: iout
66  character(len=*), intent(in) :: pkgname
67  character(len=*), intent(in) :: filtyp
68  class(disbasetype), pointer :: dis
69  ! -- local
70  integer(I4B) :: indx
71  !
72  ! Call overridden method of parent class
73  call this%ObsType%obs_df(iout, pkgname, filtyp, dis)
74  !
75  ! -- StoreObsType arguments are: (ObserveType, cumulative, indx);
76  ! indx is returned.
77  !
78  ! -- Store obs type and assign procedure pointer for head observation type
79  call this%StoreObsType('concentration', .false., indx)
80  this%obsData(indx)%ProcessIdPtr => prt_process_concentration_obs_id
81  !
82  ! -- Store obs type and assign procedure pointer for flow-ja-face observation type
83  call this%StoreObsType('flow-ja-face', .true., indx)
84  this%obsData(indx)%ProcessIdPtr => prt_process_intercell_obs_id
85 
86  end subroutine prt_obs_df
87 
88  !> @brief Save observations
89  subroutine prt_obs_bd(this)
90  ! -- dummy
91  class(prtobstype), intent(inout) :: this
92  ! -- local
93  integer(I4B) :: i, jaindex, nodenumber
94  character(len=100) :: msg
95  class(observetype), pointer :: obsrv => null()
96  !
97  call this%obs_bd_clear()
98  !
99  ! -- iterate through all PRT observations
100  if (this%npakobs > 0) then
101  do i = 1, this%npakobs
102  obsrv => this%pakobs(i)%obsrv
103  nodenumber = obsrv%NodeNumber
104  jaindex = obsrv%JaIndex
105  select case (obsrv%ObsTypeId)
106  case ('CONCENTRATION')
107  call this%SaveOneSimval(obsrv, this%x(nodenumber))
108  case ('FLOW-JA-FACE')
109  call this%SaveOneSimval(obsrv, this%flowja(jaindex))
110  case default
111  msg = 'Error: Unrecognized observation type: '//trim(obsrv%ObsTypeId)
112  call store_error(msg)
113  call store_error_unit(this%inUnitObs)
114  end select
115  end do
116  end if
117 
118  end subroutine prt_obs_bd
119 
120  !> @brief Read and prepare
121  subroutine prt_obs_rp(this)
122  class(prtobstype), intent(inout) :: this
123  ! Do PRT observations need any checking? If so, add checks here
124  end subroutine prt_obs_rp
125 
126  !> @brief Deallocate
127  subroutine prt_obs_da(this)
128  ! -- dummy
129  class(prtobstype), intent(inout) :: this
130  !
131  nullify (this%x)
132  nullify (this%flowja)
133  call this%ObsType%obs_da()
134 
135  end subroutine prt_obs_da
136 
137  !> @brief Set pointers
138  subroutine set_pointers(this, x, flowja)
139  ! -- dummy
140  class(prtobstype), intent(inout) :: this
141  real(DP), dimension(:), pointer, contiguous, intent(in) :: x
142  real(DP), dimension(:), pointer, contiguous, intent(in) :: flowja
143  !
144  this%x => x
145  this%flowja => flowja
146 
147  end subroutine set_pointers
148 
149  ! -- Procedures related to GWF observations (NOT type-bound)
150 
151  subroutine prt_process_concentration_obs_id(obsrv, dis, inunitobs, iout)
152  ! -- dummy
153  type(observetype), intent(inout) :: obsrv
154  class(disbasetype), intent(in) :: dis
155  integer(I4B), intent(in) :: inunitobs
156  integer(I4B), intent(in) :: iout
157  ! -- local
158  integer(I4B) :: nn1
159  integer(I4B) :: icol, istart, istop
160  character(len=LINELENGTH) :: ermsg, string
161  !
162  ! -- Initialize variables
163  string = obsrv%IDstring
164  icol = 1
165  !
166  ! Get node number, with option for ID string to be either node
167  ! number or lay, row, column (when dis is structured).
168  nn1 = dis%noder_from_string(icol, istart, istop, inunitobs, &
169  iout, string, .false.)
170  !
171  if (nn1 > 0) then
172  obsrv%NodeNumber = nn1
173  else
174  ermsg = 'Error reading data from ID string'
175  call store_error(ermsg)
176  call store_error_unit(inunitobs)
177  end if
178 
179  end subroutine prt_process_concentration_obs_id
180 
181  subroutine prt_process_intercell_obs_id(obsrv, dis, inunitobs, iout)
182  ! -- dummy
183  type(observetype), intent(inout) :: obsrv
184  class(disbasetype), intent(in) :: dis
185  integer(I4B), intent(in) :: inunitobs
186  integer(I4B), intent(in) :: iout
187  ! -- local
188  integer(I4B) :: nn1, nn2
189  integer(I4B) :: icol, istart, istop, jaidx
190  character(len=LINELENGTH) :: ermsg, string
191  ! formats
192 70 format('Error: No connection exists between cells identified in text: ', a)
193  !
194  ! -- Initialize variables
195  string = obsrv%IDstring
196  icol = 1
197  !
198  ! Get node number, with option for ID string to be either node
199  ! number or lay, row, column (when dis is structured).
200  nn1 = dis%noder_from_string(icol, istart, istop, inunitobs, &
201  iout, string, .false.)
202  !
203  if (nn1 > 0) then
204  obsrv%NodeNumber = nn1
205  else
206  ermsg = 'Error reading data from ID string: '//string(istart:istop)
207  call store_error(ermsg)
208  end if
209  !
210  ! Get node number, with option for ID string to be either node
211  ! number or lay, row, column (when dis is structured).
212  nn2 = dis%noder_from_string(icol, istart, istop, inunitobs, &
213  iout, string, .false.)
214  if (nn2 > 0) then
215  obsrv%NodeNumber2 = nn2
216  else
217  ermsg = 'Error reading data from ID string: '//string(istart:istop)
218  call store_error(ermsg)
219  end if
220  !
221  ! -- store JA index
222  jaidx = dis%con%getjaindex(nn1, nn2)
223  if (jaidx == 0) then
224  write (ermsg, 70) trim(string)
225  call store_error(ermsg)
226  end if
227  obsrv%JaIndex = jaidx
228  !
229  if (count_errors() > 0) then
230  call store_error_unit(inunitobs)
231  end if
232 
233  end subroutine prt_process_intercell_obs_id
234 
235 end module prtobsmodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
integer(i4b), parameter maxobstypes
maximum number of observation types
Definition: Constants.f90:47
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine prt_obs_df(this, iout, pkgname, filtyp, dis)
Define package.
Definition: prt-obs.f90:63
subroutine prt_process_intercell_obs_id(obsrv, dis, inunitobs, iout)
Definition: prt-obs.f90:182
subroutine prt_process_concentration_obs_id(obsrv, dis, inunitobs, iout)
Definition: prt-obs.f90:152
subroutine set_pointers(this, x, flowja)
Set pointers.
Definition: prt-obs.f90:139
subroutine prt_obs_da(this)
Deallocate.
Definition: prt-obs.f90:128
subroutine prt_obs_rp(this)
Read and prepare.
Definition: prt-obs.f90:122
subroutine prt_obs_ar(this, x, flowja)
Allocate and read.
Definition: prt-obs.f90:48
subroutine prt_obs_bd(this)
Save observations.
Definition: prt-obs.f90:90
subroutine, public prt_obs_cr(obs, inobs)
Create a new PrtObsType object.
Definition: prt-obs.f90:34
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168