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  ! Overridden methods
26  procedure :: model_ad
27  procedure :: model_solve
28  procedure :: model_cq
29  procedure :: model_bd
30  procedure :: model_da
31  ! Utility methods
32  procedure :: allocate_scalars
33  procedure :: allocate_arrays
34  procedure :: set_idsoln
35  end type explicitmodeltype
36 
37 contains
38 
39  !> @ brief Advance the model
40  subroutine model_ad(this)
41  class(explicitmodeltype) :: this
42  end subroutine model_ad
43 
44  !> @ brief Solve the model
45  subroutine model_solve(this)
46  class(explicitmodeltype) :: this
47  end subroutine model_solve
48 
49  !> @ brief Calculate model flows
50  subroutine model_cq(this, icnvg, isuppress_output)
51  class(explicitmodeltype) :: this
52  integer(I4B), intent(in) :: icnvg
53  integer(I4B), intent(in) :: isuppress_output
54  end subroutine model_cq
55 
56  !> @ brief Calculate model budget
57  subroutine model_bd(this, icnvg, isuppress_output)
58  class(explicitmodeltype) :: this
59  integer(I4B), intent(in) :: icnvg
60  integer(I4B), intent(in) :: isuppress_output
61  end subroutine model_bd
62 
63  !> @ brief Deallocate the model
64  subroutine model_da(this)
65  class(explicitmodeltype) :: this
66 
67  ! deallocate scalars
68  deallocate (this%filename)
69  call mem_deallocate(this%nja)
70 
71  ! deallocate arrays
72  call mem_deallocate(this%ibound)
73  call mem_deallocate(this%flowja, 'FLOWJA', this%memoryPath)
74 
75  ! nullify pointers
76  if (associated(this%ibound)) &
77  call mem_deallocate(this%ibound, 'IBOUND', this%memoryPath)
78 
79  ! deallocate derived types
80  call this%bndlist%Clear()
81  deallocate (this%bndlist)
82 
83  ! deallocate base type
84  call this%BaseModelType%model_da()
85  end subroutine model_da
86 
87  !> @ brief Allocate scalar variables
88  subroutine allocate_scalars(this, modelname)
89  class(explicitmodeltype) :: this
90  character(len=*), intent(in) :: modelname
91 
92  call this%BaseModelType%allocate_scalars(modelname)
93  allocate (this%bndlist)
94  allocate (this%filename)
95  this%filename = ''
96  end subroutine allocate_scalars
97 
98  !> @brief Allocate array variables
99  subroutine allocate_arrays(this)
100  class(explicitmodeltype) :: this
101  integer(I4B) :: i
102 
103  call mem_allocate(this%nja, 'NJA', this%memoryPath)
104  this%nja = this%dis%nja
105 
106  call mem_allocate(this%ibound, this%dis%nodes, 'IBOUND', this%memoryPath)
107  do i = 1, this%dis%nodes
108  this%ibound(i) = 1 ! active by default
109  end do
110 
111  call mem_allocate(this%flowja, this%nja, 'FLOWJA', this%memoryPath)
112  do i = 1, this%nja
113  this%flowja(i) = dzero
114  end do
115  end subroutine allocate_arrays
116 
117  !> @ brief Cast a generic object into an explicit model
118  function castasexplicitmodelclass(obj) result(res)
119  class(*), pointer, intent(inout) :: obj
120  class(explicitmodeltype), pointer :: res
121 
122  res => null()
123  if (.not. associated(obj)) return
124 
125  select type (obj)
126  class is (explicitmodeltype)
127  res => obj
128  end select
129  end function castasexplicitmodelclass
130 
131  !> @ brief Add explicit model to a generic list
132  subroutine addexplicitmodeltolist(list, model)
133  ! dummy
134  type(listtype), intent(inout) :: list
135  class(explicitmodeltype), pointer, intent(inout) :: model
136  ! local
137  class(*), pointer :: obj
138 
139  obj => model
140  call list%Add(obj)
141  end subroutine addexplicitmodeltolist
142 
143  !> @ brief Get generic object from list and return as explicit model
144  function getexplicitmodelfromlist(list, idx) result(res)
145  ! dummy
146  type(listtype), intent(inout) :: list
147  integer(I4B), intent(in) :: idx
148  class(explicitmodeltype), pointer :: res
149  ! local
150  class(*), pointer :: obj
151 
152  obj => list%GetItem(idx)
153  res => castasexplicitmodelclass(obj)
154  end function getexplicitmodelfromlist
155 
156  !> @brief Set the solution ID
157  subroutine set_idsoln(this, id)
158  class(explicitmodeltype) :: this
159  integer(I4B), intent(in) :: id
160  this%idsoln = id
161  end subroutine set_idsoln
162 
163 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.
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