MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
explicitmodelmodule Module Reference

Models that solve themselves.

Data Types

type  explicitmodeltype
 Base type for models that solve themselves. More...
 

Functions/Subroutines

subroutine model_ad (this)
 @ brief Advance the model More...
 
subroutine model_solve (this)
 @ brief Solve the model More...
 
subroutine model_cq (this, icnvg, isuppress_output)
 @ brief Calculate model flows More...
 
subroutine model_bd (this, icnvg, isuppress_output)
 @ brief Calculate model budget More...
 
subroutine model_da (this)
 @ brief Deallocate the model More...
 
subroutine allocate_scalars (this, modelname)
 @ brief Allocate scalar variables More...
 
subroutine allocate_arrays (this)
 Allocate array variables. More...
 
class(explicitmodeltype) function, pointer, public castasexplicitmodelclass (obj)
 @ brief Cast a generic object into an explicit model More...
 
subroutine, public addexplicitmodeltolist (list, model)
 @ brief Add explicit model to a generic list More...
 
class(explicitmodeltype) function, pointer, public getexplicitmodelfromlist (list, idx)
 @ brief Get generic object from list and return as explicit model More...
 
subroutine set_idsoln (this, id)
 Set the solution ID. More...
 
logical function has_pending (this)
 Whether the model has any pending work. More...
 

Function/Subroutine Documentation

◆ addexplicitmodeltolist()

subroutine, public explicitmodelmodule::addexplicitmodeltolist ( type(listtype), intent(inout)  list,
class(explicitmodeltype), intent(inout), pointer  model 
)

Definition at line 131 of file ExplicitModel.f90.

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)
Here is the caller graph for this function:

◆ allocate_arrays()

subroutine explicitmodelmodule::allocate_arrays ( class(explicitmodeltype this)
private

Definition at line 98 of file ExplicitModel.f90.

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

◆ allocate_scalars()

subroutine explicitmodelmodule::allocate_scalars ( class(explicitmodeltype this,
character(len=*), intent(in)  modelname 
)
private

Definition at line 87 of file ExplicitModel.f90.

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 = ''

◆ castasexplicitmodelclass()

class(explicitmodeltype) function, pointer, public explicitmodelmodule::castasexplicitmodelclass ( class(*), intent(inout), pointer  obj)

Definition at line 117 of file ExplicitModel.f90.

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
Here is the caller graph for this function:

◆ getexplicitmodelfromlist()

class(explicitmodeltype) function, pointer, public explicitmodelmodule::getexplicitmodelfromlist ( type(listtype), intent(inout)  list,
integer(i4b), intent(in)  idx 
)

Definition at line 143 of file ExplicitModel.f90.

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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ has_pending()

logical function explicitmodelmodule::has_pending ( class(explicitmodeltype), intent(in)  this)
private

Definition at line 163 of file ExplicitModel.f90.

164  class(ExplicitModelType), intent(in) :: this
165  has_pending = .false.

◆ model_ad()

subroutine explicitmodelmodule::model_ad ( class(explicitmodeltype this)
private

Definition at line 39 of file ExplicitModel.f90.

40  class(ExplicitModelType) :: this

◆ model_bd()

subroutine explicitmodelmodule::model_bd ( class(explicitmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)
private

Definition at line 56 of file ExplicitModel.f90.

57  class(ExplicitModelType) :: this
58  integer(I4B), intent(in) :: icnvg
59  integer(I4B), intent(in) :: isuppress_output

◆ model_cq()

subroutine explicitmodelmodule::model_cq ( class(explicitmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)
private

Definition at line 49 of file ExplicitModel.f90.

50  class(ExplicitModelType) :: this
51  integer(I4B), intent(in) :: icnvg
52  integer(I4B), intent(in) :: isuppress_output

◆ model_da()

subroutine explicitmodelmodule::model_da ( class(explicitmodeltype this)
private

Definition at line 63 of file ExplicitModel.f90.

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()

◆ model_solve()

subroutine explicitmodelmodule::model_solve ( class(explicitmodeltype this)
private

Definition at line 44 of file ExplicitModel.f90.

45  class(ExplicitModelType) :: this

◆ set_idsoln()

subroutine explicitmodelmodule::set_idsoln ( class(explicitmodeltype this,
integer(i4b), intent(in)  id 
)
private

Definition at line 156 of file ExplicitModel.f90.

157  class(ExplicitModelType) :: this
158  integer(I4B), intent(in) :: id
159  this%idsoln = id