1 |
C $Header: /u/gcmpack/MITgcm/pkg/thsice/thsice_get_bulkf.F,v 1.6 2010/10/10 20:15:06 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "THSICE_OPTIONS.h" |
5 |
#ifdef ALLOW_BULK_FORCE |
6 |
#include "BULK_FORCE_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
CBOP |
10 |
C !ROUTINE: THSICE_GET_BULKF |
11 |
C !INTERFACE: |
12 |
SUBROUTINE THSICE_GET_BULKF( |
13 |
I bi, bj, |
14 |
I iMin,iMax, jMin,jMax, |
15 |
I iceFlag, hSnow, Tsf, |
16 |
O flxExcSw, dFlxdT, evap, dEvdT, |
17 |
I myTime, myIter, myThid ) |
18 |
C !DESCRIPTION: \bv |
19 |
C *==========================================================* |
20 |
C | S/R THSICE_GET_BULKF |
21 |
C *==========================================================* |
22 |
C | Interface S/R : get Surface Fluxes from pkg BULK_FORCE |
23 |
C *==========================================================* |
24 |
C \ev |
25 |
|
26 |
C !USES: |
27 |
IMPLICIT NONE |
28 |
|
29 |
C == Global data == |
30 |
#include "SIZE.h" |
31 |
#ifdef ALLOW_BULK_FORCE |
32 |
#include "EEPARAMS.h" |
33 |
#include "BULKF_PARAMS.h" |
34 |
#include "BULKF.h" |
35 |
#endif |
36 |
|
37 |
C !INPUT/OUTPUT PARAMETERS: |
38 |
C === Routine arguments === |
39 |
C bi,bj :: tile indices |
40 |
C iMin,iMax :: computation domain: 1rst index range |
41 |
C jMin,jMax :: computation domain: 2nd index range |
42 |
C iceFlag :: sea-ice fractional mask [0-1] |
43 |
C iceFlag :: True= get fluxes at this location ; False= do nothing |
44 |
C hSnow :: snow height [m] |
45 |
C Tsf :: surface (ice or snow) temperature (oC) |
46 |
C flxExcSw :: net (downward) surface heat flux, except short-wave [W/m2] |
47 |
C dFlxdT :: deriv of flx with respect to Tsf [W/m/K] |
48 |
C evap :: surface evaporation (>0 if evaporate) [kg/m2/s] |
49 |
C dEvdT :: deriv of evap. with respect to Tsf [kg/m2/s/K] |
50 |
C myThid :: Thread no. that called this routine. |
51 |
INTEGER bi, bj |
52 |
INTEGER iMin, iMax |
53 |
INTEGER jMin, jMax |
54 |
LOGICAL iceFlag (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
55 |
_RL hSnow (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
56 |
_RL Tsf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
57 |
_RL flxExcSw(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
58 |
_RL dFlxdT (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
59 |
_RL evap (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
60 |
_RL dEvdT (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
61 |
_RL myTime |
62 |
INTEGER myIter |
63 |
INTEGER myThid |
64 |
CEOP |
65 |
|
66 |
#ifdef ALLOW_THSICE |
67 |
#ifdef ALLOW_BULK_FORCE |
68 |
|
69 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
70 |
C === Local variables === |
71 |
C iceornot :: 0=open water, 1=ice cover, 2=ice+snow |
72 |
INTEGER iceornot |
73 |
INTEGER i, j |
74 |
_RL flwup ! upward LW at surface (W m-2) |
75 |
_RL flwNet_dwn ! net (downward) LW at surface (W m-2) |
76 |
_RL fsh ! surface downward sensible heat (W m-2) |
77 |
_RL flh ! surface downward latent heat (W m-2) |
78 |
_RL ust, vst, ssq |
79 |
#ifdef ALLOW_FORMULA_AIM |
80 |
_RL Tsurf(1), SHF(1), EVPloc(1), SLRU(1) |
81 |
_RL dEvp(1), sFlx(0:2) |
82 |
#endif |
83 |
|
84 |
DO j=jMin,jMax |
85 |
DO i=iMin,iMax |
86 |
IF ( iceFlag(i,j) ) THEN |
87 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
88 |
|
89 |
IF ( hSnow(i,j).GT.3. _d -1 ) THEN |
90 |
iceornot=2 |
91 |
ELSE |
92 |
iceornot=1 |
93 |
ENDIF |
94 |
|
95 |
#ifdef ALLOW_FORMULA_AIM |
96 |
IF ( useFluxFormula_AIM ) THEN |
97 |
|
98 |
Tsurf(1) = Tsf(i,j) |
99 |
CALL BULKF_FORMULA_AIM( |
100 |
I Tsurf, flwdwn(i,j,bi,bj), |
101 |
I ThAir(i,j,bi,bj), Tair(i,j,bi,bj), |
102 |
I Qair(i,j,bi,bj), wspeed(i,j,bi,bj), |
103 |
O SHF, EVPloc, SLRU, |
104 |
O dEvp, sFlx, |
105 |
I iceornot, myThid ) |
106 |
|
107 |
flxExcSw(i,j) = sFlx(1) |
108 |
dFlxdT(i,j) = sFlx(2) |
109 |
C- convert from [g/m2/s] to [kg/m2/s] |
110 |
evap(i,j) = EVPloc(1) * 1. _d -3 |
111 |
dEvdT(i,j) = dEvp(1) * 1. _d -3 |
112 |
|
113 |
ELSE |
114 |
#else /* ALLOW_FORMULA_AIM */ |
115 |
IF ( .TRUE. ) THEN |
116 |
#endif /* ALLOW_FORMULA_AIM */ |
117 |
|
118 |
ust = 0. |
119 |
vst = 0. |
120 |
ssq = 0. |
121 |
|
122 |
IF ( blk_nIter.EQ.0 ) THEN |
123 |
CALL BULKF_FORMULA_LANL( |
124 |
I uwind(i,j,bi,bj), vwind(i,j,bi,bj), wspeed(i,j,bi,bj), |
125 |
I Tair(i,j,bi,bj), Qair(i,j,bi,bj), |
126 |
I cloud(i,j,bi,bj), Tsf(i,j), |
127 |
O flwup, flh, fsh, dFlxdT(i,j), ust, vst, |
128 |
O evap(i,j), ssq, dEvdT(i,j), |
129 |
I iceornot, myThid ) |
130 |
ELSE |
131 |
CALL BULKF_FORMULA_LAY( |
132 |
I uwind(i,j,bi,bj), vwind(i,j,bi,bj), wspeed(i,j,bi,bj), |
133 |
I Tair(i,j,bi,bj), Qair(i,j,bi,bj), Tsf(i,j), |
134 |
O flwup, flh, fsh, dFlxdT(i,j), ust, vst, |
135 |
O evap(i,j), ssq, dEvdT(i,j), |
136 |
I iceornot, i,j,bi,bj,myThid ) |
137 |
ENDIF |
138 |
|
139 |
flwNet_dwn = flwdwn(i,j,bi,bj) - flwup |
140 |
flxExcSw(i,j) = flwNet_dwn + fsh + flh |
141 |
|
142 |
ENDIF |
143 |
|
144 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
145 |
ENDIF |
146 |
ENDDO |
147 |
ENDDO |
148 |
|
149 |
#endif /* ALLOW_BULK_FORCE */ |
150 |
#endif /* ALLOW_THSICE */ |
151 |
|
152 |
RETURN |
153 |
END |