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