1 |
C $Header: /u/gcmpack/MITgcm/pkg/profiles/profiles_init_fixed.F,v 1.24 2012/09/04 14:32:29 gforget Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PROFILES_OPTIONS.h" |
5 |
#include "AD_CONFIG.h" |
6 |
|
7 |
C *==========================================================* |
8 |
C | subroutine profiles_init_fixed |
9 |
C | o initialization for netcdf profiles data |
10 |
C | started: Gael Forget 15-March-2006 |
11 |
C | extended: Gael Forget 14-June-2007 |
12 |
C *==========================================================* |
13 |
|
14 |
SUBROUTINE profiles_init_fixed( myThid ) |
15 |
|
16 |
implicit none |
17 |
|
18 |
C ==================== Global Variables =========================== |
19 |
#include "SIZE.h" |
20 |
#include "EEPARAMS.h" |
21 |
#include "PARAMS.h" |
22 |
#include "GRID.h" |
23 |
#include "DYNVARS.h" |
24 |
#ifdef ALLOW_CAL |
25 |
#include "cal.h" |
26 |
#endif |
27 |
#ifdef ALLOW_PROFILES |
28 |
# include "profiles.h" |
29 |
# include "netcdf.inc" |
30 |
#endif |
31 |
C ==================== Routine Variables ========================== |
32 |
|
33 |
integer k,l,m,bi,bj,iG,jG, myThid,num_file,length_for_tile |
34 |
_RL stopProfiles |
35 |
integer fid, dimid, varid1, varid1a, varid1b |
36 |
integer varid2,varid3 |
37 |
_RL tmpyymmdd(1000),tmphhmmss(1000),diffsecs |
38 |
_RL yymmddMin,yymmddMax |
39 |
_RL hhmmssMin,hhmmssMax |
40 |
|
41 |
integer tmpdate(4),tmpdiff(4),profIsInRunTime |
42 |
_RL tmp_lon, tmp_lon2(1000), tmp_lat2(1000) |
43 |
integer vec_start(2), vec_count(2), profno_div1000, kk |
44 |
character*(80) profilesfile, fnamedatanc |
45 |
character*(80) fnameequinc, adfnameequinc |
46 |
integer IL, JL, err |
47 |
logical exst |
48 |
|
49 |
#ifdef ALLOW_PROFILES |
50 |
|
51 |
#ifdef ALLOW_PROFILES_GENERICGRID |
52 |
integer varid_intp1, varid_intp2, varid_intp11 , varid_intp22 |
53 |
integer varid_intp3, varid_intp4, varid_intp5, q |
54 |
_RL tmp_i(1000,NUM_INTERP_POINTS) |
55 |
_RL tmp_j(1000,NUM_INTERP_POINTS) |
56 |
_RL tmp_weights(1000,NUM_INTERP_POINTS),tmp_sum_weights |
57 |
_RL tmp_xC11(1000),tmp_yC11(1000) |
58 |
_RL tmp_xCNINJ(1000),tmp_yCNINJ(1000) |
59 |
_RL stopGenericGrid |
60 |
Real*8 xy_buffer_r8(0:sNx+1,0:sNy+1) |
61 |
integer vec_start2(2), vec_count2(2) |
62 |
#endif |
63 |
|
64 |
c == external functions == |
65 |
integer ILNBLNK |
66 |
EXTERNAL ILNBLNK |
67 |
integer MDS_RECLEN |
68 |
EXTERNAL MDS_RECLEN |
69 |
character*(max_len_mbuf) msgbuf |
70 |
|
71 |
c-- == end of interface == |
72 |
|
73 |
|
74 |
write(msgbuf,'(a)') ' ' |
75 |
call print_message( msgbuf, standardmessageunit, |
76 |
& SQUEEZE_RIGHT , mythid) |
77 |
write(msgbuf,'(a)') |
78 |
&'// =======================================================' |
79 |
call print_message( msgbuf, standardmessageunit, |
80 |
& SQUEEZE_RIGHT , mythid) |
81 |
write(msgbuf,'(a)') |
82 |
&'// insitu profiles model sampling >>> START <<<' |
83 |
call print_message( msgbuf, standardmessageunit, |
84 |
& SQUEEZE_RIGHT , mythid) |
85 |
write(msgbuf,'(a)') |
86 |
&'// =======================================================' |
87 |
call print_message( msgbuf, standardmessageunit, |
88 |
& SQUEEZE_RIGHT , mythid) |
89 |
write(msgbuf,'(a)') ' ' |
90 |
call print_message( msgbuf, standardmessageunit, |
91 |
& SQUEEZE_RIGHT , mythid) |
92 |
|
93 |
|
94 |
stopProfiles=0. _d 0 |
95 |
#ifdef ALLOW_PROFILES_GENERICGRID |
96 |
stopGenericGrid=0. _d 0 |
97 |
#endif |
98 |
|
99 |
#ifndef ALLOW_PROFILES_GENERICGRID |
100 |
IF ( profilesDoGenGrid ) THEN |
101 |
WRITE(msgBuf,'(2A)') 'PROFILES_INIT_FIXED: ', |
102 |
& 'profilesDoGenGrid=.true. requires' |
103 |
CALL PRINT_ERROR( msgBuf , myThid ) |
104 |
WRITE(msgBuf,'(2A)') 'PROFILES_INIT_FIXED: ', |
105 |
& 'that ALLOW_PROFILES_GENERICGRID is defined' |
106 |
CALL PRINT_ERROR( msgBuf , myThid ) |
107 |
STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED' |
108 |
ENDIF |
109 |
#endif |
110 |
|
111 |
IF ( (.NOT.profilesDoGenGrid).AND. |
112 |
& (.NOT.usingSphericalPolarGrid .OR. rotateGrid) ) THEN |
113 |
WRITE(msgBuf,'(2A)') 'PROFILES_INIT_FIXED: ', |
114 |
& 'profilesDoGenGrid=.true. is required' |
115 |
CALL PRINT_ERROR( msgBuf , myThid ) |
116 |
WRITE(msgBuf,'(2A)') 'PROFILES_INIT_FIXED: ', |
117 |
& 'unless usingSphericalGrid=.TRUE. and rotateGrid=.FALSE.' |
118 |
CALL PRINT_ERROR( msgBuf , myThid ) |
119 |
STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED' |
120 |
ENDIF |
121 |
|
122 |
|
123 |
write(msgbuf,'(a)') ' ' |
124 |
call print_message( msgbuf, standardmessageunit, |
125 |
& SQUEEZE_RIGHT , mythid) |
126 |
write(msgbuf,'(a)') 'general packages parameters :' |
127 |
JL = ILNBLNK( profilesDir ) |
128 |
if (JL.NE.0) then |
129 |
write(msgbuf,'(a,a)') ' profilesDir ',profilesDir(1:JL) |
130 |
else |
131 |
write(msgbuf,'(a,a)') ' profilesDir ','./' |
132 |
endif |
133 |
call print_message( |
134 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
135 |
write(msgbuf,'(a,a)') ' ALLOW_PROFILES_GENERICGRID ', |
136 |
#ifdef ALLOW_PROFILES_GENERICGRID |
137 |
& 'was compiled' |
138 |
#else |
139 |
& 'was NOT compiled' |
140 |
#endif |
141 |
call print_message( |
142 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
143 |
write(msgbuf,'(a,l5)') ' profilesDoGenGrid ',profilesDoGenGrid |
144 |
call print_message( |
145 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
146 |
write(msgbuf,'(a,l5)') ' profilesDoNcOutput ',profilesDoNcOutput |
147 |
call print_message( |
148 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
149 |
write(msgbuf,'(a)') ' ' |
150 |
call print_message( msgbuf, standardmessageunit, |
151 |
& SQUEEZE_RIGHT , mythid) |
152 |
|
153 |
|
154 |
_BEGIN_MASTER( mythid ) |
155 |
DO bj=1,nSy |
156 |
DO bi=1,nSx |
157 |
|
158 |
profiles_curfile_buff(bi,bj)=0 |
159 |
yymmddMin=modelstartdate(1) |
160 |
yymmddMax=modelenddate(1) |
161 |
hhmmssMin=modelstartdate(2) |
162 |
hhmmssMax=modelenddate(2) |
163 |
|
164 |
do m=1,NLEVELMAX |
165 |
do l=1,1000 |
166 |
do k=1,NVARMAX |
167 |
profiles_data_buff(m,l,k,bi,bj)=0 |
168 |
profiles_weight_buff(m,l,k,bi,bj)=0 |
169 |
enddo |
170 |
enddo |
171 |
enddo |
172 |
|
173 |
do num_file=1,NFILESPROFMAX |
174 |
|
175 |
write(msgbuf,'(a)') ' ' |
176 |
call print_message( msgbuf, standardmessageunit, |
177 |
& SQUEEZE_RIGHT , mythid) |
178 |
|
179 |
IL = ILNBLNK( profilesfiles(num_file) ) |
180 |
if (IL.NE.0) then |
181 |
write(profilesfile(1:80),'(1a)') |
182 |
& profilesfiles(num_file)(1:IL) |
183 |
write(msgbuf,'(a,i3,a,a)') |
184 |
& 'profiles file ',num_file,' is ', profilesfile(1:80) |
185 |
call print_message( |
186 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
187 |
else |
188 |
write(profilesfile(1:80),'(1a)') ' ' |
189 |
write(msgbuf,'(a,i3,a,a)') |
190 |
& 'profiles file ',num_file,' is ',' (empty) ' |
191 |
call print_message( |
192 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
193 |
endif |
194 |
|
195 |
IL = ILNBLNK( profilesfile ) |
196 |
if (IL.NE.0) then |
197 |
|
198 |
C=========================================================== |
199 |
c open data files and read information |
200 |
C=========================================================== |
201 |
|
202 |
write(fnamedatanc(1:80),'(2a)') profilesfile(1:IL),'.nc' |
203 |
write(msgbuf,'(a,a)') ' opening ', fnamedatanc(1:80) |
204 |
call print_message( |
205 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
206 |
err = NF_OPEN(fnamedatanc, 0, fiddata(num_file,bi,bj)) |
207 |
|
208 |
c1) read the number of profiles : |
209 |
fid=fiddata(num_file,bi,bj) |
210 |
err = NF_INQ_DIMID(fid,'iPROF', dimid ) |
211 |
err = NF_INQ_DIMLEN(fid, dimid, ProfNo(num_file,bi,bj) ) |
212 |
err = NF_INQ_DIMID(fid,'iDEPTH', dimid ) |
213 |
if (err.NE.NF_NOERR) then |
214 |
err = NF_INQ_DIMID(fid,'Z', dimid ) |
215 |
endif |
216 |
err = NF_INQ_DIMLEN(fid, dimid, ProfDepthNo(num_file,bi,bj) ) |
217 |
|
218 |
write(msgbuf,'(a,i9)') |
219 |
& ' file ID is ', fid |
220 |
call print_message( |
221 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
222 |
write(msgbuf,'(a,i5)') |
223 |
& ' no. of depth levels ',ProfDepthNo(num_file,bi,bj) |
224 |
call print_message( |
225 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
226 |
write(msgbuf,'(a,i9)') |
227 |
& ' no. of profiles ', ProfNo(num_file,bi,bj) |
228 |
call print_message( |
229 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
230 |
|
231 |
c2) read the dates and positions : |
232 |
err = NF_INQ_VARID(fid,'prof_depth', varid1a ) |
233 |
if (err.NE.NF_NOERR) then |
234 |
c if no prof_depth is found, then try old variable name: |
235 |
err = NF_INQ_VARID(fid,'depth', varid1a ) |
236 |
endif |
237 |
if (err.NE.NF_NOERR) then |
238 |
c if neither is found, then stop |
239 |
IL = ILNBLNK( profilesfile ) |
240 |
WRITE(msgBuf,'(3A)') |
241 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
242 |
& '.nc is not in the pkg/profiles format (no prof_depth etc.)' |
243 |
CALL PRINT_ERROR( msgBuf, myThid) |
244 |
stopProfiles=1. _d 0 |
245 |
endif |
246 |
|
247 |
do k=1,ProfDepthNo(num_file,bi,bj) |
248 |
err = NF_GET_VAR1_DOUBLE(fid,varid1a,k, |
249 |
& prof_depth(num_file,k,bi,bj)) |
250 |
enddo |
251 |
|
252 |
err = NF_INQ_VARID(fid,'prof_YYYYMMDD', varid1a ) |
253 |
err = NF_INQ_VARID(fid,'prof_HHMMSS', varid1b ) |
254 |
err = NF_INQ_VARID(fid,'prof_lon', varid2 ) |
255 |
err = NF_INQ_VARID(fid,'prof_lat', varid3 ) |
256 |
|
257 |
if (err.NE.NF_NOERR) then |
258 |
IL = ILNBLNK( profilesfile ) |
259 |
WRITE(msgBuf,'(3A)') |
260 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
261 |
& '.nc is not in the pkg/profiles format (no prof_YYYYMMDD etc.)' |
262 |
CALL PRINT_ERROR( msgBuf, myThid) |
263 |
stopProfiles=1. _d 0 |
264 |
endif |
265 |
|
266 |
#ifdef ALLOW_PROFILES_GENERICGRID |
267 |
if (profilesDoGenGrid) then |
268 |
c3) read interpolattion information (grid points, coeffs, etc.) |
269 |
err = NF_INQ_VARID(fid,'prof_interp_XC11',varid_intp1) |
270 |
err = NF_INQ_VARID(fid,'prof_interp_YC11',varid_intp2) |
271 |
err = NF_INQ_VARID(fid,'prof_interp_XCNINJ',varid_intp11) |
272 |
err = NF_INQ_VARID(fid,'prof_interp_YCNINJ',varid_intp22) |
273 |
err = NF_INQ_VARID(fid,'prof_interp_weights',varid_intp3) |
274 |
err = NF_INQ_VARID(fid,'prof_interp_i',varid_intp4) |
275 |
err = NF_INQ_VARID(fid,'prof_interp_j',varid_intp5) |
276 |
if (err.NE.NF_NOERR) then |
277 |
IL = ILNBLNK( profilesfile ) |
278 |
WRITE(msgBuf,'(3A)') |
279 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
280 |
& '.nc is missing interpolation information (profilesDoGenGrid)' |
281 |
CALL PRINT_ERROR( msgBuf, myThid) |
282 |
stopGenericGrid=2. _d 0 |
283 |
endif |
284 |
endif |
285 |
#endif |
286 |
|
287 |
|
288 |
c4) default values |
289 |
do k=1,NOBSGLOB |
290 |
prof_time(num_file,k,bi,bj)=-999 |
291 |
prof_lon(num_file,k,bi,bj)=-999 |
292 |
prof_lat(num_file,k,bi,bj)=-999 |
293 |
prof_ind_glob(num_file,k,bi,bj)=-999 |
294 |
#ifdef ALLOW_PROFILES_GENERICGRID |
295 |
do q = 1,NUM_INTERP_POINTS |
296 |
prof_interp_i(num_file,k,q,bi,bj) = -999 |
297 |
prof_interp_j(num_file,k,q,bi,bj) = -999 |
298 |
prof_interp_weights(num_file,k,q,bi,bj) = -999 |
299 |
enddo |
300 |
prof_interp_xC11(num_file,k,bi,bj)=-999 |
301 |
prof_interp_yC11(num_file,k,bi,bj)=-999 |
302 |
prof_interp_xCNINJ(num_file,k,bi,bj)=-999 |
303 |
prof_interp_yCNINJ(num_file,k,bi,bj)=-999 |
304 |
#endif |
305 |
enddo |
306 |
|
307 |
|
308 |
c5) main loop: look for profiles in this tile |
309 |
length_for_tile=0 |
310 |
profno_div1000=max(0,int(ProfNo(num_file,bi,bj)/1000)) |
311 |
|
312 |
do kk=1,profno_div1000+1 |
313 |
|
314 |
if (min(ProfNo(num_file,bi,bj), 1000*kk).GE. |
315 |
& 1+1000*(kk-1)) then |
316 |
|
317 |
c5.1) read a chunk |
318 |
vec_start(1)=1 |
319 |
vec_start(2)=1+1000*(kk-1) |
320 |
vec_count(1)=1 |
321 |
vec_count(2)=min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
322 |
|
323 |
if ( (vec_count(2).LE.0).OR.(vec_count(2).GT.1000).OR. |
324 |
& (vec_start(2).LE.0).OR. |
325 |
& (vec_count(2)+vec_start(2)-1.GT.ProfNo(num_file,bi,bj)) ) |
326 |
& then |
327 |
IL = ILNBLNK( profilesfile ) |
328 |
WRITE(msgBuf,'(3A)') |
329 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
330 |
& '.nc was not read properly (case 1).' |
331 |
CALL PRINT_ERROR( msgBuf, myThid) |
332 |
stopProfiles=1. _d 0 |
333 |
endif |
334 |
|
335 |
err = NF_GET_VARA_DOUBLE(fid,varid1a,vec_start(2), |
336 |
& vec_count(2), tmpyymmdd) |
337 |
err = NF_GET_VARA_DOUBLE(fid,varid1b,vec_start(2), |
338 |
& vec_count(2), tmphhmmss) |
339 |
err = NF_GET_VARA_DOUBLE(fid,varid2,vec_start(2), |
340 |
& vec_count(2), tmp_lon2) |
341 |
err = NF_GET_VARA_DOUBLE(fid,varid3,vec_start(2), |
342 |
& vec_count(2), tmp_lat2) |
343 |
|
344 |
if (err.NE.NF_NOERR) then |
345 |
WRITE(msgBuf,'(3A)') |
346 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
347 |
& '.nc was not read properly (case 2).' |
348 |
CALL PRINT_ERROR( msgBuf, myThid) |
349 |
stopProfiles=1. _d 0 |
350 |
endif |
351 |
|
352 |
#ifdef ALLOW_PROFILES_GENERICGRID |
353 |
if (profilesDoGenGrid) then |
354 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp1,vec_start(2), |
355 |
& vec_count(2), tmp_xC11) |
356 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp2,vec_start(2), |
357 |
& vec_count(2), tmp_yC11) |
358 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp11,vec_start(2), |
359 |
& vec_count(2), tmp_xCNINJ) |
360 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp22,vec_start(2), |
361 |
& vec_count(2), tmp_yCNINJ) |
362 |
do q=1,NUM_INTERP_POINTS |
363 |
vec_start2(1)=q |
364 |
vec_start2(2)=1+1000*(kk-1) |
365 |
vec_count2(1)=1 |
366 |
vec_count2(2)=min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
367 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp3,vec_start2, |
368 |
& vec_count2, tmp_weights(1,q)) |
369 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp4,vec_start2, |
370 |
& vec_count2, tmp_i(1,q)) |
371 |
err = NF_GET_VARA_DOUBLE(fid,varid_intp5,vec_start2, |
372 |
& vec_count2, tmp_j(1,q)) |
373 |
enddo |
374 |
endif |
375 |
#endif |
376 |
|
377 |
c5.2) loop through this chunk |
378 |
do k=1,min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
379 |
|
380 |
if ( stopProfiles .EQ. 0.) then |
381 |
|
382 |
profIsInRunTime=1 |
383 |
|
384 |
if (( ( tmpyymmdd(k).GT.yymmddMin ).OR.(( tmpyymmdd(k).EQ. |
385 |
& yymmddMin ).AND.( tmphhmmss(k).GT.hhmmssMin ))).AND. |
386 |
& ( ( tmpyymmdd(k).LT.yymmddMax ).OR.(( tmpyymmdd(k).EQ. |
387 |
& yymmddMax ).AND.( tmphhmmss(k).LT.hhmmssMax ))) ) then |
388 |
profIsInRunTime=1 |
389 |
call cal_FullDate( int(tmpyymmdd(k)),int(tmphhmmss(k)), |
390 |
& tmpdate,mythid ) |
391 |
call cal_TimePassed( modelstartdate,tmpdate,tmpdiff,mythid ) |
392 |
call cal_ToSeconds (tmpdiff,diffsecs,mythid) |
393 |
diffsecs=diffsecs+nIter0*deltaTclock |
394 |
else |
395 |
profIsInRunTime=0 |
396 |
diffsecs=-deltaTclock |
397 |
endif |
398 |
|
399 |
#ifdef ALLOW_PROFILES_GENERICGRID |
400 |
if (.NOT.profilesDoGenGrid) then |
401 |
#endif |
402 |
if (xC(sNx+1,1,bi,bj).LT.xC(1,1,bi,bj)) then |
403 |
tmp_lon=xC(sNx+1,1,bi,bj)+360 |
404 |
else |
405 |
tmp_lon=xC(sNx+1,1,bi,bj) |
406 |
endif |
407 |
if ((xC(1,1,bi,bj).LE.tmp_lon2(k)).AND. |
408 |
& (tmp_lon.GT.tmp_lon2(k)).AND. |
409 |
& (yC(1,1,bi,bj).LE.tmp_lat2(k)).AND. |
410 |
& (yC(1,sNy+1,bi,bj).GT.tmp_lat2(k)) |
411 |
& .AND.(profIsInRunTime.EQ.1)) then |
412 |
length_for_tile=length_for_tile+1 |
413 |
prof_time(num_file,length_for_tile,bi,bj)=diffsecs |
414 |
prof_lon(num_file,length_for_tile,bi,bj)=tmp_lon2(k) |
415 |
prof_lat(num_file,length_for_tile,bi,bj)=tmp_lat2(k) |
416 |
prof_ind_glob(num_file,length_for_tile,bi,bj)=k+1000*(kk-1) |
417 |
if (length_for_tile.EQ.NOBSGLOB) then |
418 |
WRITE(msgBuf,'(3A)') |
419 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
420 |
& '.nc was not read properly (increase NOBSGLOB).' |
421 |
CALL PRINT_ERROR( msgBuf, myThid) |
422 |
stopProfiles=1. _d 0 |
423 |
endif |
424 |
elseif (xC(sNx+1,1,bi,bj).LT.xC(1,1,bi,bj)) then |
425 |
if ((xC(1,1,bi,bj).LE.tmp_lon2(k)+360).AND. |
426 |
& (tmp_lon.GT.tmp_lon2(k)+360).AND. |
427 |
& (yC(1,1,bi,bj).LE.tmp_lat2(k)).AND. |
428 |
& (yC(1,sNy+1,bi,bj).GT.tmp_lat2(k)) |
429 |
& ) then |
430 |
length_for_tile=length_for_tile+1 |
431 |
prof_time(num_file,length_for_tile,bi,bj)=diffsecs |
432 |
prof_lon(num_file,length_for_tile,bi,bj)=tmp_lon2(k)+360 |
433 |
prof_lat(num_file,length_for_tile,bi,bj)=tmp_lat2(k) |
434 |
prof_ind_glob(num_file,length_for_tile,bi,bj)=k+1000*(kk-1) |
435 |
if (length_for_tile.EQ.NOBSGLOB) then |
436 |
WRITE(msgBuf,'(3A)') |
437 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
438 |
& '.nc was not read properly (increase NOBSGLOB).' |
439 |
CALL PRINT_ERROR( msgBuf, myThid) |
440 |
stopProfiles=1. _d 0 |
441 |
endif |
442 |
endif |
443 |
endif |
444 |
#ifdef ALLOW_PROFILES_GENERICGRID |
445 |
else |
446 |
if (stopGenericGrid.EQ.0.) then |
447 |
|
448 |
if ( ( abs( tmp_xC11(k) - xC(1,1,bi,bj) ).LT.0.0001 ) .AND. |
449 |
& ( abs( tmp_yC11(k) - yC(1,1,bi,bj) ).LT.0.0001 ) .AND. |
450 |
& ( abs( tmp_xCNINJ(k) - xC(sNx,sNy,bi,bj) ).LT.0.0001 ) .AND. |
451 |
& ( abs( tmp_yCNINJ(k) - yC(sNx,sNy,bi,bj) ).LT.0.0001 ) |
452 |
& .AND.(profIsInRunTime.EQ.1)) then |
453 |
|
454 |
length_for_tile=length_for_tile+1 |
455 |
prof_time(num_file,length_for_tile,bi,bj)=diffsecs |
456 |
prof_interp_xC11(num_file,length_for_tile,bi,bj)=tmp_xC11(k) |
457 |
prof_interp_yC11(num_file,length_for_tile,bi,bj)=tmp_yC11(k) |
458 |
prof_interp_xCNINJ(num_file,length_for_tile,bi,bj)=tmp_xCNINJ(k) |
459 |
prof_interp_yCNINJ(num_file,length_for_tile,bi,bj)=tmp_yCNINJ(k) |
460 |
tmp_sum_weights=0. _d 0 |
461 |
do q = 1,NUM_INTERP_POINTS |
462 |
prof_interp_weights(num_file,length_for_tile,q,bi,bj) |
463 |
& =tmp_weights(k,q) |
464 |
prof_interp_i(num_file,length_for_tile,q,bi,bj) |
465 |
& =tmp_i(k,q) |
466 |
prof_interp_j(num_file,length_for_tile,q,bi,bj) |
467 |
& =tmp_j(k,q) |
468 |
tmp_sum_weights=tmp_sum_weights+tmp_weights(k,q) |
469 |
c more test of the inputs: is the offline-computed |
470 |
c interpolation information consistent (self and with grid) |
471 |
if ( (tmp_i(k,q).LT.0).OR.(tmp_j(k,q).LT.0) |
472 |
& .OR.(tmp_i(k,q).GT.sNx+1).OR.(tmp_j(k,q).GT.sNy+1) ) then |
473 |
WRITE(msgBuf,'(4A)') |
474 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
475 |
& '.nc includes inconsistent interpolation ', |
476 |
& 'points (profilesDoGenGrid; out of tile)' |
477 |
CALL PRINT_ERROR( msgBuf, myThid) |
478 |
stopGenericGrid=1. _d 0 |
479 |
endif |
480 |
if ( tmp_weights(k,q) .NE. 0. ) then |
481 |
if ( ((tmp_i(k,q).EQ.0).AND.(tmp_j(k,q).EQ.0)) |
482 |
& .OR.((tmp_i(k,q).EQ.sNx+1).AND.(tmp_j(k,q).EQ.sNy+1)) |
483 |
& .OR.((tmp_i(k,q).EQ.0).AND.(tmp_j(k,q).EQ.sNy+1)) |
484 |
& .OR.((tmp_i(k,q).EQ.sNx+1).AND.(tmp_j(k,q).EQ.0)) ) then |
485 |
WRITE(msgBuf,'(4A)') |
486 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
487 |
& '.nc includes inconsistent interpolation ', |
488 |
& 'points (profilesDoGenGrid; using overlap corners)' |
489 |
CALL PRINT_ERROR( msgBuf, myThid) |
490 |
stopGenericGrid=1. _d 0 |
491 |
endif |
492 |
endif |
493 |
if ( (tmp_weights(k,q).LT.0).OR.(tmp_weights(k,q).GT.1) ) then |
494 |
WRITE(msgBuf,'(4A)') |
495 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
496 |
& '.nc includes inconsistent interpolation ', |
497 |
& 'weights (profilesDoGenGrid; sum oustide 0-1)' |
498 |
CALL PRINT_ERROR( msgBuf, myThid) |
499 |
stopGenericGrid=1. _d 0 |
500 |
endif |
501 |
|
502 |
enddo |
503 |
|
504 |
if ( abs(tmp_sum_weights -1. ) .GT. 0.0001 ) then |
505 |
WRITE(msgBuf,'(4A)') |
506 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
507 |
& '.nc includes inconsistent interpolation ', |
508 |
& 'weights (profilesDoGenGrid; dont add up to 1)' |
509 |
CALL PRINT_ERROR( msgBuf, myThid) |
510 |
stopGenericGrid=1. _d 0 |
511 |
endif |
512 |
|
513 |
prof_ind_glob(num_file,length_for_tile,bi,bj)=k+1000*(kk-1) |
514 |
if (length_for_tile.EQ.NOBSGLOB) then |
515 |
WRITE(msgBuf,'(3A)') |
516 |
& 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL), |
517 |
& '.nc was not read properly (increase NOBSGLOB).' |
518 |
CALL PRINT_ERROR( msgBuf, myThid) |
519 |
stopProfiles=1. _d 0 |
520 |
endif |
521 |
|
522 |
endif |
523 |
endif |
524 |
endif !if (.NOT.profilesDoGenGrid) then |
525 |
#endif /* ALLOW_PROFILES_GENERICGRID */ |
526 |
endif !if ( stopProfiles .EQ. 0.) then |
527 |
enddo !do k=1,min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
528 |
endif !if (min(ProfNo(num_file,bi,bj), 1000... |
529 |
enddo !do kk=1,profno_div1000+1 |
530 |
|
531 |
ProfNo(num_file,bi,bj)=length_for_tile |
532 |
|
533 |
write(msgbuf,'(a,i9)') |
534 |
& ' within tile & run ', ProfNo(num_file,bi,bj) |
535 |
call print_message( |
536 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
537 |
|
538 |
c6) available variablesin the data set |
539 |
|
540 |
do k=1,NVARMAX |
541 |
prof_num_var_cur(num_file,k,bi,bj)=0 |
542 |
enddo |
543 |
prof_num_var_tot(num_file,bi,bj)=0 |
544 |
|
545 |
err = NF_INQ_VARID(fid,'prof_T', varid1 ) |
546 |
if (err.EQ.NF_NOERR) then |
547 |
vec_quantities(num_file,1,bi,bj)=.TRUE. |
548 |
prof_num_var_tot(num_file,bi,bj)= |
549 |
& prof_num_var_tot(num_file,bi,bj)+1 |
550 |
prof_num_var_cur(num_file,1,bi,bj)= |
551 |
& prof_num_var_tot(num_file,bi,bj) |
552 |
else |
553 |
vec_quantities(num_file,1,bi,bj)=.FALSE. |
554 |
endif |
555 |
err = NF_INQ_VARID(fid,'prof_S', varid1 ) |
556 |
if (err.EQ.NF_NOERR) then |
557 |
vec_quantities(num_file,2,bi,bj)=.TRUE. |
558 |
prof_num_var_tot(num_file,bi,bj)= |
559 |
& prof_num_var_tot(num_file,bi,bj)+1 |
560 |
prof_num_var_cur(num_file,2,bi,bj)= |
561 |
& prof_num_var_tot(num_file,bi,bj) |
562 |
else |
563 |
vec_quantities(num_file,2,bi,bj)=.FALSE. |
564 |
endif |
565 |
#ifndef ALLOW_PROFILES_GENERICGRID |
566 |
cgf This bloc wont work when model u/v are not |
567 |
cgf zonal/meridional components, while prof_U/V are. |
568 |
err = NF_INQ_VARID(fid,'prof_U', varid1 ) |
569 |
if (err.EQ.NF_NOERR) then |
570 |
vec_quantities(num_file,3,bi,bj)=.TRUE. |
571 |
prof_num_var_tot(num_file,bi,bj)= |
572 |
& prof_num_var_tot(num_file,bi,bj)+1 |
573 |
prof_num_var_cur(num_file,3,bi,bj)= |
574 |
& prof_num_var_tot(num_file,bi,bj) |
575 |
else |
576 |
vec_quantities(num_file,3,bi,bj)=.FALSE. |
577 |
endif |
578 |
err = NF_INQ_VARID(fid,'prof_V', varid1 ) |
579 |
if (err.EQ.NF_NOERR) then |
580 |
vec_quantities(num_file,4,bi,bj)=.TRUE. |
581 |
prof_num_var_tot(num_file,bi,bj)= |
582 |
& prof_num_var_tot(num_file,bi,bj)+1 |
583 |
prof_num_var_cur(num_file,4,bi,bj)= |
584 |
& prof_num_var_tot(num_file,bi,bj) |
585 |
else |
586 |
vec_quantities(num_file,4,bi,bj)=.FALSE. |
587 |
endif |
588 |
#endif |
589 |
err = NF_INQ_VARID(fid,'prof_ptr', varid1 ) |
590 |
if (err.EQ.NF_NOERR) then |
591 |
vec_quantities(num_file,5,bi,bj)=.TRUE. |
592 |
prof_num_var_tot(num_file,bi,bj)= |
593 |
& prof_num_var_tot(num_file,bi,bj)+1 |
594 |
prof_num_var_cur(num_file,5,bi,bj)= |
595 |
& prof_num_var_tot(num_file,bi,bj) |
596 |
else |
597 |
vec_quantities(num_file,5,bi,bj)=.FALSE. |
598 |
endif |
599 |
err = NF_INQ_VARID(fid,'prof_ssh', varid1 ) |
600 |
if (err.EQ.NF_NOERR) then |
601 |
vec_quantities(num_file,6,bi,bj)=.TRUE. |
602 |
prof_num_var_tot(num_file,bi,bj)= |
603 |
& prof_num_var_tot(num_file,bi,bj)+1 |
604 |
prof_num_var_cur(num_file,6,bi,bj)= |
605 |
& prof_num_var_tot(num_file,bi,bj) |
606 |
else |
607 |
vec_quantities(num_file,6,bi,bj)=.FALSE. |
608 |
endif |
609 |
|
610 |
write(msgbuf,'(a,6L5)') ' incl. variables (T/F flags)' |
611 |
&,vec_quantities(num_file,1,bi,bj),vec_quantities(num_file,2,bi,bj) |
612 |
&,vec_quantities(num_file,3,bi,bj),vec_quantities(num_file,4,bi,bj) |
613 |
&,vec_quantities(num_file,5,bi,bj),vec_quantities(num_file,6,bi,bj) |
614 |
call print_message( |
615 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
616 |
|
617 |
C=========================================================== |
618 |
c create files for model counterparts to observations |
619 |
C=========================================================== |
620 |
|
621 |
if (ProfNo(num_file,bi,bj).GT.0) then |
622 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
623 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
624 |
|
625 |
JL = ILNBLNK( profilesDir ) |
626 |
|
627 |
if (profilesDoNcOutput) then |
628 |
|
629 |
write(fnameequinc(1:80),'(3a,i3.3,a,i3.3,a)') |
630 |
& profilesDir(1:JL),profilesfile(1:IL),'.',iG,'.',jG,'.equi.nc' |
631 |
write(adfnameequinc(1:80),'(4a,i3.3,a,i3.3,a)') |
632 |
& profilesDir(1:JL),'ad', |
633 |
& profilesfile(1:IL),'.',iG,'.',jG,'.equi.nc' |
634 |
|
635 |
inquire( file=fnameequinc, exist=exst ) |
636 |
if (.NOT.exst) then |
637 |
call profiles_init_ncfile(num_file, |
638 |
& fiddata(num_file,bi,bj),fnameequinc, |
639 |
& fidforward(num_file,bi,bj),ProfNo(num_file,bi,bj), |
640 |
& ProfDepthNo(num_file,bi,bj), |
641 |
& bi,bj,myThid) |
642 |
else |
643 |
err = NF_OPEN(fnameequinc,NF_WRITE,fidforward(num_file,bi,bj)) |
644 |
endif |
645 |
#ifdef ALLOW_ADJOINT_RUN |
646 |
inquire( file=adfnameequinc, exist=exst ) |
647 |
if (.NOT.exst) then |
648 |
call profiles_init_ncfile(num_file,fiddata(num_file,bi,bj), |
649 |
& adfnameequinc, fidadjoint(num_file,bi,bj), |
650 |
& ProfNo(num_file,bi,bj),ProfDepthNo(num_file,bi,bj), |
651 |
& bi,bj, myThid) |
652 |
else |
653 |
err = NF_OPEN(adfnameequinc,NF_WRITE,fidadjoint(num_file,bi,bj)) |
654 |
endif |
655 |
#endif |
656 |
else |
657 |
|
658 |
write(fnameequinc(1:80),'(3a,i3.3,a,i3.3,a)') |
659 |
& profilesDir(1:JL),profilesfile(1:IL),'.',iG,'.',jG,'.equi.data' |
660 |
write(adfnameequinc(1:80),'(4a,i3.3,a,i3.3,a)') |
661 |
& profilesDir(1:JL),'ad', |
662 |
& profilesfile(1:IL),'.',iG,'.',jG,'.equi.data' |
663 |
|
664 |
inquire( file=fnameequinc, exist=exst ) |
665 |
fidforward(num_file,bi,bj)= |
666 |
& 2000+(num_file-1)*nSx*nSx+(bi-1)*nSy+bj |
667 |
if (.NOT.exst) then |
668 |
call profiles_init_ncfile(num_file,fiddata(num_file,bi,bj), |
669 |
& fnameequinc,fidforward(num_file,bi,bj), |
670 |
& ProfNo(num_file,bi,bj),ProfDepthNo(num_file,bi,bj), |
671 |
& bi,bj,myThid) |
672 |
else |
673 |
c call MDSFINDUNIT( fidforward(num_file,bi,bj) , mythid ) |
674 |
open( fidforward(num_file,bi,bj),file=fnameequinc, |
675 |
& form ='unformatted',status='unknown', access='direct', |
676 |
& recl= (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 ) |
677 |
endif |
678 |
#ifdef ALLOW_ADJOINT_RUN |
679 |
inquire( file=adfnameequinc, exist=exst ) |
680 |
fidadjoint(num_file,bi,bj)= |
681 |
& 4000+(num_file-1)*nSx*nSx+(bi-1)*nSy+bj |
682 |
if (.NOT.exst) then |
683 |
call profiles_init_ncfile(num_file,fiddata(num_file,bi,bj), |
684 |
& adfnameequinc, fidadjoint(num_file,bi,bj), |
685 |
& ProfNo(num_file,bi,bj),ProfDepthNo(num_file,bi,bj), |
686 |
& bi,bj, myThid) |
687 |
else |
688 |
c call MDSFINDUNIT( fidadjoint(num_file,bi,bj) , mythid ) |
689 |
open( fidadjoint(num_file,bi,bj),file=adfnameequinc, |
690 |
& form ='unformatted',status='unknown', access='direct', |
691 |
& recl= (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 ) |
692 |
endif |
693 |
#endif |
694 |
|
695 |
endif |
696 |
|
697 |
endif |
698 |
|
699 |
|
700 |
C=========================================================== |
701 |
else |
702 |
ProfNo(num_file,bi,bj)=0 |
703 |
do k=1,NVARMAX |
704 |
prof_num_var_cur(num_file,k,bi,bj)=0 |
705 |
vec_quantities(num_file,k,bi,bj)=.FALSE. |
706 |
enddo |
707 |
prof_num_var_tot(num_file,bi,bj)=0 |
708 |
do k=1,NOBSGLOB |
709 |
prof_time(num_file,k,bi,bj)=-999 |
710 |
prof_lon(num_file,k,bi,bj)=-999 |
711 |
prof_lat(num_file,k,bi,bj)=-999 |
712 |
prof_ind_glob(num_file,k,bi,bj)=-999 |
713 |
#ifdef ALLOW_PROFILES_GENERICGRID |
714 |
do q = 1,NUM_INTERP_POINTS |
715 |
prof_interp_i(num_file,k,q,bi,bj) = -999 |
716 |
prof_interp_j(num_file,k,q,bi,bj) = -999 |
717 |
prof_interp_weights(num_file,k,q,bi,bj) = -999 |
718 |
enddo |
719 |
prof_interp_xC11(num_file,k,bi,bj)=-999 |
720 |
prof_interp_yC11(num_file,k,bi,bj)=-999 |
721 |
prof_interp_xCNINJ(num_file,k,bi,bj)=-999 |
722 |
prof_interp_yCNINJ(num_file,k,bi,bj)=-999 |
723 |
#endif |
724 |
enddo |
725 |
|
726 |
endif !if (IL.NE.0) then |
727 |
enddo ! do num_file=1,NFILESPROFMAX |
728 |
|
729 |
C=========================================================== |
730 |
C error cases: |
731 |
C=========================================================== |
732 |
|
733 |
#ifdef ALLOW_PROFILES_GENERICGRID |
734 |
|
735 |
c1) you want to provide interpolation information |
736 |
|
737 |
if ( stopGenericGrid.EQ.2.) then |
738 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
739 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
740 |
cgf XC grid |
741 |
call MDSFINDUNIT( fid , mythid ) |
742 |
write(fnameequinc(1:80),'(a,i3.3,a,i3.3,a,i4.4,a,i4.4,a)') |
743 |
& 'profilesXCincl1PointOverlap.',iG,'.',jG,'.',sNx,'.',sNy,'.data' |
744 |
k=MDS_RECLEN(64,(sNx+2)*(sNy+2),mythid) |
745 |
WRITE(standardMessageUnit,'(A,/,2A)') |
746 |
& 'PROFILES_INIT_FIXED: creating grid from profiles; file:', |
747 |
& fnameequinc |
748 |
open( fid, file= fnameequinc, form ='unformatted', |
749 |
& status='unknown',access='direct', recl= k) |
750 |
DO m=0,sNy+1 |
751 |
DO l=0,sNx+1 |
752 |
xy_buffer_r8(l,m)=xC(l,m,bi,bj) |
753 |
ENDDO |
754 |
ENDDO |
755 |
#ifdef _BYTESWAPIO |
756 |
call MDS_BYTESWAPR8((sNx+2)*(sNy+2),xy_buffer_r8) |
757 |
#endif |
758 |
write(fid,rec=1) xy_buffer_r8 |
759 |
close(fid) |
760 |
cgf YC grid |
761 |
call MDSFINDUNIT( fid , mythid ) |
762 |
write(fnameequinc(1:80),'(a,i3.3,a,i3.3,a,i4.4,a,i4.4,a)') |
763 |
& 'profilesYCincl1PointOverlap.',iG,'.',jG,'.',sNx,'.',sNy,'.data' |
764 |
k=MDS_RECLEN(64,(sNx+2)*(sNy+2),mythid) |
765 |
WRITE(standardMessageUnit,'(A,/,A)') |
766 |
& 'PROFILES_INIT_FIXED: creating grid from profiles; file:', |
767 |
& fnameequinc |
768 |
open( fid, file= fnameequinc, form ='unformatted', |
769 |
& status='unknown', access='direct', recl= k) |
770 |
DO m=0,sNy+1 |
771 |
DO l=0,sNx+1 |
772 |
xy_buffer_r8(l,m)=yC(l,m,bi,bj) |
773 |
ENDDO |
774 |
ENDDO |
775 |
#ifdef _BYTESWAPIO |
776 |
call MDS_BYTESWAPR8((sNx+2)*(sNy+2),xy_buffer_r8) |
777 |
#endif |
778 |
write(fid,rec=1) xy_buffer_r8 |
779 |
close(fid) |
780 |
|
781 |
WRITE(msgBuf,'(3A)') |
782 |
& 'PROFILES_INIT_FIXED : ', |
783 |
& 'when using ALLOW_PROFILES_GENERICGRID ', |
784 |
& 'you have to provide interpolation coeffs etc. ' |
785 |
CALL PRINT_ERROR( msgBuf, myThid) |
786 |
WRITE(msgBuf,'(2A)') |
787 |
& 'and some of your nc files dont have them. ', |
788 |
& 'You could use profiles_prep_mygrid.m and/or' |
789 |
CALL PRINT_ERROR( msgBuf, myThid) |
790 |
WRITE(msgBuf,'(A)') |
791 |
& 'use the grid info in profiles*incl1PointOverlap*data' |
792 |
CALL PRINT_ERROR( msgBuf, myThid) |
793 |
stopProfiles=1. _d 0 |
794 |
|
795 |
endif |
796 |
|
797 |
#endif |
798 |
|
799 |
ENDDO |
800 |
ENDDO |
801 |
|
802 |
_END_MASTER( mythid ) |
803 |
_BARRIER |
804 |
|
805 |
c2) stop after other kind of errors |
806 |
_GLOBAL_SUM_RL( stopProfiles , myThid ) |
807 |
if ( stopProfiles.GE.1.) then |
808 |
STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED' |
809 |
endif |
810 |
#ifdef ALLOW_PROFILES_GENERICGRID |
811 |
_GLOBAL_SUM_RL( stopGenericGrid , myThid ) |
812 |
if ( stopGenericGrid.GE.1.) then |
813 |
STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED' |
814 |
endif |
815 |
#endif |
816 |
|
817 |
|
818 |
write(msgbuf,'(a)') ' ' |
819 |
call print_message( msgbuf, standardmessageunit, |
820 |
& SQUEEZE_RIGHT , mythid) |
821 |
write(msgbuf,'(a)') |
822 |
&'// =======================================================' |
823 |
call print_message( msgbuf, standardmessageunit, |
824 |
& SQUEEZE_RIGHT , mythid) |
825 |
write(msgbuf,'(a)') |
826 |
&'// insitu profiles model sampling >>> END <<<' |
827 |
call print_message( msgbuf, standardmessageunit, |
828 |
& SQUEEZE_RIGHT , mythid) |
829 |
write(msgbuf,'(a)') |
830 |
&'// =======================================================' |
831 |
call print_message( msgbuf, standardmessageunit, |
832 |
& SQUEEZE_RIGHT , mythid) |
833 |
write(msgbuf,'(a)') ' ' |
834 |
call print_message( msgbuf, standardmessageunit, |
835 |
& SQUEEZE_RIGHT , mythid) |
836 |
|
837 |
|
838 |
#endif |
839 |
|
840 |
RETURN |
841 |
END |