1 |
C $Header: /u/gcmpack/MITgcm/pkg/bulk_force/bulkf_readparms.F,v 1.7 2014/05/27 23:41:28 jmc Exp $ |
2 |
C $Name: BASE $ |
3 |
|
4 |
#include "BULK_FORCE_OPTIONS.h" |
5 |
|
6 |
SUBROUTINE BULKF_READPARMS( myThid ) |
7 |
C *==========================================================* |
8 |
C | SUBROUTINE BULKF_READPARMS | |
9 |
C | o Routine to initialize BULKF variables and constants. | |
10 |
C *==========================================================* |
11 |
C | Initialize BULKF parameters, read in data.blk | |
12 |
C *==========================================================* |
13 |
IMPLICIT NONE |
14 |
|
15 |
C === Global variables === |
16 |
#include "SIZE.h" |
17 |
#include "EEPARAMS.h" |
18 |
#include "PARAMS.h" |
19 |
#include "BULKF_PARAMS.h" |
20 |
#include "BULKF.h" |
21 |
#ifdef CONSERV_BULKF |
22 |
#include "BULKF_CONSERV.h" |
23 |
#endif |
24 |
|
25 |
C === Routine arguments === |
26 |
INTEGER myThid |
27 |
|
28 |
#ifdef ALLOW_BULK_FORCE |
29 |
C === Local variables === |
30 |
C msgBuf :: Informational/error message buffer |
31 |
C iUnit :: Work variable for IO unit number |
32 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
33 |
INTEGER iUnit |
34 |
|
35 |
C-- Bulk Formula parameter |
36 |
NAMELIST /BULKF_CONST/ |
37 |
& rhoA, rhoFW, |
38 |
& cpAir, Lvap, Lfresh, |
39 |
& Tf0kel, Rgas, xkar, stefan, |
40 |
& zref, zwd, zth, |
41 |
& cDrag_1, cDrag_2, cDrag_3, |
42 |
& cStantonS, cStantonU, cDalton, |
43 |
& umin, humid_fac, saltQsFac, gamma_blk, |
44 |
& atm_emissivity, ocean_emissivity, |
45 |
& snow_emissivity, ice_emissivity, |
46 |
#ifdef ALLOW_FORMULA_AIM |
47 |
& FWIND0, CHS, VGUST, DTHETA, dTstab, FSTAB, |
48 |
#endif |
49 |
& ocean_albedo |
50 |
|
51 |
NAMELIST /BULKF_PARM01/ |
52 |
& useFluxFormula_AIM, blk_nIter, calcWindStress, |
53 |
& blk_taveFreq, |
54 |
& AirTempFile, AirHumidityFile, RainFile, |
55 |
& SolarFile, LongwaveFile, UWindFile, |
56 |
& VWindFile, RunoffFile, WSpeedFile, QnetFile, |
57 |
& EmPFile, CloudFile, airPotTempFile |
58 |
|
59 |
#ifdef CONSERV_BULKF |
60 |
c- conserving qnet, empmr |
61 |
NAMELIST /BULKF_PARM02/ |
62 |
& qnet_off, empmr_off, conservcycle |
63 |
#endif |
64 |
|
65 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
66 |
|
67 |
IF ( .NOT.useBulkForce ) THEN |
68 |
C- pkg BULK_FORCE is not used |
69 |
_BEGIN_MASTER(myThid) |
70 |
C- Track pkg activation status: |
71 |
C print a (weak) warning if data.blk is found |
72 |
CALL PACKAGES_UNUSED_MSG( |
73 |
I 'useBulkForce', 'BULKF_READPARMS', 'blk' ) |
74 |
_END_MASTER(myThid) |
75 |
RETURN |
76 |
ENDIF |
77 |
|
78 |
_BEGIN_MASTER(myThid) |
79 |
|
80 |
CALL OPEN_COPY_DATA_FILE( |
81 |
I 'data.blk', 'BULKF_READPARMS', |
82 |
O iUnit, |
83 |
I myThid ) |
84 |
|
85 |
C-- Default values |
86 |
C- Physical constant : |
87 |
c slp0 = atm_Po / 100. ! reference sea-level atmospheric pressure [mb] |
88 |
rhoA = 1.3 _d 0 |
89 |
rhoFW = rhoConstFresh |
90 |
cpAir = atm_Cp |
91 |
c cpwv = 1.81 _d 3 |
92 |
Lvap = 2.5 _d 6 |
93 |
Lfresh = 3.34 _d 5 |
94 |
c Lvap_ice = 2.83 _d 6 |
95 |
Tf0kel = celsius2K |
96 |
Rgas = atm_Rd |
97 |
c Rvap = 461. _d 0 |
98 |
xkar = 0.4 _d 0 |
99 |
stefan = 5.67 _d -8 |
100 |
zref = 10.0 _d 0 |
101 |
zwd = zref |
102 |
zth = zref |
103 |
cDrag_1 = 2.70 _d -3 |
104 |
cDrag_2 = 0.142 _d -3 |
105 |
cDrag_3 = 0.0764 _d -3 |
106 |
cStantonS = 18.0 _d -3 |
107 |
cStantonU = 32.7 _d -3 |
108 |
cDalton = 34.6 _d -3 |
109 |
umin = 1.0 _d 0 |
110 |
humid_fac = 0.606 _d 0 |
111 |
saltQsFac = 0.980 _d 0 |
112 |
gamma_blk = 0.010 _d 0 |
113 |
atm_emissivity = .90 _d 0 |
114 |
ocean_emissivity= .985 _d 0 |
115 |
snow_emissivity = .98 _d 0 |
116 |
ice_emissivity = .98 _d 0 |
117 |
ocean_albedo = .10 _d 0 |
118 |
#ifdef ALLOW_FORMULA_AIM |
119 |
FWIND0 = 0.6 _d 0 |
120 |
CHS = 0.8 _d -3 |
121 |
VGUST = 5. _d 0 |
122 |
DTHETA = 3. _d 0 |
123 |
dTstab = 1. _d 0 |
124 |
FSTAB = 0.67 _d 0 |
125 |
#endif |
126 |
|
127 |
C- bulk-forcing parameters: |
128 |
useFluxFormula_AIM = .FALSE. |
129 |
blk_nIter = 5 |
130 |
calcWindStress = zonalWindFile .EQ. ' ' |
131 |
& .AND. meridWindFile .EQ. ' ' |
132 |
blk_taveFreq = taveFreq |
133 |
|
134 |
C- Input data files names : |
135 |
AirTempFile=' ' |
136 |
AirHumidityFile=' ' |
137 |
RainFile=' ' |
138 |
SolarFile=' ' |
139 |
LongwaveFile=' ' |
140 |
UWindFile=' ' |
141 |
VWindFile=' ' |
142 |
WspeedFile=' ' |
143 |
RunoffFile=' ' |
144 |
QnetFile=' ' |
145 |
EmPFile=' ' |
146 |
CloudFile=' ' |
147 |
SnowFile=' ' |
148 |
airPotTempFile=' ' |
149 |
|
150 |
C-- Read parameters from open data file |
151 |
WRITE(msgBuf,'(A)')' BULKF_READPARMS: starts to read BULKF_CONST' |
152 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
153 |
& SQUEEZE_RIGHT , myThid) |
154 |
READ(UNIT=iUnit,NML=BULKF_CONST) |
155 |
WRITE(msgBuf,'(A)') ' BULKF_READPARMS: read BULKF_CONST : OK' |
156 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
157 |
& SQUEEZE_RIGHT , myThid) |
158 |
|
159 |
WRITE(msgBuf,'(A)')' BULKF_READPARMS: starts to read BULKF_PARM01' |
160 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
161 |
& SQUEEZE_RIGHT , myThid) |
162 |
READ(UNIT=iUnit,NML=BULKF_PARM01) |
163 |
WRITE(msgBuf,'(A)') ' BULKF_READPARMS: read BULKF_PARM01 : OK' |
164 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
165 |
& SQUEEZE_RIGHT , myThid) |
166 |
|
167 |
#ifdef CONSERV_BULKF |
168 |
c -- default |
169 |
qnet_off=0.d0 |
170 |
empmr_off=0.d0 |
171 |
WRITE(msgBuf,'(A)')' BULKF_READPARMS: starts reading BULKF_PARM02' |
172 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
173 |
& SQUEEZE_RIGHT , myThid) |
174 |
READ(UNIT=iUnit,NML=BULKF_PARM02) |
175 |
WRITE(msgBuf,'(A)') ' BULKF_READPARMS: read BULKF_PARM02 : OK' |
176 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
177 |
& SQUEEZE_RIGHT , myThid) |
178 |
|
179 |
#endif /* CONSERV_BULKF */ |
180 |
|
181 |
C-- Close the open data file |
182 |
#ifdef SINGLE_DISK_IO |
183 |
CLOSE(iUnit) |
184 |
#else |
185 |
CLOSE(iUnit,STATUS='DELETE') |
186 |
#endif /* SINGLE_DISK_IO */ |
187 |
|
188 |
C- check that CPP option is "defined" when running-flag parameter is on: |
189 |
#ifndef ALLOW_FORMULA_AIM |
190 |
IF ( useFluxFormula_AIM ) THEN |
191 |
WRITE(msgBuf,'(2A)') ' BULKF_READPARMS: ', |
192 |
& 'useFluxFormula_AIM is TRUE and #undef ALLOW_FORMULA_AIM' |
193 |
CALL PRINT_ERROR( msgBuf , myThid) |
194 |
WRITE(msgBuf,'(2A)') ' BULKF_READPARMS: => recompile with', |
195 |
& ' #define ALLOW_FORMULA_AIM in BULK_FORCE_OPTIONS.h' |
196 |
CALL PRINT_ERROR( msgBuf , myThid) |
197 |
STOP 'ABNORMAL END: S/R CONFIG_CHECK' |
198 |
ENDIF |
199 |
#endif |
200 |
|
201 |
C- Define other constants (from previous ones): |
202 |
c Qcoef = 6.11 _d 0 * 0.622 _d 0 / p0 |
203 |
c Sha = Rgas / .286 _d 0 |
204 |
|
205 |
useQnetch = QnetFile .NE. ' ' |
206 |
useEmPch = EmPFile .NE. ' ' |
207 |
|
208 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
209 |
iUnit = standardMessageUnit |
210 |
c iUnit=88 |
211 |
c OPEN(iUnit,file='bulkf_check_params',status='unknown') |
212 |
c WRITE(iUnit,*) 'BlkF: slp0 =',slp0 |
213 |
WRITE(iUnit,*) 'BlkF: rhoA =',rhoA |
214 |
WRITE(iUnit,*) 'BlkF: rhoFW =',rhoFW |
215 |
WRITE(iUnit,*) 'BlkF: cpAir =',cpAir |
216 |
c WRITE(iUnit,*) 'BlkF: cpwv =',cpwv |
217 |
WRITE(iUnit,*) 'BlkF: Lvap =',Lvap |
218 |
WRITE(iUnit,*) 'BlkF: Lfresh =',Lfresh |
219 |
c WRITE(iUnit,*) 'BlkF: Lvap_ice =',Lvap_ice |
220 |
c WRITE(iUnit,*) 'BlkF: Sha =',Sha |
221 |
WRITE(iUnit,*) 'BlkF: Tf0kel =',Tf0kel |
222 |
WRITE(iUnit,*) 'BlkF: Rgas =',Rgas |
223 |
c WRITE(iUnit,*) 'BlkF: Rvap =',Rvap |
224 |
WRITE(iUnit,*) 'BlkF: xkar =',xkar |
225 |
WRITE(iUnit,*) 'BlkF: stefan =',stefan |
226 |
WRITE(iUnit,*) 'BlkF: zref =',zref |
227 |
WRITE(iUnit,*) 'BlkF: zwd =',zwd |
228 |
WRITE(iUnit,*) 'BlkF: zth =',zth |
229 |
WRITE(iUnit,*) 'BlkF: cDrag_1 =',cDrag_1 |
230 |
WRITE(iUnit,*) 'BlkF: cDrag_2 =',cDrag_2 |
231 |
WRITE(iUnit,*) 'BlkF: cDrag_3 =',cDrag_3 |
232 |
WRITE(iUnit,*) 'BlkF: cStantonS=',cStantonS |
233 |
WRITE(iUnit,*) 'BlkF: cStantonU=',cStantonU |
234 |
WRITE(iUnit,*) 'BlkF: cDalton =',cDalton |
235 |
WRITE(iUnit,*) 'BlkF: umin =',umin |
236 |
WRITE(iUnit,*) 'BlkF: humid_fac=',humid_fac |
237 |
WRITE(iUnit,*) 'BlkF: saltQsFac=',saltQsFac |
238 |
WRITE(iUnit,*) 'BlkF: gamma_blk=',gamma_blk |
239 |
WRITE(iUnit,*) 'BlkF: atm_emissivity =',atm_emissivity |
240 |
WRITE(iUnit,*) 'BlkF: ocean_emissivity=',ocean_emissivity |
241 |
WRITE(iUnit,*) 'BlkF: snow_emissivity =',snow_emissivity |
242 |
WRITE(iUnit,*) 'BlkF: ice_emissivity =',ice_emissivity |
243 |
WRITE(iUnit,*) 'BlkF: ocean_albedo =',ocean_albedo |
244 |
#ifdef ALLOW_FORMULA_AIM |
245 |
WRITE(iUnit,*) 'BlkF: FWIND0 =', FWIND0 |
246 |
WRITE(iUnit,*) 'BlkF: CHS =', CHS |
247 |
WRITE(iUnit,*) 'BlkF: VGUST =', VGUST |
248 |
WRITE(iUnit,*) 'BlkF: DTHETA =', DTHETA |
249 |
WRITE(iUnit,*) 'BlkF: dTstab =', dTstab |
250 |
WRITE(iUnit,*) 'BlkF: FSTAB =', FSTAB |
251 |
#endif |
252 |
WRITE(iUnit,*) 'BlkF: useFluxFormula_AIM=',useFluxFormula_AIM |
253 |
WRITE(iUnit,*) 'BlkF: calcWindStress =', calcWindStress |
254 |
WRITE(iUnit,*) 'BlkF: useQnetch =', useQnetch |
255 |
WRITE(iUnit,*) 'BlkF: useEmPch =', useEmPch |
256 |
WRITE(iUnit,*) 'BlkF: blk_nIter =',blk_nIter |
257 |
WRITE(iUnit,*) 'BlkF: blk_taveFreq=', blk_taveFreq |
258 |
IF (iUnit.EQ.88) CLOSE(iUnit) |
259 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
260 |
|
261 |
_END_MASTER(myThid) |
262 |
|
263 |
C-- Everyone else must wait for the parameters to be loaded |
264 |
_BARRIER |
265 |
|
266 |
#endif /* ALLOW_BULK_FORCE */ |
267 |
|
268 |
RETURN |
269 |
END |