/[MITgcm]/MITgcm/pkg/bulk_force/bulkf_forcing.F
ViewVC logotype

Contents of /MITgcm/pkg/bulk_force/bulkf_forcing.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.4 - (show annotations) (download)
Thu Oct 9 04:19:19 2003 UTC (20 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint51o_pre, checkpoint51l_post, checkpoint52, checkpoint51t_post, checkpoint51n_post, checkpoint51s_post, checkpoint51n_pre, checkpoint52b_pre, checkpoint51l_pre, checkpoint51q_post, checkpoint51r_post, checkpoint51i_post, checkpoint52a_pre, checkpoint51o_post, checkpoint52a_post, ecco_c52_e35, checkpoint51m_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh, tg2-branch, checkpoint51n_branch
Changes since 1.3: +3 -1 lines
 o first check-in for the "branch-genmake2" merge
 o verification suite as run on shelley (gcc 3.2.2):

Wed Oct  8 23:42:29 EDT 2003
                T           S           U           V
G D M    c        m  s        m  s        m  s        m  s
E p a R  g  m  m  e  .  m  m  e  .  m  m  e  .  m  m  e  .
N n k u  2  i  a  a  d  i  a  a  d  i  a  a  d  i  a  a  d
2 d e n  d  n  x  n  .  n  x  n  .  n  x  n  .  n  x  n  .

OPTFILE=NONE

Y Y Y Y 13 16 16 16  0 16 16 16 16 16 16 16 16 13 12  0  0 pass  adjustment.128x64x1
Y Y Y Y 16 16 16 16  0 16 16 16 16 16 16  0  0 16 16  0  0 pass  adjustment.cs-32x32x1
Y Y Y Y 16 16 16 16  0 16 16 16 16 16 16 22  0 16 16 22  0 pass  adjust_nlfs.cs-32x32x1
Y Y Y Y -- 13 13 16 16 13 13 13 13 16 16 16 16 16 16 16 16 N/O   advect_cs
Y Y Y Y -- 22 16 16 16 16 16 16 13 16 16 16 16 16 16 16 16 N/O   advect_xy
Y Y Y Y -- 13 16 13 16 16 16 16 16 16 16 22 16 16 16 16 16 N/O   advect_xz
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  aim.5l_cs
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 16 16 16 16 16 13 16 pass  aim.5l_Equatorial_Channel
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 13 16 16 13 13 16 pass  aim.5l_LatLon
Y Y Y Y 13 16 16 16 16 16 16 16 16 16 13 12 13 13 16 13 16 pass  exp0
Y Y Y Y 14 16 16 16 16 16 16 16 22 16 16 16 13 16 16 22 16 pass  exp1
Y Y Y Y 13 13 16 13 16 16 16 16 16 13 13 16 16 13 13 13 13 pass  exp2
Y Y Y Y 16 16 16 16 16 16 16 16 22 16 16 16 16 16 16 16 16 pass  exp4
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 22 16 16 16 22 16 pass  exp5
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  front_relax
Y Y Y Y 14 16 16 13 13 16 16 13 13 16 13 13 16 12 13 13 16 pass  global_ocean.90x40x15
Y Y Y Y 10 16 16 13 13 16 13 16 16 13 13 13 13 16 16 13 16 FAIL  global_ocean.cs32x15
Y Y Y Y  6 11 12 13 13 12 13 16 13  9  9  9  9 10  9  9 11 FAIL  global_ocean_pressure
Y Y Y Y 14 16 16 13 16 16 16 13 13 13 13 13 16 12 16 13 16 pass  global_with_exf
Y Y Y Y 14 16 16 16 16 16 16 16 16 11 13 22 13 16 16  9 16 pass  hs94.128x64x5
Y Y Y Y 13 16 16 16 16 16 16 16 16 11 16 16 16 13 16 22 13 pass  hs94.1x64x5
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 16 13 13 16 16 22 13 pass  hs94.cs-32x32x5
Y Y Y Y 10 10 16 13 13 16 16 16 22 16 13 13 13 13 13 22 13 FAIL  ideal_2D_oce
Y Y Y Y  8 16 16 16 16 16 16 16 16 13 13  8 16 16 16 16 16 FAIL  internal_wave
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 13 22 13 13 13 22 16 pass  inverted_barometer
Y Y Y Y 12 16 16 16 16 16 16 16 16 16 13 12 13 13 13 13 13 FAIL  lab_sea
Y Y Y Y 11 16 16 16 16 16 16 16 13 13 13 12 13 16 13 12 13 FAIL  natl_box
Y Y Y Y 16 16 16 16 16 16 16 16 22 16 16 16 16 16 16 16 16 pass  plume_on_slope
Y Y Y Y 13 16 16 16 16 13 16 16 16 16 16 16 16 13 16 16 16 pass  solid-body.cs-32x32x1

1 C $Header: /u/u3/gcmpack/MITgcm/pkg/bulk_force/bulkf_forcing.F,v 1.3.4.1 2003/10/02 18:18:33 adcroft Exp $
2 C $Name: $
3 c swd -- bulkf formula pkg
4
5 #include "BULK_FORCE_OPTIONS.h"
6
7 subroutine bulkf_forcing( bi,bj,
8 I mycurrenttime,
9 I mycurrentiter,
10 I mythid
11 I )
12
13 c ==================================================================
14 c SUBROUTINE BULKF_FORCING
15 c ==================================================================
16 c
17 c o Get the surface fluxes used to force ocean model
18 c Output:
19 c ------
20 c ustress, vstress - wind stress
21 c qnet - net heat flux
22 c empmr - freshwater flux
23 c ---------
24 c
25 c Input:
26 c ------
27 c uwind, vwind - mean wind speed (m/s) at height hu (m)
28 c Tair - mean air temperature (K) at height ht (m)
29 c Qair - mean air humidity (kg/kg) at height hq (m)
30 c theta(k=1) - sea surface temperature (K)
31 c rain - precipitation
32 c runoff - river(ice) runoff
33 c
34 c ==================================================================
35 c SUBROUTINE bulkf_forcing
36 c ==================================================================
37
38 implicit none
39
40 c == global variables ==
41
42 #include "EEPARAMS.h"
43 #include "SIZE.h"
44 #include "PARAMS.h"
45 #include "DYNVARS.h"
46 #include "GRID.h"
47 #include "FFIELDS.h"
48 #ifdef ALLOW_BULK_FORCE
49 #include "BULKF.h"
50 #include "BULKF_INT.h"
51 #include "BULKF_DIAG.h"
52 #endif
53 #ifdef ALLOW_THERM_SEAICE
54 #include "ICE.h"
55 #endif
56 c == routine arguments ==
57
58 integer mythid
59 integer mycurrentiter
60 _RL mycurrenttime
61 integer bi,bj
62
63 C == Local variables ==
64 integer i,j,k
65
66 #ifdef ALLOW_BULK_FORCE
67
68
69 _RL df0dT,hfl
70
71 c variables to include seaice effect
72 _RL tmp
73 _RL albedo
74 integer iceornot
75
76 c == external functions ==
77
78 #endif /* ALLOW_BULK_FORCE */
79
80 c determine forcing field records
81
82 #ifdef ALLOW_BULK_FORCE
83
84 k = 1
85 do j = 1-oly,sny+oly
86 do i = 1-olx,snx+olx
87 if (hFacC(i,j,1,bi,bj).ne.0.0) then
88 #ifdef ALLOW_THERM_SEAICE
89 if (ICEMASK(i,j,bi,bj).gt.0.d0) then
90 tmp=Tsrf(i,j,bi,bj)
91 if (snowheight(i,j,bi,bj).gt.3.d-1) then
92 iceornot=2
93 else
94 iceornot=1
95 endif
96 else
97 tmp=theta(i,j,1,bi,bj)
98 iceornot=0
99 endif
100 #else
101 tmp=theta(i,j,1,bi,bj)
102 iceornot=0
103 #endif
104
105
106 call BULKF_FORMULA_LANL(
107 & uwind(i,j,bi,bj), vwind(i,j,bi,bj),
108 & wspeed(i,j,bi,bj), Tair(i,j,bi,bj), Qair(i,j,bi,bj),
109 & cloud(i,j,bi,bj), tmp,
110 & flwup(i,j,bi,bj), flh(i,j,bi,bj),
111 & fsh(i,j,bi,bj), df0dT,
112 & ustress(i,j,bi,bj), vstress(i,j,bi,bj),
113 & evap(i,j,bi,bj), savssq(i,j,bi,bj),
114 & iceornot, readwindstress)
115
116 C Note that the LANL flux conventions are opposite
117 C of what they are in the model.
118
119 cQQ use down long wave data
120 flwupnet(i,j,bi,bj)=flwup(i,j,bi,bj)-flw(i,j,bi,bj)
121 cQQ using down solar, need to have water albedo -- .1
122 #ifdef ALLOW_THERM_SEAICE
123 if (ICEMASK(i,j,bi,bj).gt.0.d0) then
124 call sfc_albedo(ICEHEIGHT(i,j,bi,bj),
125 & SNOWHEIGHT(i,j,bi,bj),
126 & Tsrf(i,j,bi,bj),
127 & sage(i,j,bi,bj), albedo)
128 else
129 albedo=.1
130 endif
131 #else
132 albedo=.1
133 #endif
134 fswnet(i,j,bi,bj)=solar(i,j,bi,bj)*(1.d0-albedo)
135 else
136 ustress(i,j,bi,bj) = 0. _d 0
137 vstress(i,j,bi,bj) = 0. _d 0
138 fsh(i,j,bi,bj) = 0. _d 0
139 flh(i,j,bi,bj) = 0. _d 0
140 flwup(i,j,bi,bj) =0. _d 0
141 evap(i,j,bi,bj) =0. _d 0
142 fswnet(i,j,bi,bj) =0. _d 0
143 savssq(i,j,bi,bj) =0. _d 0
144 endif
145 enddo
146 enddo
147
148
149 if ( .NOT.readwindstress) then
150 cswd move wind stresses to u and v points
151 do j = 1,sny
152 do i = 1,snx
153 fu(i,j,bi,bj)=
154 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))/2*
155 & maskW(i,j,1,bi,bj)
156 fv(i,j,bi,bj)=
157 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))/2*
158 & maskS(i,j,1,bi,bj)
159 enddo
160 enddo
161 endif
162 c
163 c
164 c Add all contributions.
165 k = 1
166 do j = 1,sny
167 do i = 1,snx
168 if (hFacC(i,j,1,bi,bj).ne.0.0) then
169 c Net surface heat flux.
170 hfl = 0. _d 0
171 hfl = hfl + fsh(i,j,bi,bj)
172 hfl = hfl + flh(i,j,bi,bj)
173 hfl = hfl + flwupnet(i,j,bi,bj)
174 c QQ minus solar for ncep data
175 c QQ plus solar for dasilva
176 c hfl = hfl - solar(i,j,bi,bj)
177 cQQ for calculated sw net
178 hfl = hfl - fswnet(i,j,bi,bj)
179 c Heat flux:
180 Qnet(i,j,bi,bj) = -hfl*maskc(i,j,k,bi,bj)
181 cswdcou -- add ---
182 #ifdef COUPLE_MODEL
183 dFdT(i,j,bi,bj) = df0dT
184 #endif
185 cswdcou -- end add ---
186 c Salt flux from Precipitation and Evaporation.
187 EmPmR(i,j,bi,bj) = (-evap(i,j,bi,bj)+rain(i,j,bi,bj)
188 & - runoff(i,j,bi,bj))*maskc(i,j,k,bi,bj)
189
190 cccccccccccCHEATccccccccccccccccccccccccc
191 c Qnet(i,j,bi,bj) = Qnetch(i,j,bi,bj)
192 c EmPmR(i,j,bi,bj) = EmPch(i,j,bi,bj)
193 cccccccccccccccccccccccccccccccccccccccccc
194 else
195 Qnet(i,j,bi,bj) =0.d0
196 EmPmR(i,j,bi,bj)=0.d0
197 cswdcou -- add ---
198 #ifdef COUPLE_MODEL
199 dFdT(i,j,bi,bj) = 0.d0
200 #endif
201 cswdcou -- end add ---
202 endif
203 enddo
204 enddo
205
206
207 cc Update the tile edges.
208 c _EXCH_XY_R8(Qnet, mythid)
209 c _EXCH_XY_R8(EmPmR, mythid)
210 c _EXCH_XY_R8(fu , mythid)
211 c _EXCH_XY_R8(fv , mythid)
212
213
214 #ifdef ALLOW_TIMEAVE
215 call BULKF_AVE(bi,bj,mythid)
216 #endif /*ALLOW_TIMEAVE*/
217
218 #endif /*ALLOW_BULK_FORCE*/
219
220 return
221 end

  ViewVC Help
Powered by ViewVC 1.1.22