MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
observemodule Module Reference

This module contains the derived types ObserveType and ObsDataType. More...

Data Types

type  observetype
 
type  obsdatatype
 
interface  ProcessIdSub
 @ brief Process user-provided IDstring More...
 

Functions/Subroutines

subroutine resetcurrentvalue (this)
 @ brief Reset current observation value More...
 
subroutine writeto (this, obstab, btagfound, fnamein)
 @ brief Write observation input data More...
 
subroutine resetobsindex (this)
 @ brief Reset a observation index More...
 
subroutine addobsindex (this, indx)
 @ brief Add a observation index More...
 
subroutine da (this)
 @ brief Deallocate a observation More...
 
subroutine, public constructobservation (newObservation, defLine, numunit, formatted, indx, obsData, inunit)
 @ brief Construct a new ObserveType More...
 
type(observetype) function, pointer castasobservetype (obj)
 @ brief Cast a object as a ObserveType More...
 
subroutine, public addobstolist (list, obs)
 @ brief Add a ObserveType to a list More...
 
type(observetype) function, pointer, public getobsfromlist (list, idx)
 @ brief Get an ObserveType from a list More...
 

Detailed Description

This module contains the derived types ObserveType and ObsDataType.

  • ObserveType – is designed to contain all information and functionality needed for one observation. ObserveType contains a pointer to an ObsDataType object.
  • ObsDataType – is for storing package ID, observation type, and a pointer to a subroutine that will be called to process the IDstring provided in Obs input. The ProcessIdPtr member of ObsDataType requires a pointer to an ObserveType object.

Function/Subroutine Documentation

◆ addobsindex()

subroutine observemodule::addobsindex ( class(observetype), intent(inout)  this,
integer(i4b), intent(in)  indx 
)
private

Subroutine to add the observation index to the observation index array (indxbnds). The observation index count (indxbnds_count) is also incremented by one and the observation index array is expanded, if necessary.

Parameters
[in]indxobservation index

Definition at line 192 of file Observe.f90.

193  ! -- dummy
194  class(ObserveType), intent(inout) :: this
195  integer(I4B), intent(in) :: indx !< observation index
196  !
197  ! -- Increment the index count
198  this%indxbnds_count = this%indxbnds_count + 1
199  !
200  ! -- Expand the observation index array, if necessary
201  call expandarraywrapper(this%indxbnds_count, this%indxbnds, loginc=.true.)
202  !
203  ! -- add index to observation index
204  this%indxbnds(this%indxbnds_count) = indx

◆ addobstolist()

subroutine, public observemodule::addobstolist ( type(listtype), intent(inout)  list,
type(observetype), intent(inout), pointer  obs 
)

Subroutine to add a ObserveType to a list.

Parameters
[in,out]listObserveType list
[in,out]obsObserveType

Definition at line 318 of file Observe.f90.

319  ! -- dummy
320  type(ListType), intent(inout) :: list !< ObserveType list
321  type(ObserveType), pointer, intent(inout) :: obs !< ObserveType
322  ! -- local
323  class(*), pointer :: obj
324  !
325  obj => obs
326  call list%Add(obj)
Here is the caller graph for this function:

◆ castasobservetype()

type(observetype) function, pointer observemodule::castasobservetype ( class(*), intent(inout), pointer  obj)
private

Function to cast an object as a ObserveType object.

Parameters
[in,out]objobject
Returns
returned ObserveType object

Definition at line 298 of file Observe.f90.

299  ! -- dummy
300  class(*), pointer, intent(inout) :: obj !< object
301  ! -- return
302  type(ObserveType), pointer :: res !< returned ObserveType object
303  !
304  res => null()
305  if (.not. associated(obj)) return
306  !
307  select type (obj)
308  type is (observetype)
309  res => obj
310  end select
Here is the caller graph for this function:

◆ constructobservation()

subroutine, public observemodule::constructobservation ( type(observetype), pointer  newObservation,
character(len=*), intent(in)  defLine,
integer(i4b), intent(in)  numunit,
logical, intent(in)  formatted,
integer(i4b), intent(in)  indx,
type(obsdatatype), dimension(:), intent(in), pointer  obsData,
integer(i4b), intent(in)  inunit 
)

Subroutine to construct and return an ObserveType object based on the contents of defLine.

Parameters
newobservationnew ObserveType
[in]deflinestring with observation data
[in]numunitOutput unit number
[in]formattedlogical indicating if formatted output will be written
[in]indxIndex in ObsOutput array
[in]obsdataobsData type
[in]inunitobservation input file unit

Definition at line 228 of file Observe.f90.

230  ! -- dummy variables
231  type(ObserveType), pointer :: newObservation !< new ObserveType
232  character(len=*), intent(in) :: defLine !< string with observation data
233  integer(I4B), intent(in) :: numunit !< Output unit number
234  logical, intent(in) :: formatted !< logical indicating if formatted output will be written
235  integer(I4B), intent(in) :: indx !< Index in ObsOutput array
236  type(ObsDataType), dimension(:), pointer, intent(in) :: obsData !< obsData type
237  integer(I4B), intent(in) :: inunit !< observation input file unit
238  ! -- local
239  real(DP) :: r
240  integer(I4B) :: i
241  integer(I4B) :: icol
242  integer(I4B) :: iout
243  integer(I4B) :: istart
244  integer(I4B) :: istop
245  integer(I4B) :: n
246  !
247  ! -- initialize
248  iout = 0
249  icol = 1
250  !
251  ! -- Allocate an ObserveType object.
252  allocate (newobservation)
253  allocate (newobservation%indxbnds(0))
254  !
255  ! -- Set indxbnds_count to 0
256  newobservation%indxbnds_count = 0
257  !
258  ! -- Define the contents of the ObservationSingleType object based on the
259  ! contents of defLine.
260  !
261  ! -- Get observation name and store it
262  call urword(defline, icol, istart, istop, 1, n, r, iout, inunit)
263  newobservation%Name = defline(istart:istop)
264  !
265  ! -- Get observation type, convert it to uppercase, and store it.
266  call urword(defline, icol, istart, istop, 1, n, r, iout, inunit)
267  newobservation%ObsTypeId = defline(istart:istop)
268  !
269  ! -- Look up package ID for this observation type and store it
270  do i = 1, maxobstypes
271  if (obsdata(i)%ObsTypeID == newobservation%ObsTypeId) then
272  newobservation%obsDatum => obsdata(i)
273  exit
274  elseif (obsdata(i)%ObsTypeID == '') then
275  exit
276  end if
277  end do
278  !
279  ! -- Remaining text is ID [and ID2]; store the remainder of the string
280  istart = istop + 1
281  istop = len_trim(defline)
282  if (istart > istop) then
283  istart = istop
284  end if
285  newobservation%IDstring = defline(istart:istop)
286  !
287  ! Store UnitNumber, FormattedOutput, and IndxObsOutput
288  newobservation%UnitNumber = numunit
289  newobservation%FormattedOutput = formatted
290  newobservation%IndxObsOutput = indx
Here is the call graph for this function:
Here is the caller graph for this function:

◆ da()

subroutine observemodule::da ( class(observetype), intent(inout)  this)
private

Subroutine to deallocated a observation (ObserveType).

Definition at line 212 of file Observe.f90.

213  ! -- dummy
214  class(ObserveType), intent(inout) :: this
215  if (allocated(this%indxbnds)) then
216  deallocate (this%indxbnds)
217  end if

◆ getobsfromlist()

type(observetype) function, pointer, public observemodule::getobsfromlist ( type(listtype), intent(inout)  list,
integer(i4b), intent(in)  idx 
)

Function to get an ObserveType from a list.

Parameters
[in,out]listObserveType list
[in]idxObserveType list index
Returns
returned ObserveType

Definition at line 334 of file Observe.f90.

335  ! -- dummy
336  type(ListType), intent(inout) :: list !< ObserveType list
337  integer(I4B), intent(in) :: idx !< ObserveType list index
338  ! -- return
339  type(ObserveType), pointer :: res !< returned ObserveType
340  ! -- local
341  class(*), pointer :: obj
342  !
343  obj => list%GetItem(idx)
344  res => castasobservetype(obj)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ resetcurrentvalue()

subroutine observemodule::resetcurrentvalue ( class(observetype), intent(inout)  this)

Subroutine to reset the current observation value.

Definition at line 117 of file Observe.f90.

118  ! -- dummy
119  class(ObserveType), intent(inout) :: this
120  !
121  ! -- Reset current value to zero.
122  this%CurrentTimeStepEndValue = dzero

◆ resetobsindex()

subroutine observemodule::resetobsindex ( class(observetype), intent(inout)  this)
private

Subroutine to reset the observation index count and array.

Definition at line 168 of file Observe.f90.

169  ! -- dummy
170  class(ObserveType), intent(inout) :: this
171  !
172  ! -- Reset the index count
173  this%indxbnds_count = 0
174  !
175  ! -- Deallocate observation index array, if necessary
176  if (allocated(this%indxbnds)) then
177  deallocate (this%indxbnds)
178  end if
179  !
180  ! -- Allocate observation index array to size 0
181  allocate (this%indxbnds(0))

◆ writeto()

subroutine observemodule::writeto ( class(observetype), intent(inout)  this,
type(tabletype), intent(inout)  obstab,
character(len=*), intent(in)  btagfound,
character(len=*), intent(in)  fnamein 
)
private

Subroutine to write observation input data to a table in the model list file.

Parameters
[in,out]obstabobservation table
[in]btagfoundlogical indicating if boundname was found
[in]fnameinobservation input file name

Definition at line 131 of file Observe.f90.

132  ! -- dummy
133  class(ObserveType), intent(inout) :: this
134  type(TableType), intent(inout) :: obstab !< observation table
135  character(len=*), intent(in) :: btagfound !< logical indicating if boundname was found
136  character(len=*), intent(in) :: fnamein !< observation input file name
137  ! -- local
138  character(len=12) :: tag
139  character(len=80) :: fnameout
140  !
141  ! -- write btagfound to tag
142  if (len_trim(btagfound) > 12) then
143  tag = btagfound(1:12)
144  else
145  write (tag, '(a12)') btagfound
146  end if
147  !
148  ! -- write fnamein to fnameout
149  if (len_trim(fnamein) > 80) then
150  fnameout = fnamein(1:80)
151  else
152  write (fnameout, '(a80)') fnamein
153  end if
154  !
155  ! -- write data to observation table
156  call obstab%add_term(this%Name)
157  call obstab%add_term(tag//trim(this%ObsTypeId))
158  call obstab%add_term('ALL TIMES')
159  call obstab%add_term('"'//trim(this%IDstring)//'"')
160  call obstab%add_term(fnameout)