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