1+ ! > @brief This module contains methods for calculating the shortwave radiation heat flux
2+ ! !
3+ ! ! This module contains the methods used to calculate the shortwave radiation heat flux
4+ ! ! for surface-water boundaries, like streams and lakes. In its current form,
5+ ! ! this class acts like a package to a package, similar to the TVK package that
6+ ! ! can be invoked from the NPF package. Once this package is completed in its
7+ ! ! prototyped form, it will likely be moved around.
8+ ! <
9+
10+ module ShortwaveModule
11+ use ConstantsModule, only: LINELENGTH, LENMEMPATH, DZERO, LENVARNAME
12+ use KindModule, only: I4B, DP
13+ use MemoryManagerModule, only: mem_setptr
14+ use MemoryHelperModule, only: create_mem_path
15+ use SimModule, only: store_error
16+ use SimVariablesModule, only: errmsg
17+ use PbstBaseModule, only: PbstBaseType, pbstbase_da
18+
19+ implicit none
20+
21+ private
22+
23+ public :: SwrType
24+ public :: swr_cr
25+
26+ character (len= 16 ) :: text = ' SWR'
27+
28+ type, extends(PbstBaseType) :: SwrType
29+
30+ real (DP), dimension (:), pointer , contiguous :: solr = > null () ! < solar radiation
31+ real (DP), dimension (:), pointer , contiguous :: shd = > null () ! < shade fraction
32+ real (DP), dimension (:), pointer , contiguous :: swrefl = > null () ! < shortwave reflectance of water surface
33+
34+ contains
35+
36+ procedure :: da = > swr_da
37+ procedure :: read_option = > swr_read_option
38+ ! procedure :: pbst_options => swr_options
39+ procedure :: subpck_set_stressperiod = > swr_set_stressperiod
40+ procedure :: pbst_allocate_arrays = > swr_allocate_arrays
41+ ! procedure, private :: swr_allocate_scalars
42+ procedure , public :: swr_cq
43+
44+ end type SwrType
45+
46+ contains
47+
48+ ! > @brief Create a new SwrType object
49+ ! !
50+ ! ! Create a new shortwave radiation flux (SwrType) object. Initially for use with
51+ ! ! the SFE package.
52+ ! <
53+ subroutine swr_cr (swr , name_model , inunit , iout , ncv )
54+ ! -- dummy
55+ type (SwrType), pointer , intent (out ) :: swr
56+ character (len=* ), intent (in ) :: name_model
57+ integer (I4B), intent (in ) :: inunit
58+ integer (I4B), intent (in ) :: iout
59+ integer (I4B), target , intent (in ) :: ncv
60+ !
61+ allocate (swr)
62+ call swr% init(name_model, ' SWR' , ' SWR' , inunit, iout, ncv)
63+ swr% text = text
64+ !
65+ ! -- allocate scalars
66+ ! call swr%swr_allocate_scalars()
67+ end subroutine swr_cr
68+
69+ ! > @brief Allocate scalars specific to the streamflow energy transport (SFE)
70+ ! ! package.
71+ ! ! NO SCALARS IN THIS YET
72+ ! <
73+ ! subroutine swr_allocate_scalars(this)
74+ ! -- modules
75+ ! use MemoryManagerModule, only: mem_allocate
76+ ! -- dummy
77+ ! class(SwrType) :: this
78+ !
79+ ! -- allocate
80+ ! call mem_allocate(this%solr, 'SOLR', this%memoryPath)
81+ ! call mem_allocate(this%shd, 'SHD', this%memoryPath)
82+ ! call mem_allocate(this%swrefl, 'SWREFL', this%memoryPath)
83+ !
84+ ! -- initialize to default values
85+ ! this%solr = 600 ! W/m3
86+ ! this%shd = 0.25 ! dimensionless fraction
87+ ! this%swrefl = 0.03 ! dimensionless fraction
88+ ! end subroutine swr_allocate_scalars
89+
90+ ! > @brief Allocate arrays specific to the sensible heat flux (SWR) package
91+ ! <
92+ subroutine swr_allocate_arrays (this )
93+ ! -- modules
94+ use MemoryManagerModule, only: mem_allocate
95+ ! -- dummy
96+ class(SwrType), intent (inout ) :: this
97+ ! -- local
98+ integer (I4B) :: n
99+ !
100+ ! -- time series
101+ call mem_allocate(this% solr, this% ncv, ' SOLR' , this% memoryPath)
102+ call mem_allocate(this% shd, this% ncv, ' SHD' , this% memoryPath)
103+ call mem_allocate(this% swrefl, this% ncv, ' SWREFL' , this% memoryPath)
104+ !
105+ ! -- initialize
106+ do n = 1 , this% ncv
107+ this% solr(n) = DZERO
108+ this% shd(n) = DZERO
109+ this% swrefl(n) = DZERO
110+ end do
111+ end subroutine
112+
113+ ! > @brief Set options specific to the SwrType
114+ ! !
115+ ! NOT USED RIGHT NOW BECAUSE NO CONSTANTS
116+ ! ! This routine overrides PbstBaseType%bnd_options
117+ ! <
118+ ! subroutine swr_options(this, option, found)
119+ ! -- dummy
120+ ! class(SwrType), intent(inout) :: this
121+ ! character(len=*), intent(inout) :: option
122+ ! logical, intent(inout) :: found
123+ !
124+ ! found = .true.
125+ ! select case (option)
126+ ! case ('DENSITY_AIR')
127+ ! this%rhoa = this%parser%GetDouble()
128+ ! if (this%rhoa <= 0.0) then
129+ ! write (errmsg, '(a)') 'Specified value for the density of &
130+ ! &the atmosphere must be greater than 0.0.'
131+ ! call store_error(errmsg)
132+ ! call this%parser%StoreErrorUnit()
133+ ! else
134+ ! write (this%iout, '(4x,a,1pg15.6)') &
135+ ! "The density of the atmosphere has been set to: ", this%rhoa
136+ ! end if
137+ ! case ('HEAT_CAPACITY_AIR')
138+ ! this%cpa = this%parser%GetDouble()
139+ ! if (this%cpa <= 0.0) then
140+ ! write (errmsg, '(a)') 'Specified value for the heat capacity of &
141+ ! &the atmosphere must be greater than 0.0.'
142+ ! call store_error(errmsg)
143+ ! call this%parser%StoreErrorUnit()
144+ ! else
145+ ! write (this%iout, '(4x,a,1pg15.6)') &
146+ ! "The heat capacity of the atmosphere has been set to: ", this%cpa
147+ ! end if
148+ ! case ('DRAG_COEFFICIENT')
149+ ! this%cd = this%parser%GetDouble()
150+ ! if (this%cd <= 0.0) then
151+ ! write (errmsg, '(a)') 'Specified value for the drag coefficient &
152+ ! &must be greater than 0.0.'
153+ ! call store_error(errmsg)
154+ ! call this%parser%StoreErrorUnit()
155+ ! else
156+ ! write (this%iout, '(4x,a,1pg15.6)') &
157+ ! "The heat capacity of the atmosphere has been set to: ", this%cpa
158+ ! end if
159+ ! case default
160+ ! write (errmsg, '(a,a)') 'Unknown SHF option: ', trim(option)
161+ ! call store_error(errmsg)
162+ ! call this%parser%StoreErrorUnit()
163+ ! end select
164+ ! end subroutine swr_options
165+
166+ ! > @brief Calculate Shortwave Radiation Heat Flux
167+ ! !
168+ ! ! Calculate and return the shortwave radiation heat flux for one reach
169+ ! <
170+ subroutine swr_cq (this , ifno , swrflx )
171+ ! -- dummy
172+ class(SwrType), intent (inout ) :: this
173+ integer (I4B), intent (in ) :: ifno ! < stream reach integer id
174+ ! real(DP), intent(in) :: tstrm !< temperature of the stream reach
175+ real (DP), intent (inout ) :: swrflx ! < calculated shortwave radiation heat flux amount
176+ ! -- local
177+ ! real(DP) :: swr_const
178+ !
179+ ! -- calculate shortwave radiation heat flux (version: user input of sol rad data)
180+ swrflx = (1 - this% shd(ifno)) * (1 - this% swrefl(ifno)) * this% solr(ifno)
181+ end subroutine swr_cq
182+
183+ ! > @brief Deallocate package memory
184+ ! !
185+ ! ! Deallocate TVK package scalars and arrays.
186+ ! <
187+ subroutine swr_da (this )
188+ ! -- modules
189+ use MemoryManagerModule, only: mem_deallocate
190+ ! -- dummy
191+ class(SwrType) :: this
192+ !
193+ ! -- Nullify pointers to other package variables
194+ ! call mem_deallocate(this%rhoa)
195+ ! call mem_deallocate(this%cpa)
196+ ! call mem_deallocate(this%cd)
197+ !
198+ ! -- Deallocate time series
199+ call mem_deallocate(this% shd)
200+ call mem_deallocate(this% swrefl)
201+ call mem_deallocate(this% solr)
202+ !
203+ ! -- Deallocate parent
204+ call pbstbase_da(this)
205+ end subroutine swr_da
206+
207+ ! > @brief Read a SWR-specific option from the OPTIONS block
208+ ! !
209+ ! ! Process a single SWR-specific option. Used when reading the OPTIONS block
210+ ! ! of the SWR package input file.
211+ ! <
212+ function swr_read_option (this , keyword ) result(success)
213+ ! -- dummy
214+ class(SwrType) :: this
215+ character (len=* ), intent (in ) :: keyword
216+ ! -- return
217+ logical :: success
218+ !
219+ ! -- There are no SWR-specific options, so just return false
220+ success = .false.
221+ end function swr_read_option
222+
223+ ! > @brief Set the stress period attributes based on the keyword
224+ ! <
225+ subroutine swr_set_stressperiod (this , itemno , keyword , found )
226+ ! -- module
227+ use TimeSeriesManagerModule, only: read_value_or_time_series_adv
228+ ! -- dummy
229+ class(SwrType), intent (inout ) :: this
230+ integer (I4B), intent (in ) :: itemno
231+ character (len=* ), intent (in ) :: keyword
232+ logical , intent (inout ) :: found
233+ ! -- local
234+ character (len= LINELENGTH) :: text
235+ integer (I4B) :: ierr
236+ integer (I4B) :: jj
237+ real (DP), pointer :: bndElem = > null ()
238+ !
239+ ! <shd> SHADE
240+ ! <swrefl> REFLECTANCE OF SHORTWAVE RADIATION OFF WATER SURFACE
241+ ! <solr> SOLAR RADIATION
242+ !
243+ found = .true.
244+ select case (keyword)
245+ case (' SHD' )
246+ ierr = this% pbst_check_valid(itemno)
247+ if (ierr /= 0 ) then
248+ goto 999
249+ end if
250+ call this% parser% GetString(text)
251+ jj = 1
252+ bndElem = > this% shd(itemno)
253+ call read_value_or_time_series_adv(text, itemno, jj, bndElem, &
254+ this% packName, ' BND' , this% tsManager, &
255+ this% iprpak, ' SHD' )
256+ case (' SWREFL' )
257+ ierr = this% pbst_check_valid(itemno)
258+ if (ierr /= 0 ) then
259+ goto 999
260+ end if
261+ call this% parser% GetString(text)
262+ jj = 1
263+ bndElem = > this% swrefl(itemno)
264+ call read_value_or_time_series_adv(text, itemno, jj, bndElem, &
265+ this% packName, ' BND' , this% tsManager, &
266+ this% iprpak, ' SWREFL' )
267+ case (' SOLR' )
268+ ierr = this% pbst_check_valid(itemno)
269+ if (ierr /= 0 ) then
270+ goto 999
271+ end if
272+ call this% parser% GetString(text)
273+ jj = 1
274+ bndElem = > this% solr(itemno)
275+ call read_value_or_time_series_adv(text, itemno, jj, bndElem, &
276+ this% packName, ' BND' , this% tsManager, &
277+ this% iprpak, ' SOLR' )
278+ case default
279+ !
280+ ! -- Keyword not recognized so return to caller with found = .false.
281+ found = .false.
282+ end select
283+ !
284+ 999 continue
285+ end subroutine swr_set_stressperiod
286+
287+ end module ShortwaveModule
0 commit comments