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