ROMS
Loading...
Searching...
No Matches
ana_winds.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_winds (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets surface wind components using an analytical !
12! expression. !
13! !
14!=======================================================================
15!
16 USE mod_param
17 USE mod_forces
18 USE mod_grid
19 USE mod_ncparam
20!
21! Imported variable declarations.
22!
23 integer, intent(in) :: ng, tile, model
24!
25! Local variable declarations.
26!
27 character (len=*), parameter :: MyFile = &
28 & __FILE__
29!
30#include "tile.h"
31!
32 CALL ana_winds_tile (ng, tile, model, &
33 & lbi, ubi, lbj, ubj, &
34 & imins, imaxs, jmins, jmaxs, &
35#ifdef SPHERICAL
36 & grid(ng) % lonr, &
37 & grid(ng) % latr, &
38#else
39 & grid(ng) % xr, &
40 & grid(ng) % yr, &
41#endif
42 & forces(ng) % Uwind, &
43 & forces(ng) % Vwind)
44!
45! Set analytical header file name used.
46!
47#ifdef DISTRIBUTE
48 IF (lanafile) THEN
49#else
50 IF (lanafile.and.(tile.eq.0)) THEN
51#endif
52 ananame(36)=myfile
53 END IF
54!
55 RETURN
56 END SUBROUTINE ana_winds
57!
58!***********************************************************************
59 SUBROUTINE ana_winds_tile (ng, tile, model, &
60 & LBi, UBi, LBj, UBj, &
61 & IminS, ImaxS, JminS, JmaxS, &
62#ifdef SPHERICAL
63 & lonr, latr, &
64#else
65 & xr, yr, &
66#endif
67 & Uwind, Vwind)
68!***********************************************************************
69!
70 USE mod_param
71 USE mod_scalars
72!
74#ifdef DISTRIBUTE
76#endif
77!
78! Imported variable declarations.
79!
80 integer, intent(in) :: ng, tile, model
81 integer, intent(in) :: LBi, UBi, LBj, UBj
82 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
83!
84#ifdef ASSUMED_SHAPE
85# ifdef SPHERICAL
86 real(r8), intent(in) :: lonr(LBi:,LBj:)
87 real(r8), intent(in) :: latr(LBi:,LBj:)
88# else
89 real(r8), intent(in) :: xr(LBi:,LBj:)
90 real(r8), intent(in) :: yr(LBi:,LBj:)
91# endif
92 real(r8), intent(out) :: Uwind(LBi:,LBj:)
93 real(r8), intent(out) :: Vwind(LBi:,LBj:)
94#else
95# ifdef SPHERICAL
96 real(r8), intent(in) :: lonr(LBi:UBi,LBj:UBj)
97 real(r8), intent(in) :: latr(LBi:UBi,LBj:UBj)
98# else
99 real(r8), intent(in) :: xr(LBi:UBi,LBj:UBj)
100 real(r8), intent(in) :: yr(LBi:UBi,LBj:UBj)
101# endif
102 real(r8), intent(out) :: Uwind(LBi:UBi,LBj:UBj)
103 real(r8), intent(out) :: Vwind(LBi:UBi,LBj:UBj)
104#endif
105!
106! Local variable declarations.
107!
108 integer :: i, j
109!
110 real(r8) :: Wdir, Wmag, cff, u_wind, v_wind
111
112#include "set_bounds.h"
113!
114!-----------------------------------------------------------------------
115! Set surface wind components (m/s) at RHO-points.
116!-----------------------------------------------------------------------
117!
118#if defined BENCHMARK
119 wmag=15.0_r8
120 DO j=jstrt,jendt
121 DO i=istrt,iendt
122 cff=0.2_r8*(60.0_r8+latr(i,j))
123 uwind(i,j)=wmag*exp(-cff*cff)
124 vwind(i,j)=0.0_r8
125 END DO
126 END DO
127#elif defined BL_TEST
128 IF ((tdays(ng)-dstart).le.6.0_r8) THEN
129 u_wind=0.0_r8
130!! v_wind=4.7936_r8
131 v_wind=10.0_r8
132 END IF
133 DO j=jstrt,jendt
134 DO i=istrt,iendt
135 uwind(i,j)=u_wind
136 vwind(i,j)=v_wind
137 END DO
138 END DO
139#else
140 ana_winds.h: no values provided for uwind and vwind.
141#endif
142!
143! Exchange boundary data.
144!
145 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
146 CALL exchange_r2d_tile (ng, tile, &
147 & lbi, ubi, lbj, ubj, &
148 & uwind)
149 CALL exchange_r2d_tile (ng, tile, &
150 & lbi, ubi, lbj, ubj, &
151 & vwind)
152 END IF
153
154#ifdef DISTRIBUTE
155 CALL mp_exchange2d (ng, tile, model, 2, &
156 & lbi, ubi, lbj, ubj, &
157 & nghostpoints, &
158 & ewperiodic(ng), nsperiodic(ng), &
159 & uwind, vwind)
160#endif
161!
162 RETURN
163 END SUBROUTINE ana_winds_tile
subroutine ana_winds(ng, tile, model)
Definition ana_winds.h:3
subroutine ana_winds_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, lonr, latr, xr, yr, uwind, vwind)
Definition ana_winds.h:68
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
logical lanafile
character(len=256), dimension(39) ananame
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:), allocatable tdays
real(dp) dstart
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)