MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
BaseModel.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
5  use listmodule, only: listtype
6  implicit none
7 
8  private
11 
12  !> @brief Highest level model type. All models extend this parent type.
13  type :: basemodeltype
14  character(len=LENMEMPATH) :: memorypath !< the location in the memory manager where the variables are stored
15  character(len=LENMODELNAME), pointer :: name => null() !< name of the model
16  character(len=3), pointer :: macronym => null() !< 3 letter model acronym (GWF, GWT, ...)
17  integer(I4B), pointer :: idsoln => null() !< id of the solution model is in
18  integer(I4B), pointer :: id => null() !< model id
19  integer(I4B), pointer :: iout => null() !< output unit number
20  integer(I4B), pointer :: inewton => null() !< newton-raphson flag
21  integer(I4B), pointer :: iprpak => null() !< integer flag to echo input
22  integer(I4B), pointer :: iprflow => null() !< flag to print simulated flows
23  integer(I4B), pointer :: ipakcb => null() !< save_flows flag
24  contains
25  procedure :: model_df
26  procedure :: model_ar
27  procedure :: model_rp
28  procedure :: model_dt
29  procedure :: model_ot
30  procedure :: model_fp
31  procedure :: model_da
32  procedure :: allocate_scalars
33  procedure :: model_message
34  end type basemodeltype
35 
36 contains
37 
38  !> @brief Define the model
39  !<
40  subroutine model_df(this)
41  class(basemodeltype) :: this
42  end subroutine model_df
43 
44  !> @brief Allocate and read
45  !<
46  subroutine model_ar(this)
47  class(basemodeltype) :: this
48  end subroutine model_ar
49 
50  !> @brief Read and prepare
51  !<
52  subroutine model_rp(this)
53  class(basemodeltype) :: this
54  end subroutine model_rp
55 
56  !> @brief Calculate time step length
57  !<
58  subroutine model_dt(this)
59  class(basemodeltype) :: this
60  end subroutine model_dt
61 
62  !> @brief Output results
63  !<
64  subroutine model_ot(this)
65  class(basemodeltype) :: this
66  end subroutine model_ot
67 
68  !> @brief Write line to model iout
69  !<
70  subroutine model_message(this, line, fmt)
71  ! -- dummy
72  class(basemodeltype) :: this
73  character(len=*), intent(in) :: line
74  character(len=*), intent(in), optional :: fmt
75  ! -- local
76  character(len=LINELENGTH) :: cfmt
77  !
78  ! -- process optional variables
79  if (present(fmt)) then
80  cfmt = fmt
81  else
82  cfmt = '(1x,a)'
83  end if
84  !
85  ! -- write line
86  write (this%iout, trim(cfmt)) trim(line)
87  end subroutine model_message
88 
89  !> @brief Final processing
90  !<
91  subroutine model_fp(this)
92  class(basemodeltype) :: this
93  end subroutine model_fp
94 
95  !> @brief Allocate scalar variables
96  !<
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  !<
127  subroutine model_da(this)
128  ! -- modules
130  ! -- dummy
131  class(basemodeltype) :: this
132  !
133  ! -- Strings
134  call mem_deallocate(this%name, 'NAME', this%memoryPath)
135  call mem_deallocate(this%macronym, 'MACRONYM', this%memoryPath)
136  !
137  ! -- Scalars
138  call mem_deallocate(this%id)
139  call mem_deallocate(this%iout)
140  call mem_deallocate(this%inewton)
141  call mem_deallocate(this%iprpak)
142  call mem_deallocate(this%iprflow)
143  call mem_deallocate(this%ipakcb)
144  call mem_deallocate(this%idsoln)
145  end subroutine model_da
146 
147  function castasbasemodelclass(obj) result(res)
148  class(*), pointer, intent(inout) :: obj
149  class(basemodeltype), pointer :: res
150  !
151  res => null()
152  if (.not. associated(obj)) return
153  !
154  select type (obj)
155  class is (basemodeltype)
156  res => obj
157  end select
158  end function castasbasemodelclass
159 
160  subroutine addbasemodeltolist(list, model)
161  ! -- dummy
162  type(listtype), intent(inout) :: list
163  class(basemodeltype), pointer, intent(inout) :: model
164  ! -- local
165  class(*), pointer :: obj
166  !
167  obj => model
168  call list%Add(obj)
169  end subroutine addbasemodeltolist
170 
171  function getbasemodelfromlist(list, idx) result(res)
172  ! -- dummy
173  type(listtype), intent(inout) :: list
174  integer(I4B), intent(in) :: idx
175  class(basemodeltype), pointer :: res
176  ! -- local
177  class(*), pointer :: obj
178  !
179  obj => list%GetItem(idx)
180  res => castasbasemodelclass(obj)
181  end function getbasemodelfromlist
182 
183 end module basemodelmodule
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:161
subroutine model_ar(this)
Allocate and read.
Definition: BaseModel.f90:47
subroutine model_ot(this)
Output results.
Definition: BaseModel.f90:65
subroutine model_message(this, line, fmt)
Write line to model iout.
Definition: BaseModel.f90:71
class(basemodeltype) function, pointer, public castasbasemodelclass(obj)
Definition: BaseModel.f90:148
class(basemodeltype) function, pointer, public getbasemodelfromlist(list, idx)
Definition: BaseModel.f90:172
subroutine model_dt(this)
Calculate time step length.
Definition: BaseModel.f90:59
subroutine allocate_scalars(this, modelname)
Allocate scalar variables.
Definition: BaseModel.f90:98
subroutine model_rp(this)
Read and prepare.
Definition: BaseModel.f90:53
subroutine model_df(this)
Define the model.
Definition: BaseModel.f90:41
subroutine model_fp(this)
Final processing.
Definition: BaseModel.f90:92
subroutine model_da(this)
Deallocate.
Definition: BaseModel.f90:128
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
This module defines variable data types.
Definition: kind.f90:8
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
A generic heterogeneous doubly-linked list.
Definition: List.f90:14