MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
SolutionGroup.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
6  use listmodule, only: listtype
7 
8  implicit none
9  private
12  private :: castassolutiongroupclass
13 
15  integer(I4B), pointer :: id
16  integer(I4B), pointer :: mxiter
17  integer(I4B), pointer :: nsolutions
18  integer(I4B), dimension(:), allocatable :: idsolutions !array of solution ids in basesolutionlist
19 
20  contains
21 
22  procedure :: sgp_ca
23  procedure :: sgp_da
24  procedure, private :: allocate_scalars
25  procedure :: add_solution
26  end type solutiongrouptype
27 
28 contains
29 
30  !> @brief Create a new solution group
31  !<
32  subroutine solutiongroup_create(sgp, id)
33  ! -- dummy
34  type(solutiongrouptype), pointer :: sgp
35  integer(I4B), intent(in) :: id
36  !
37  allocate (sgp)
38  call sgp%allocate_scalars()
39  sgp%id = id
40  end subroutine solutiongroup_create
41 
42  !> @brief Calculate the solution group
43  !!
44  !! Solve each solution group and each solution. Start with converge
45  !! flag equal true and reset to zero if any non-convergence triggers
46  !! are encountered.
47  !<
48  subroutine sgp_ca(this)
49  ! -- modules
50  use constantsmodule, only: linelength
52  use tdismodule, only: kstp, kper
53  ! -- dummy
54  class(solutiongrouptype) :: this
55  ! -- local
56  class(basesolutiontype), pointer :: sp
57  integer(I4B) :: kpicard, isgcnvg, isuppress_output
58  integer(I4B) :: is, isoln
59  ! -- formats
60  character(len=*), parameter :: fmtnocnvg = &
61  "(1X,'Solution Group ', i0, ' did not converge for stress period ', i0, &
62  &' and time step ', i0)"
63  !
64  ! -- Suppress output during picard iterations
65  if (this%mxiter > 1) then
66  isuppress_output = 1
67  else
68  isuppress_output = 0
69  end if
70  !
71  ! -- set failed flag
72  laststepfailed = 0
73  !
74  ! -- Picard loop
75  picardloop: do kpicard = 1, this%mxiter
76  if (this%mxiter > 1) then
77  write (iout, '(/a,i6/)') 'SOLUTION GROUP PICARD ITERATION: ', kpicard
78  end if
79  isgcnvg = 1
80  do is = 1, this%nsolutions
81  isoln = this%idsolutions(is)
83  call sp%sln_ca(isgcnvg, isuppress_output)
84  end do
85  if (isgcnvg == 1) exit picardloop
86  end do picardloop
87  !
88  ! -- if a picard loop was used and the solution group converged
89  ! then rerun the timestep and save the output. Or if there
90  ! is only one picard iteration, then do nothing as models
91  ! are assumed to be explicitly coupled.
92  if (isgcnvg == 1) then
93  if (this%mxiter > 1) then
94  isuppress_output = 0
95  do is = 1, this%nsolutions
96  isoln = this%idsolutions(is)
98  call sp%sln_ca(isgcnvg, isuppress_output)
99  end do
100  end if
101  else
102  isimcnvg = 0
103  laststepfailed = 1
104  write (iout, fmtnocnvg) this%id, kper, kstp
105  end if
106  end subroutine sgp_ca
107 
108  !> @brief Deallocate
109  !<
110  subroutine sgp_da(this)
111  ! -- dummy
112  class(solutiongrouptype) :: this
113  !
114  deallocate (this%id)
115  deallocate (this%mxiter)
116  deallocate (this%nsolutions)
117  deallocate (this%idsolutions)
118  end subroutine sgp_da
119 
120  !> @brief Allocate scalars
121  !<
122  subroutine allocate_scalars(this)
123  ! -- dummy
124  class(solutiongrouptype) :: this
125  !
126  allocate (this%id)
127  allocate (this%mxiter)
128  allocate (this%nsolutions)
129  this%id = 0
130  this%mxiter = 1
131  this%nsolutions = 0
132  end subroutine allocate_scalars
133 
134  !> @brief Add solution
135  !<
136  subroutine add_solution(this, isoln, sp)
137  ! -- modules
139  ! -- dummy
140  class(solutiongrouptype) :: this
141  integer(I4B), intent(in) :: isoln
142  class(basesolutiontype), pointer, intent(in) :: sp
143  ! -- local
144  integer(I4B) :: ipos
145  !
146  call expandarray(this%idsolutions)
147  ipos = size(this%idsolutions)
148  this%idsolutions(ipos) = isoln
149  this%nsolutions = this%nsolutions + 1
150  end subroutine add_solution
151 
152  function castassolutiongroupclass(obj) result(res)
153  implicit none
154  ! -- dummy
155  class(*), pointer, intent(inout) :: obj
156  ! -- return
157  class(solutiongrouptype), pointer :: res
158  !
159  res => null()
160  if (.not. associated(obj)) return
161  !
162  select type (obj)
163  class is (solutiongrouptype)
164  res => obj
165  end select
166  end function castassolutiongroupclass
167 
168  subroutine addsolutiongrouptolist(list, solutiongroup)
169  implicit none
170  ! -- dummy
171  type(listtype), intent(inout) :: list
172  type(solutiongrouptype), pointer, intent(inout) :: solutiongroup
173  ! -- local
174  class(*), pointer :: obj
175  !
176  obj => solutiongroup
177  call list%Add(obj)
178  end subroutine addsolutiongrouptolist
179 
180  function getsolutiongroupfromlist(list, idx) result(res)
181  implicit none
182  ! -- dummy
183  type(listtype), intent(inout) :: list
184  integer(I4B), intent(in) :: idx
185  class(solutiongrouptype), pointer :: res
186  ! -- local
187  class(*), pointer :: obj
188  !
189  obj => list%GetItem(idx)
190  res => castassolutiongroupclass(obj)
191  end function getsolutiongroupfromlist
192 
193 end module solutiongroupmodule
subroutine, public addbasesolutiontolist(list, solution)
class(basesolutiontype) function, pointer, public getbasesolutionfromlist(list, idx)
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
This module defines variable data types.
Definition: kind.f90:8
type(listtype), public basesolutionlist
Definition: mf6lists.f90:19
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) laststepfailed
flag indicating if the last step failed (1) if last step failed; (0) otherwise (set in converge_check...
integer(i4b) iout
file unit number for simulation output
integer(i4b) isimcnvg
simulation convergence flag (1) if all objects have converged, (0) otherwise
class(solutiongrouptype) function, pointer, private castassolutiongroupclass(obj)
subroutine add_solution(this, isoln, sp)
Add solution.
class(solutiongrouptype) function, pointer, public getsolutiongroupfromlist(list, idx)
subroutine, public solutiongroup_create(sgp, id)
Create a new solution group.
subroutine sgp_da(this)
Deallocate.
subroutine sgp_ca(this)
Calculate the solution group.
subroutine, public addsolutiongrouptolist(list, solutiongroup)
subroutine allocate_scalars(this)
Allocate scalars.
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
A generic heterogeneous doubly-linked list.
Definition: List.f90:14