ROMS
Loading...
Searching...
No Matches
ice_spdiw.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#ifdef ICE_MODEL
5!
6!git $Id$
7!=======================================================================
8! Copyright (c) 2002-2025 The ROMS Group !
9! Licensed under a MIT/X style license W. Paul Budgell !
10! See License_ROMS.md Katherine Hedstrom !
11!================================================== Hernan G. Arango ===
12! !
13! This module computes the magnitude of the shear between the ice !
14! and the surface water. In this case, the surface water is defined !
15! as the water in a surface mixed layer, so that velocity must be !
16! computed first. !
17! !
18!=======================================================================
19!
20 USE mod_param
21 USE mod_coupling
22 USE mod_forces
23 USE mod_grid
24 USE mod_ice
25# ifdef LMD_SKPP
26 USE mod_mixing
27# endif
28 USE mod_ocean
29 USE mod_scalars
30!
31 USE bc_2d_mod
32# ifdef DISTRIBUTE
34# endif
35!
36 implicit none
37!
38 PRIVATE
39 PUBLIC :: ice_spdiw
40!
41 CONTAINS
42!
43!***********************************************************************
44 SUBROUTINE ice_spdiw (ng, tile, model)
45!***********************************************************************
46!
47 USE mod_stepping
48!
49! Imported variable declarations.
50!
51 integer, intent(in) :: ng, tile, model
52!
53! Local variable declarations.
54!
55 character (len=*), parameter :: MyFile = &
56 & __FILE__
57!
58# include "tile.h"
59!
60# ifdef PROFILE
61 CALL wclock_on (ng, model, 42, __line__, myfile)
62# endif
63 CALL ice_spdiw_tile (ng, tile, model, &
64 & lbi, ubi, lbj, ubj, &
65 & imins, imaxs, jmins, jmaxs, &
66 & nrhs(ng), &
67# ifdef ICE_MODEL
68 & liuol(ng), &
69# endif
70 & grid(ng) % z_r, &
71 & grid(ng) % z_w, &
72 & ocean(ng) % u, &
73 & ocean(ng) % v, &
74# ifdef LMD_SKPP
75 & mixing(ng) % hsbl, &
76# endif
77 & ice(ng) % Fi, &
78 & ice(ng) % Si)
79# ifdef PROFILE
80 CALL wclock_off (ng, model, 42, __line__, myfile)
81# endif
82!
83 RETURN
84 END SUBROUTINE ice_spdiw
85!
86!***********************************************************************
87 SUBROUTINE ice_spdiw_tile (ng, tile, model, &
88 & LBi, UBi, LBj, UBj, &
89 & IminS, ImaxS, JminS, JmaxS, &
90 & nrhs, &
91# ifdef ICE_MODEL
92 & liuol, &
93# endif
94 & z_r, z_w, &
95 & u, v, &
96# ifdef LMD_SKPP
97 & hsbl, &
98# endif
99 & Fi, Si)
100!***********************************************************************
101!
102! Imported variable declarations.
103!
104 integer, intent(in) :: ng, tile, model
105 integer, intent(in) :: LBi, UBi, LBj, UBj
106 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
107 integer, intent(in) :: nrhs
108 integer, intent(in) :: liuol
109!
110# ifdef ASSUMED_SHAPE
111 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
112 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
113 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
114 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
115# ifdef LMD_SKPP
116 real(r8), intent(in) :: hsbl(LBi:,LBj:)
117# endif
118 real(r8), intent(in) :: Si(LBi:,LBj:,:,:)
119 real(r8), intent(inout) :: Fi(LBi:,LBj:,:)
120# else
121 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
122 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
123 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
124 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
125# ifdef LMD_SKPP
126 real(r8), intent(in) :: hsbl(LBi:UBi,LBj:UBj)
127# endif
128 real(r8), intent(in) :: Si(LBi:UBi,LBj:UBj,2,nIceS)
129 real(r8), intent(inout) :: Fi(LBi:UBi,LBj:UBj,nIceF)
130# endif
131!
132! Local variable declarations.
133!
134 integer :: i, j
135 integer :: nlio, nbotu, nbotv, k
136!
137 integer, dimension(IminS:ImaxS,JminS:JmaxS) :: nbot
138!
139 real(r8) :: dml, mlio, totml
140!
141 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: uw
142 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vw
143
144# include "set_bounds.h"
145!
146!-----------------------------------------------------------------------
147! Compute magnitude of the shear between the ice and the surface water.
148!-----------------------------------------------------------------------
149!
150! Determine model level associated with the depth of the mixed layer.
151!
152 DO j=jstrm2,jendp2
153 DO i=istrm2,iendp2
154# ifdef LMD_SKPP
155 mlio=min(-abs(hsbl(i,j)),-10.0_r8)
156# else
157 mlio=-10.0_r8
158# endif
159 nbot(i,j)=1
160 DO k=n(ng),1,-1
161 IF (z_r(i,j,k).lt.mlio) THEN
162 nbot(i,j)=min(k,n(ng))
163 nbot(i,j)=max(nbot(i,j),1)
164 EXIT
165 END IF
166 END DO
167 END DO
168 END DO
169!
170! Compute verticaly averaged U-velocity over mixed-layer thickness.
171!
172 DO j=jstr,jend
173 DO i=istru-1,iend+1
174 nlio=0
175 nbotu=nint(0.5_r8*(nbot(i-1,j)+nbot(i,j)))
176 nbotu=max(min(nbotu,n(ng)),1)
177 uw(i,j)=0.0_r8
178 totml=0.0_r8
179 DO k=n(ng),nbotu,-1
180 nlio=nlio+1
181 dml=0.5_r8*(z_w(i-1,j,k)-z_w(i-1,j,k-1)+ &
182 & z_w(i ,j,k)-z_w(i ,j,k-1))
183 uw(i,j)=uw(i,j)+u(i,j,k,nrhs)*dml
184 totml=totml+dml
185 END DO
186 uw(i,j)=uw(i,j)/totml
187 END DO
188 END DO
189!
190! Compute verticaly averaged U-velocity over mixed-layer thickness.
191!
192 DO j=jstrv-1,jend+1
193 DO i=istr,iend
194 nlio=0
195 nbotv=nint(0.5_r8*(nbot(i,j-1)+nbot(i,j)))
196 nbotv=max(min(nbotv,n(ng)),1)
197 vw(i,j)=0.0_r8
198 totml=0.0_r8
199 DO k=n(ng),nbotv,-1
200 nlio=nlio+1
201 dml=0.5_r8*(z_w(i,j-1,k)-z_w(i,j-1,k-1)+ &
202 & z_w(i,j ,k)-z_w(i,j ,k-1))
203 vw(i,j)=vw(i,j)+v(i,j,k,nrhs)*dml
204 totml=totml+dml
205 END DO
206 vw(i,j)=vw(i,j)/totml
207 END DO
208 END DO
209!
210! Compute magnitude of the shear between ice velocities and vertically
211! averaged mixed layer velocities.
212!
213 DO j=jstr,jend
214 DO i=istr,iend
215 fi(i,j,iciovs)=0.5_r8* &
216 & sqrt((uw(i ,j)-si(i ,j,liuol,isuice)+ &
217 & uw(i+1,j)-si(i+1,j,liuol,isuice))**2+ &
218 & (vw(i,j )-si(i,j ,liuol,isvice)+ &
219 & vw(i,j+1)-si(i,j+1,liuol,isvice))**2)
220 END DO
221 END DO
222!
223! Load vertically averaged mixed layer velocity components.
224!
225 DO j=jstr,jend
226 DO i=istrp,iend
227 fi(i,j,icuavg)=uw(i,j)
228 END DO
229 END DO
230 DO j=jstrp,jend
231 DO i=istr,iend
232 fi(i,j,icvavg)=vw(i,j)
233 END DO
234 END DO
235!
236! Set lateral boundary conditions.
237!
238 CALL bc_r2d_tile (ng, tile, &
239 & lbi, ubi, lbj, ubj, &
240 & fi(:,:,iciovs))
241
242 CALL bc_u2d_tile (ng, tile, &
243 & lbi, ubi, lbj, ubj, &
244 & fi(:,:,icuavg))
245
246 CALL bc_v2d_tile (ng, tile, &
247 & lbi, ubi, lbj, ubj, &
248 & fi(:,:,icvavg))
249
250# ifdef DISTRIBUTE
251 CALL mp_exchange2d (ng, tile, model, 3, &
252 & lbi, ubi, lbj, ubj, &
253 & nghostpoints, ewperiodic(ng), nsperiodic(ng), &
254 & fi(:,:,iciovs), &
255 & fi(:,:,icuavg), &
256 & fi(:,:,icvavg))
257# endif
258!
259 RETURN
260 END SUBROUTINE ice_spdiw_tile
261#endif
262 END MODULE ice_spdiw_mod
subroutine bc_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition bc_2d.F:44
subroutine bc_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition bc_2d.F:345
subroutine bc_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition bc_2d.F:167
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer, parameter isvice
Definition ice_mod.h:147
integer, parameter icuavg
Definition ice_mod.h:187
type(t_ice), dimension(:), allocatable ice
Definition ice_mod.h:283
integer, parameter icvavg
Definition ice_mod.h:188
integer, parameter iciovs
Definition ice_mod.h:178
integer, parameter isuice
Definition ice_mod.h:146
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
integer, dimension(:), allocatable nrhs
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3