ROMS
Loading...
Searching...
No Matches
mod_bbl Module Reference

Data Types

type  t_bbl
 

Functions/Subroutines

subroutine, public allocate_bbl (ng, lbi, ubi, lbj, ubj)
 
subroutine, public deallocate_bbl (ng)
 
subroutine, public initialize_bbl (ng, tile)
 

Variables

type(t_bbl), dimension(:), allocatable bbl
 

Function/Subroutine Documentation

◆ allocate_bbl()

subroutine, public mod_bbl::allocate_bbl ( integer, intent(in) ng,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )

Definition at line 66 of file mod_bbl.F.

67!
68!=======================================================================
69! !
70! This routine allocates all variables in the module for all nested !
71! grids. !
72! !
73!=======================================================================
74!
75 USE mod_param
76!
77! Imported variable declarations.
78!
79 integer, intent(in) :: ng, LBi, UBi, LBj, UBj
80!
81! Local variable declarations.
82!
83 real(r8) :: size2d
84!
85!-----------------------------------------------------------------------
86! Allocate module variables.
87!-----------------------------------------------------------------------
88!
89 IF (ng.eq.1) allocate ( bbl(ngrids) )
90!
91! Set horizontal array size.
92!
93 size2d=real((ubi-lbi+1)*(ubj-lbj+1),r8)
94!
95! Bottom boundary layer arrays.
96!
97 allocate ( bbl(ng) % Iconv(lbi:ubi,lbj:ubj) )
98 dmem(ng)=dmem(ng)+size2d
99
100 allocate ( bbl(ng) % Ubot(lbi:ubi,lbj:ubj) )
101 dmem(ng)=dmem(ng)+size2d
102
103 allocate ( bbl(ng) % Ur(lbi:ubi,lbj:ubj) )
104 dmem(ng)=dmem(ng)+size2d
105
106 allocate ( bbl(ng) % Vbot(lbi:ubi,lbj:ubj) )
107 dmem(ng)=dmem(ng)+size2d
108
109 allocate ( bbl(ng) % Vr(lbi:ubi,lbj:ubj) )
110 dmem(ng)=dmem(ng)+size2d
111
112 allocate ( bbl(ng) % bustrc(lbi:ubi,lbj:ubj) )
113 dmem(ng)=dmem(ng)+size2d
114
115 allocate ( bbl(ng) % bvstrc(lbi:ubi,lbj:ubj) )
116 dmem(ng)=dmem(ng)+size2d
117
118 allocate ( bbl(ng) % bustrw(lbi:ubi,lbj:ubj) )
119 dmem(ng)=dmem(ng)+size2d
120
121 allocate ( bbl(ng) % bvstrw(lbi:ubi,lbj:ubj) )
122 dmem(ng)=dmem(ng)+size2d
123
124 allocate ( bbl(ng) % bustrcwmax(lbi:ubi,lbj:ubj) )
125 dmem(ng)=dmem(ng)+size2d
126
127 allocate ( bbl(ng) % bvstrcwmax(lbi:ubi,lbj:ubj) )
128 dmem(ng)=dmem(ng)+size2d
129
130 RETURN
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
integer ngrids
Definition mod_param.F:113

References bbl, mod_param::dmem, mod_param::ngrids, and mod_kinds::r8.

Referenced by mod_arrays::roms_allocate_arrays().

Here is the caller graph for this function:

◆ deallocate_bbl()

subroutine, public mod_bbl::deallocate_bbl ( integer, intent(in) ng)

Definition at line 133 of file mod_bbl.F.

134!
135!=======================================================================
136! !
137! This routine deallocates all variables in the module for all nested !
138! grids. !
139! !
140!=======================================================================
141!
142 USE mod_param, ONLY : ngrids
143# ifdef SUBOBJECT_DEALLOCATION
144 USE destroy_mod, ONLY : destroy
145# endif
146!
147! Imported variable declarations.
148!
149 integer, intent(in) :: ng
150!
151! Local variable declarations.
152!
153 character (len=*), parameter :: MyFile = &
154 & __FILE__//", deallocate_bbl"
155
156# ifdef SUBOBJECT_DEALLOCATION
157!
158!-----------------------------------------------------------------------
159! Deallocate each variable in the derived-type T_BBL structure
160! separately.
161!-----------------------------------------------------------------------
162!
163 IF (.not.destroy(ng, bbl(ng)%Iconv, myfile, &
164 & __line__, 'BBL(ng)%Iconv')) RETURN
165
166 IF (.not.destroy(ng, bbl(ng)%Ubot, myfile, &
167 & __line__, 'BBL(ng)%Ubot')) RETURN
168
169 IF (.not.destroy(ng, bbl(ng)%Ur, myfile, &
170 & __line__, 'BBL(ng)%Ur')) RETURN
171
172 IF (.not.destroy(ng, bbl(ng)%Vbot, myfile, &
173 & __line__, 'BBL(ng)%Vbot')) RETURN
174
175 IF (.not.destroy(ng, bbl(ng)%Vr, myfile, &
176 & __line__, 'BBL(ng)%Vr')) RETURN
177
178 IF (.not.destroy(ng, bbl(ng)%bustrc, myfile, &
179 & __line__, 'BBL(ng)%bustrc')) RETURN
180
181 IF (.not.destroy(ng, bbl(ng)%bvstrc, myfile, &
182 & __line__, 'BBL(ng)%bvstrc')) RETURN
183
184 IF (.not.destroy(ng, bbl(ng)%bustrw, myfile, &
185 & __line__, 'BBL(ng)%bustrw')) RETURN
186
187 IF (.not.destroy(ng, bbl(ng)%bvstrw, myfile, &
188 & __line__, 'BBL(ng)%bvstrw')) RETURN
189
190 IF (.not.destroy(ng, bbl(ng)%bustrcwmax, myfile, &
191 & __line__, 'BBL(ng)%bustrcwmax')) RETURN
192
193 IF (.not.destroy(ng, bbl(ng)%bvstrcwmax, myfile, &
194 & __line__, 'BBL(ng)%bvstrcwmax')) RETURN
195# endif
196!
197!-----------------------------------------------------------------------
198! Deallocate derived-type BBL structure.
199!-----------------------------------------------------------------------
200!
201 IF (ng.eq.ngrids) THEN
202 IF (allocated(bbl)) deallocate ( bbl )
203 END IF
204!
205 RETURN

References bbl, and mod_param::ngrids.

Referenced by mod_arrays::roms_deallocate_arrays().

Here is the caller graph for this function:

◆ initialize_bbl()

subroutine, public mod_bbl::initialize_bbl ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 208 of file mod_bbl.F.

209!
210!=======================================================================
211! !
212! This routine initialize all variables in the module using first !
213! touch distribution policy. In shared-memory configuration, this !
214! operation actually performs propagation of the "shared arrays" !
215! across the cluster, unless another policy is specified to !
216! override the default. !
217! !
218!=======================================================================
219!
220 USE mod_param
221!
222! Imported variable declarations.
223!
224 integer, intent(in) :: ng, tile
225!
226! Local variable declarations.
227!
228 integer :: Imin, Imax, Jmin, Jmax
229 integer :: i, j
230
231 real(r8), parameter :: IniVal = 0.0_r8
232
233# include "set_bounds.h"
234!
235! Set array initialization range.
236
237# ifdef DISTRIBUTE
238 imin=bounds(ng)%LBi(tile)
239 imax=bounds(ng)%UBi(tile)
240 jmin=bounds(ng)%LBj(tile)
241 jmax=bounds(ng)%UBj(tile)
242# else
243 IF (domain(ng)%Western_Edge(tile)) THEN
244 imin=bounds(ng)%LBi(tile)
245 ELSE
246 imin=istr
247 END IF
248 IF (domain(ng)%Eastern_Edge(tile)) THEN
249 imax=bounds(ng)%UBi(tile)
250 ELSE
251 imax=iend
252 END IF
253 IF (domain(ng)%Southern_Edge(tile)) THEN
254 jmin=bounds(ng)%LBj(tile)
255 ELSE
256 jmin=jstr
257 END IF
258 IF (domain(ng)%Northern_Edge(tile)) THEN
259 jmax=bounds(ng)%UBj(tile)
260 ELSE
261 jmax=jend
262 END IF
263# endif
264!
265!-----------------------------------------------------------------------
266! Initialize module variables.
267!-----------------------------------------------------------------------
268!
269 DO j=jmin,jmax
270 DO i=imin,imax
271 bbl(ng) % Iconv(i,j) = 0
272
273 bbl(ng) % Ubot(i,j) = inival
274 bbl(ng) % Ur(i,j) = inival
275
276 bbl(ng) % Vbot(i,j) = inival
277 bbl(ng) % Vr(i,j) = inival
278
279 bbl(ng) % bustrc(i,j) = inival
280 bbl(ng) % bvstrc(i,j) = inival
281
282 bbl(ng) % bustrw(i,j) = inival
283 bbl(ng) % bvstrw(i,j) = inival
284
285 bbl(ng) % bustrcwmax(i,j) = inival
286 bbl(ng) % bvstrcwmax(i,j) = inival
287 END DO
288 END DO
289!
290 RETURN
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329

References bbl, mod_param::bounds, and mod_param::domain.

Referenced by mod_arrays::roms_initialize_arrays().

Here is the caller graph for this function:

Variable Documentation

◆ bbl