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

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

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


Revision 1.3 - (hide annotations) (download)
Wed Dec 11 14:23:35 2002 UTC (21 years, 6 months ago) by cheisey
Branch: MAIN
CVS Tags: checkpoint47e_post, checkpoint47c_post, checkpoint50c_post, checkpoint48e_post, checkpoint50c_pre, checkpoint48i_post, checkpoint50d_pre, checkpoint51, checkpoint50, checkpoint50d_post, checkpoint50b_pre, checkpoint51f_post, checkpoint48b_post, checkpoint51d_post, checkpoint48c_pre, checkpoint47d_pre, checkpoint48d_pre, checkpoint51j_post, checkpoint47i_post, checkpoint47d_post, checkpoint48d_post, checkpoint48f_post, checkpoint48h_post, checkpoint51b_pre, checkpoint47g_post, checkpoint51h_pre, checkpoint48a_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint47j_post, branch-exfmods-tag, branchpoint-genmake2, checkpoint48c_post, checkpoint51b_post, checkpoint51c_post, checkpoint50g_post, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51i_pre, checkpoint47f_post, checkpoint50e_post, checkpoint51e_post, checkpoint48, checkpoint49, checkpoint51f_pre, checkpoint47h_post, checkpoint51g_post, checkpoint50b_post, checkpoint51a_post, checkpoint48g_post
Branch point for: branch-exfmods-curt, branch-genmake2
Changes since 1.2: +29 -33 lines
Tidying up some recent bulk_forcing changes.

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

  ViewVC Help
Powered by ViewVC 1.1.22