MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
ExplicitModel.f90
Go to the documentation of this file.
1 !> @brief Models that solve themselves
3 
4  use kindmodule, only: i4b, dp, lgp
6  use listmodule, only: listtype
8  use basedismodule, only: disbasetype
10 
11  implicit none
12  private
13  public :: explicitmodeltype, &
17 
18  !> @brief Base type for models that solve themselves.
19  !!
20  !! An explicit solution simply scrolls through a list of explicit
21  !! models and calls solution procedures in a prescribed sequence.
22  !<
23  type, extends(basemodeltype) :: explicitmodeltype
24  contains
25  procedure :: model_ad
26  procedure :: model_solve
27  procedure :: model_cq
28  procedure :: model_bd
29  procedure :: model_da
30  procedure :: allocate_scalars
31  procedure :: allocate_arrays
32  procedure :: set_idsoln
33  procedure :: has_pending
34  end type explicitmodeltype
35 
36 contains
37 
38  !> @ brief Advance the model
39  subroutine model_ad(this)
40  class(explicitmodeltype) :: this
41  end subroutine model_ad
42 
43  !> @ brief Solve the model
44  subroutine model_solve(this)
45  class(explicitmodeltype) :: this
46  end subroutine model_solve
47 
48  !> @ brief Calculate model flows
49  subroutine model_cq(this, icnvg, isuppress_output)
50  class(explicitmodeltype) :: this
51  integer(I4B), intent(in) :: icnvg
52  integer(I4B), intent(in) :: isuppress_output
53  end subroutine model_cq
54 
55  !> @ brief Calculate model budget
56  subroutine model_bd(this, icnvg, isuppress_output)
57  class(explicitmodeltype) :: this
58  integer(I4B), intent(in) :: icnvg
59  integer(I4B), intent(in) :: isuppress_output
60  end subroutine model_bd
61 
62  !> @ brief Deallocate the model
63  subroutine model_da(this)
64  class(explicitmodeltype) :: this
65 
66  ! deallocate scalars
67  deallocate (this%filename)
68  call mem_deallocate(this%nja)
69 
70  ! deallocate arrays
71  call mem_deallocate(this%ibound)
72  call mem_deallocate(this%flowja, 'FLOWJA', this%memoryPath)
73 
74  ! nullify pointers
75  if (associated(this%ibound)) &
76  call mem_deallocate(this%ibound, 'IBOUND', this%memoryPath)
77 
78  ! deallocate derived types
79  call this%bndlist%Clear()
80  deallocate (this%bndlist)
81 
82  ! deallocate base type
83  call this%BaseModelType%model_da()
84  end subroutine model_da
85 
86  !> @ brief Allocate scalar variables
87  subroutine allocate_scalars(this, modelname)
88  class(explicitmodeltype) :: this
89  character(len=*), intent(in) :: modelname
90 
91  call this%BaseModelType%allocate_scalars(modelname)
92  allocate (this%bndlist)
93  allocate (this%filename)
94  this%filename = ''
95  end subroutine allocate_scalars
96 
97  !> @brief Allocate array variables
98  subroutine allocate_arrays(this)
99  class(explicitmodeltype) :: this
100  integer(I4B) :: i
101 
102  call mem_allocate(this%nja, 'NJA', this%memoryPath)
103  this%nja = this%dis%nja
104 
105  call mem_allocate(this%ibound, this%dis%nodes, 'IBOUND', this%memoryPath)
106  do i = 1, this%dis%nodes
107  this%ibound(i) = 1 ! active by default
108  end do
109 
110  call mem_allocate(this%flowja, this%nja, 'FLOWJA', this%memoryPath)
111  do i = 1, this%nja
112  this%flowja(i) = dzero
113  end do
114  end subroutine allocate_arrays
115 
116  !> @ brief Cast a generic object into an explicit model
117  function castasexplicitmodelclass(obj) result(res)
118  class(*), pointer, intent(inout) :: obj
119  class(explicitmodeltype), pointer :: res
120 
121  res => null()
122  if (.not. associated(obj)) return
123 
124  select type (obj)
125  class is (explicitmodeltype)
126  res => obj
127  end select
128  end function castasexplicitmodelclass
129 
130  !> @ brief Add explicit model to a generic list
131  subroutine addexplicitmodeltolist(list, model)
132  ! dummy
133  type(listtype), intent(inout) :: list
134  class(explicitmodeltype), pointer, intent(inout) :: model
135  ! local
136  class(*), pointer :: obj
137 
138  obj => model
139  call list%Add(obj)
140  end subroutine addexplicitmodeltolist
141 
142  !> @ brief Get generic object from list and return as explicit model
143  function getexplicitmodelfromlist(list, idx) result(res)
144  ! dummy
145  type(listtype), intent(inout) :: list
146  integer(I4B), intent(in) :: idx
147  class(explicitmodeltype), pointer :: res
148  ! local
149  class(*), pointer :: obj
150 
151  obj => list%GetItem(idx)
152  res => castasexplicitmodelclass(obj)
153  end function getexplicitmodelfromlist
154 
155  !> @brief Set the solution ID
156  subroutine set_idsoln(this, id)
157  class(explicitmodeltype) :: this
158  integer(I4B), intent(in) :: id
159  this%idsoln = id
160  end subroutine set_idsoln
161 
162  !> @brief Whether the model has any pending work.
163  logical function has_pending(this)
164  class(explicitmodeltype), intent(in) :: this
165  has_pending = .false.
166  end function has_pending
167 
168 end module explicitmodelmodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Models that solve themselves.
logical function has_pending(this)
Whether the model has any pending work.
class(explicitmodeltype) function, pointer, public getexplicitmodelfromlist(list, idx)
@ brief Get generic object from list and return as explicit model
subroutine allocate_arrays(this)
Allocate array variables.
subroutine set_idsoln(this, id)
Set the solution ID.
subroutine model_bd(this, icnvg, isuppress_output)
@ brief Calculate model budget
subroutine, public addexplicitmodeltolist(list, model)
@ brief Add explicit model to a generic list
subroutine model_da(this)
@ brief Deallocate the model
subroutine allocate_scalars(this, modelname)
@ brief Allocate scalar variables
subroutine model_solve(this)
@ brief Solve the model
class(explicitmodeltype) function, pointer, public castasexplicitmodelclass(obj)
@ brief Cast a generic object into an explicit model
subroutine model_cq(this, icnvg, isuppress_output)
@ brief Calculate model flows
subroutine model_ad(this)
@ brief Advance the model
This module defines variable data types.
Definition: kind.f90:8
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:16
Base type for models that solve themselves.
A generic heterogeneous doubly-linked list.
Definition: List.f90:14