MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
BaseModel.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
5  use listmodule, only: listtype
6  use basedismodule, only: disbasetype
9  implicit none
10 
11  private
14 
15  !> @brief Highest level model type. All models extend this parent type.
16  type :: basemodeltype
17  character(len=LENMEMPATH) :: memorypath !< the location in the memory manager where the variables are stored
18  character(len=LENMODELNAME), pointer :: name => null() !< name of the model
19  character(len=LINELENGTH), pointer :: filename => null() !< input file name
20  character(len=3), pointer :: macronym => null() !< 3 letter model acronym (GWF, GWT, ...)
21  integer(I4B), pointer :: idsoln => null() !< id of the solution model is in
22  integer(I4B), pointer :: id => null() !< model id
23  integer(I4B), pointer :: iout => null() !< output unit number
24  integer(I4B), pointer :: inewton => null() !< newton-raphson flag
25  integer(I4B), pointer :: iprpak => null() !< integer flag to echo input
26  integer(I4B), pointer :: iprflow => null() !< flag to print simulated flows
27  integer(I4B), pointer :: ipakcb => null() !< save_flows flag
28  integer(I4B), pointer :: nja => null() !< number of connections
29  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< ibound
30  real(dp), dimension(:), pointer, contiguous :: flowja => null() !< intercell flows
31  type(listtype), pointer :: bndlist => null() !< array of boundary packages
32  class(disbasetype), pointer :: dis => null() !< discretization object
33  contains
34  procedure :: model_df
35  procedure :: model_ar
36  procedure :: model_rp
37  procedure :: model_dt
38  procedure :: model_ot
39  procedure :: model_fp
40  procedure :: model_da
41  procedure :: allocate_scalars
42  procedure :: model_message
43  procedure :: create_lstfile
44  end type basemodeltype
45 
46 contains
47 
48  !> @brief Define the model
49  subroutine model_df(this)
50  class(basemodeltype) :: this
51  end subroutine model_df
52 
53  !> @brief Allocate and read
54  subroutine model_ar(this)
55  class(basemodeltype) :: this
56  end subroutine model_ar
57 
58  !> @brief Read and prepare
59  subroutine model_rp(this)
60  class(basemodeltype) :: this
61  end subroutine model_rp
62 
63  !> @brief Calculate time step length
64  subroutine model_dt(this)
65  class(basemodeltype) :: this
66  end subroutine model_dt
67 
68  !> @brief Output results
69  subroutine model_ot(this)
70  class(basemodeltype) :: this
71  end subroutine model_ot
72 
73  !> @brief Write line to model iout
74  subroutine model_message(this, line, fmt)
75  ! dummy
76  class(basemodeltype) :: this
77  character(len=*), intent(in) :: line
78  character(len=*), intent(in), optional :: fmt
79  ! local
80  character(len=LINELENGTH) :: cfmt
81 
82  if (present(fmt)) then
83  cfmt = fmt
84  else
85  cfmt = '(1x,a)'
86  end if
87 
88  write (this%iout, trim(cfmt)) trim(line)
89  end subroutine model_message
90 
91  !> @brief Final processing
92  subroutine model_fp(this)
93  class(basemodeltype) :: this
94  end subroutine model_fp
95 
96  !> @brief Allocate scalar variables
97  subroutine allocate_scalars(this, modelname)
98  ! modules
100  ! dummy
101  class(basemodeltype) :: this
102  character(len=*), intent(in) :: modelname
103 
104  call mem_allocate(this%name, lenmodelname, 'NAME', this%memoryPath)
105  call mem_allocate(this%macronym, 3, 'MACRONYM', this%memoryPath)
106  call mem_allocate(this%id, 'ID', this%memoryPath)
107  call mem_allocate(this%iout, 'IOUT', this%memoryPath)
108  call mem_allocate(this%inewton, 'INEWTON', this%memoryPath)
109  call mem_allocate(this%iprpak, 'IPRPAK', this%memoryPath)
110  call mem_allocate(this%iprflow, 'IPRFLOW', this%memoryPath)
111  call mem_allocate(this%ipakcb, 'IPAKCB', this%memoryPath)
112  call mem_allocate(this%idsoln, 'IDSOLN', this%memoryPath)
113 
114  this%name = modelname
115  this%macronym = ''
116  this%idsoln = 0
117  this%id = 0
118  this%iout = 0
119  this%iprpak = 0
120  this%iprflow = 0
121  this%ipakcb = 0
122  this%inewton = 0 !default is standard formulation
123  end subroutine allocate_scalars
124 
125  !> @brief Deallocate
126  subroutine model_da(this)
127  ! modules
129  ! dummy
130  class(basemodeltype) :: this
131 
132  ! deallocate strings
133  call mem_deallocate(this%name, 'NAME', this%memoryPath)
134  call mem_deallocate(this%macronym, 'MACRONYM', this%memoryPath)
135  !
136  ! deallocate scalars
137  call mem_deallocate(this%id)
138  call mem_deallocate(this%iout)
139  call mem_deallocate(this%inewton)
140  call mem_deallocate(this%iprpak)
141  call mem_deallocate(this%iprflow)
142  call mem_deallocate(this%ipakcb)
143  call mem_deallocate(this%idsoln)
144  end subroutine model_da
145 
146  function castasbasemodelclass(obj) result(res)
147  class(*), pointer, intent(inout) :: obj
148  class(basemodeltype), pointer :: res
149 
150  res => null()
151  if (.not. associated(obj)) return
152 
153  select type (obj)
154  class is (basemodeltype)
155  res => obj
156  end select
157  end function castasbasemodelclass
158 
159  subroutine addbasemodeltolist(list, model)
160  ! dummy
161  type(listtype), intent(inout) :: list
162  class(basemodeltype), pointer, intent(inout) :: model
163  ! local
164  class(*), pointer :: obj
165 
166  obj => model
167  call list%Add(obj)
168  end subroutine addbasemodeltolist
169 
170  function getbasemodelfromlist(list, idx) result(res)
171  ! dummy
172  type(listtype), intent(inout) :: list
173  integer(I4B), intent(in) :: idx
174  class(basemodeltype), pointer :: res
175  ! local
176  class(*), pointer :: obj
177 
178  obj => list%GetItem(idx)
179  res => castasbasemodelclass(obj)
180  end function getbasemodelfromlist
181 
182  subroutine create_lstfile(this, lst_fname, model_fname, defined, headertxt)
183  ! dummy
184  class(basemodeltype) :: this
185  character(len=*), intent(inout) :: lst_fname
186  character(len=*), intent(in) :: model_fname
187  logical(LGP), intent(in) :: defined
188  character(len=*), intent(in) :: headertxt
189  ! local
190  integer(I4B) :: i, istart, istop
191 
192  ! set list file name if not provided
193  if (.not. defined) then
194 
195  ! initialize
196  lst_fname = ' '
197  istart = 0
198  istop = len_trim(model_fname)
199 
200  ! identify '.' character position from back of string
201  do i = istop, 1, -1
202  if (model_fname(i:i) == '.') then
203  istart = i
204  exit
205  end if
206  end do
207 
208  ! if not found start from string end
209  if (istart == 0) istart = istop + 1
210 
211  ! set list file name
212  lst_fname = model_fname(1:istart)
213  istop = istart + 3
214  lst_fname(istart:istop) = '.lst'
215  end if
216 
217  ! create the list file
218  this%iout = getunit()
219  call openfile(this%iout, 0, lst_fname, 'LIST', filstat_opt='REPLACE')
220 
221  ! write list file header
222  call write_listfile_header(this%iout, headertxt)
223  end subroutine create_lstfile
224 
225 end module basemodelmodule
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:160
subroutine model_ar(this)
Allocate and read.
Definition: BaseModel.f90:55
subroutine model_ot(this)
Output results.
Definition: BaseModel.f90:70
subroutine model_message(this, line, fmt)
Write line to model iout.
Definition: BaseModel.f90:75
class(basemodeltype) function, pointer, public castasbasemodelclass(obj)
Definition: BaseModel.f90:147
class(basemodeltype) function, pointer, public getbasemodelfromlist(list, idx)
Definition: BaseModel.f90:171
subroutine model_dt(this)
Calculate time step length.
Definition: BaseModel.f90:65
subroutine allocate_scalars(this, modelname)
Allocate scalar variables.
Definition: BaseModel.f90:98
subroutine model_rp(this)
Read and prepare.
Definition: BaseModel.f90:60
subroutine model_df(this)
Define the model.
Definition: BaseModel.f90:50
subroutine model_fp(this)
Final processing.
Definition: BaseModel.f90:93
subroutine create_lstfile(this, lst_fname, model_fname, defined, headertxt)
Definition: BaseModel.f90:183
subroutine model_da(this)
Deallocate.
Definition: BaseModel.f90:127
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 lenmodelname
maximum length of the model name
Definition: Constants.f90:22
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
This module contains version information.
Definition: version.f90:7
subroutine write_listfile_header(iout, cmodel_type, write_sys_command, write_kind_info)
@ brief Write program header
Definition: version.f90:98
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:16
A generic heterogeneous doubly-linked list.
Definition: List.f90:14