1 |
C $Header: /u/gcmpack/MITgcm/pkg/profiles/profiles_init_fixed.F,v 1.7 2006/10/25 01:15:54 gforget Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PROFILES_OPTIONS.h" |
5 |
|
6 |
C o==========================================================o |
7 |
C | subroutine profiles_init_fixed | |
8 |
C | o initialization for netcdf profiles data | |
9 |
C | started: Gael Forget 15-March-2006 | |
10 |
C o==========================================================o |
11 |
|
12 |
SUBROUTINE profiles_init_fixed( myThid ) |
13 |
|
14 |
implicit none |
15 |
|
16 |
C ==================== Global Variables =========================== |
17 |
#include "SIZE.h" |
18 |
#include "EEPARAMS.h" |
19 |
#include "PARAMS.h" |
20 |
#include "GRID.h" |
21 |
#include "DYNVARS.h" |
22 |
#ifdef ALLOW_CAL |
23 |
#include "cal.h" |
24 |
#endif |
25 |
#ifdef ALLOW_PROFILES |
26 |
# include "profiles.h" |
27 |
# include "netcdf.inc" |
28 |
#endif |
29 |
C ==================== Routine Variables ========================== |
30 |
|
31 |
integer k,l,m,bi,bj,iG,jG, myThid,num_file,length_for_tile |
32 |
integer fid, dimid, varid1, varid1a, varid1b |
33 |
integer varid2,varid3 |
34 |
_RL tmpyymmdd(1000),tmphhmmss(1000),diffsecs |
35 |
integer tmpdate(4),tmpdiff(4) |
36 |
_RL tmp_lon, tmp_lon2(1000), tmp_lat2(1000) |
37 |
integer vec_start(2), vec_count(2), profno_div1000, kk |
38 |
character*(80) profilesfile, fnamedatanc |
39 |
character*(80) fnameequinc, adfnameequinc |
40 |
integer IL, err |
41 |
logical exst |
42 |
|
43 |
#ifdef ALLOW_PROFILES |
44 |
|
45 |
c == external functions == |
46 |
integer ILNBLNK |
47 |
character*(max_len_mbuf) msgbuf |
48 |
|
49 |
c-- == end of interface == |
50 |
|
51 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
52 |
DO bj = myByLo(myThid), myByHi(myThid) |
53 |
|
54 |
profiles_curfile_buff(bi,bj)=0 |
55 |
|
56 |
do m=1,NLEVELMAX |
57 |
do l=1,1000 |
58 |
do k=1,NVARMAX |
59 |
profiles_data_buff(m,l,k,bi,bj)=0 |
60 |
profiles_weight_buff(m,l,k,bi,bj)=0 |
61 |
enddo |
62 |
enddo |
63 |
enddo |
64 |
|
65 |
c remplacer par une boucle ensuite : |
66 |
do num_file=1,NFILESPROFMAX |
67 |
|
68 |
IL = ILNBLNK( profilesfiles(num_file) ) |
69 |
if (IL.NE.0) then |
70 |
write(profilesfile(1:80),'(1a)') |
71 |
& profilesfiles(num_file)(1:IL) |
72 |
write(msgbuf,'(a,X,i3,X,a)') |
73 |
& 'Profiles num_file is ', num_file, profilesfile(1:80) |
74 |
call print_message( |
75 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
76 |
else |
77 |
write(profilesfile(1:80),'(1a)') ' ' |
78 |
write(msgbuf,'(a,X,i3,X,a)') |
79 |
& 'Profiles num_file is ', num_file, ' empty ' |
80 |
call print_message( |
81 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
82 |
endif |
83 |
|
84 |
IL = ILNBLNK( profilesfile ) |
85 |
if (IL.NE.0) then |
86 |
|
87 |
C=========================================================== |
88 |
c open data files and read the position vectors |
89 |
C=========================================================== |
90 |
|
91 |
write(fnamedatanc(1:80),'(2a)') profilesfile(1:IL),'.nc' |
92 |
write(msgbuf,'(a,X,i3,X,a)') |
93 |
& 'Opening num_file ', num_file, fnamedatanc(1:80) |
94 |
call print_message( |
95 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
96 |
err = NF_OPEN(fnamedatanc, 0, fiddata(num_file,bi,bj)) |
97 |
|
98 |
c1) read the number of profiles : |
99 |
cgf err = NF_OPEN(filename, 0, fid) |
100 |
fid=fiddata(num_file,bi,bj) |
101 |
err = NF_INQ_DIMID(fid,'iPROF', dimid ) |
102 |
err = NF_INQ_DIMLEN(fid, dimid, ProfNo(num_file,bi,bj) ) |
103 |
err = NF_INQ_DIMID(fid,'iDEPTH', dimid ) |
104 |
if (err.NE.NF_NOERR) then |
105 |
err = NF_INQ_DIMID(fid,'Z', dimid ) |
106 |
else |
107 |
STOP 'ERROR reading profile file in profiles_init_fixed' |
108 |
endif |
109 |
err = NF_INQ_DIMLEN(fid, dimid, ProfDepthNo(num_file,bi,bj) ) |
110 |
write(msgbuf,'(a,X,4i9)') |
111 |
& ' fid, num_file, ProfNo, ProfDepthNo ', |
112 |
& fid, num_file, ProfNo(num_file,bi,bj), |
113 |
& ProfDepthNo(num_file,bi,bj) |
114 |
call print_message( |
115 |
& msgbuf, standardmessageunit, SQUEEZE_RIGHT , mythid) |
116 |
|
117 |
c2) read the dates and positions : |
118 |
err = NF_INQ_VARID(fid,'depth', varid1a ) |
119 |
do k=1,ProfDepthNo(num_file,bi,bj) |
120 |
err = NF_GET_VAR1_DOUBLE(fid,varid1a,k, |
121 |
& prof_depth(num_file,k,bi,bj)) |
122 |
enddo |
123 |
|
124 |
err = NF_INQ_VARID(fid,'prof_YYYYMMDD', varid1a ) |
125 |
err = NF_INQ_VARID(fid,'prof_HHMMSS', varid1b ) |
126 |
err = NF_INQ_VARID(fid,'prof_lon', varid2 ) |
127 |
err = NF_INQ_VARID(fid,'prof_lat', varid3 ) |
128 |
|
129 |
c DO bi = myBxLo(myThid), myBxHi(myThid) |
130 |
c DO bj = myByLo(myThid), myByHi(myThid) |
131 |
|
132 |
do k=1,NOBSGLOB |
133 |
prof_time(num_file,k,bi,bj)=-999 |
134 |
prof_lon(num_file,k,bi,bj)=-999 |
135 |
prof_lat(num_file,k,bi,bj)=-999 |
136 |
prof_ind_glob(num_file,k,bi,bj)=-999 |
137 |
enddo |
138 |
|
139 |
|
140 |
length_for_tile=0 |
141 |
profno_div1000=max(0,int(ProfNo(num_file,bi,bj)/1000)) |
142 |
|
143 |
do kk=1,profno_div1000+1 |
144 |
|
145 |
if (min(ProfNo(num_file,bi,bj), 1000*kk).GE. |
146 |
& 1+1000*(kk-1)) then |
147 |
|
148 |
vec_start(1)=1 |
149 |
vec_start(2)=1+1000*(kk-1) |
150 |
vec_count(1)=1 |
151 |
vec_count(2)=min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
152 |
|
153 |
if ( (vec_count(2).LE.0).OR.(vec_count(2).GT.1000).OR. |
154 |
& (vec_start(2).LE.0).OR. |
155 |
& (vec_count(2)+vec_start(2)-1.GT.ProfNo(num_file,bi,bj)) ) |
156 |
& then |
157 |
print*,"stop 1",vec_start, vec_count |
158 |
stop |
159 |
endif |
160 |
|
161 |
err = NF_GET_VARA_DOUBLE(fid,varid1a,vec_start(2), |
162 |
& vec_count(2), tmpyymmdd) |
163 |
err = NF_GET_VARA_DOUBLE(fid,varid1b,vec_start(2), |
164 |
& vec_count(2), tmphhmmss) |
165 |
err = NF_GET_VARA_DOUBLE(fid,varid2,vec_start(2), |
166 |
& vec_count(2), tmp_lon2) |
167 |
err = NF_GET_VARA_DOUBLE(fid,varid3,vec_start(2), |
168 |
& vec_count(2), tmp_lat2) |
169 |
|
170 |
if (err.NE.NF_NOERR) then |
171 |
print*,"stop 2",vec_start(2),vec_count(2), |
172 |
& kk,min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
173 |
stop |
174 |
endif |
175 |
|
176 |
do k=1,min(1000,ProfNo(num_file,bi,bj)-1000*(kk-1)) |
177 |
|
178 |
call cal_FullDate( int(tmpyymmdd(k)),int(tmphhmmss(k)), |
179 |
& tmpdate,bi,bj,mythid ) |
180 |
call cal_TimePassed( modelstartdate,tmpdate,tmpdiff,mythid ) |
181 |
call cal_ToSeconds (tmpdiff,diffsecs,mythid) |
182 |
diffsecs=diffsecs+nIter0*deltaTclock |
183 |
|
184 |
if (xC(sNx+1,1,bi,bj).LT.xC(1,1,bi,bj)) then |
185 |
tmp_lon=xC(sNx+1,1,bi,bj)+360 |
186 |
else |
187 |
tmp_lon=xC(sNx+1,1,bi,bj) |
188 |
endif |
189 |
if ((xC(1,1,bi,bj).LE.tmp_lon2(k)).AND. |
190 |
& (tmp_lon.GT.tmp_lon2(k)).AND. |
191 |
& (yC(1,1,bi,bj).LE.tmp_lat2(k)).AND. |
192 |
& (yC(1,sNy+1,bi,bj).GT.tmp_lat2(k)) |
193 |
& ) then |
194 |
length_for_tile=length_for_tile+1 |
195 |
prof_time(num_file,length_for_tile,bi,bj)=diffsecs |
196 |
prof_lon(num_file,length_for_tile,bi,bj)=tmp_lon2(k) |
197 |
prof_lat(num_file,length_for_tile,bi,bj)=tmp_lat2(k) |
198 |
prof_ind_glob(num_file,length_for_tile,bi,bj)=k+1000*(kk-1) |
199 |
if (length_for_tile.GT.NOBSGLOB) then |
200 |
print*,"too much profiles: need to increase NOBSGLOB," |
201 |
print*," or split the data file (less memory cost)" |
202 |
stop |
203 |
endif |
204 |
elseif (xC(sNx+1,1,bi,bj).LT.xC(1,1,bi,bj)) then |
205 |
if ((xC(1,1,bi,bj).LE.tmp_lon2(k)+360).AND. |
206 |
& (tmp_lon.GT.tmp_lon2(k)+360).AND. |
207 |
& (yC(1,1,bi,bj).LE.tmp_lat2(k)).AND. |
208 |
& (yC(1,sNy+1,bi,bj).GT.tmp_lat2(k)) |
209 |
& ) then |
210 |
length_for_tile=length_for_tile+1 |
211 |
prof_time(num_file,length_for_tile,bi,bj)=diffsecs |
212 |
prof_lon(num_file,length_for_tile,bi,bj)=tmp_lon2(k)+360 |
213 |
prof_lat(num_file,length_for_tile,bi,bj)=tmp_lat2(k) |
214 |
prof_ind_glob(num_file,length_for_tile,bi,bj)=k+1000*(kk-1) |
215 |
if (length_for_tile.GT.NOBSGLOB) then |
216 |
print*,"too much profiles: need to increase NOBSGLOB," |
217 |
print*," or split the data file (less memory cost)" |
218 |
stop |
219 |
endif |
220 |
endif |
221 |
endif |
222 |
enddo |
223 |
endif |
224 |
enddo |
225 |
|
226 |
ProfNo(num_file,bi,bj)=length_for_tile |
227 |
print*,"fid dimid ProfNo(num_file,bi,bj)",fid, dimid, |
228 |
& num_file, ProfNo(num_file,bi,bj) |
229 |
|
230 |
do k=1,NVARMAX |
231 |
prof_num_var_cur(num_file,k,bi,bj)=0 |
232 |
enddo |
233 |
prof_num_var_tot(num_file,bi,bj)=0 |
234 |
|
235 |
c3) detect available data types |
236 |
err = NF_INQ_VARID(fid,'prof_T', varid1 ) |
237 |
if (err.EQ.NF_NOERR) then |
238 |
vec_quantities(num_file,1,bi,bj)=.TRUE. |
239 |
prof_num_var_tot(num_file,bi,bj)= |
240 |
& prof_num_var_tot(num_file,bi,bj)+1 |
241 |
prof_num_var_cur(num_file,1,bi,bj)= |
242 |
& prof_num_var_tot(num_file,bi,bj) |
243 |
else |
244 |
vec_quantities(num_file,1,bi,bj)=.FALSE. |
245 |
endif |
246 |
err = NF_INQ_VARID(fid,'prof_S', varid1 ) |
247 |
if (err.EQ.NF_NOERR) then |
248 |
vec_quantities(num_file,2,bi,bj)=.TRUE. |
249 |
prof_num_var_tot(num_file,bi,bj)= |
250 |
& prof_num_var_tot(num_file,bi,bj)+1 |
251 |
prof_num_var_cur(num_file,2,bi,bj)= |
252 |
& prof_num_var_tot(num_file,bi,bj) |
253 |
else |
254 |
vec_quantities(num_file,2,bi,bj)=.FALSE. |
255 |
endif |
256 |
err = NF_INQ_VARID(fid,'prof_U', varid1 ) |
257 |
if (err.EQ.NF_NOERR) then |
258 |
vec_quantities(num_file,3,bi,bj)=.TRUE. |
259 |
prof_num_var_tot(num_file,bi,bj)= |
260 |
& prof_num_var_tot(num_file,bi,bj)+1 |
261 |
prof_num_var_cur(num_file,3,bi,bj)= |
262 |
& prof_num_var_tot(num_file,bi,bj) |
263 |
else |
264 |
vec_quantities(num_file,3,bi,bj)=.FALSE. |
265 |
endif |
266 |
err = NF_INQ_VARID(fid,'prof_V', varid1 ) |
267 |
if (err.EQ.NF_NOERR) then |
268 |
vec_quantities(num_file,4,bi,bj)=.TRUE. |
269 |
prof_num_var_tot(num_file,bi,bj)= |
270 |
& prof_num_var_tot(num_file,bi,bj)+1 |
271 |
prof_num_var_cur(num_file,4,bi,bj)= |
272 |
& prof_num_var_tot(num_file,bi,bj) |
273 |
else |
274 |
vec_quantities(num_file,4,bi,bj)=.FALSE. |
275 |
endif |
276 |
err = NF_INQ_VARID(fid,'prof_ptr', varid1 ) |
277 |
if (err.EQ.NF_NOERR) then |
278 |
vec_quantities(num_file,5,bi,bj)=.TRUE. |
279 |
prof_num_var_tot(num_file,bi,bj)= |
280 |
& prof_num_var_tot(num_file,bi,bj)+1 |
281 |
prof_num_var_cur(num_file,5,bi,bj)= |
282 |
& prof_num_var_tot(num_file,bi,bj) |
283 |
else |
284 |
vec_quantities(num_file,5,bi,bj)=.FALSE. |
285 |
endif |
286 |
err = NF_INQ_VARID(fid,'prof_ssh', varid1 ) |
287 |
if (err.EQ.NF_NOERR) then |
288 |
vec_quantities(num_file,6,bi,bj)=.TRUE. |
289 |
prof_num_var_tot(num_file,bi,bj)= |
290 |
& prof_num_var_tot(num_file,bi,bj)+1 |
291 |
prof_num_var_cur(num_file,6,bi,bj)= |
292 |
& prof_num_var_tot(num_file,bi,bj) |
293 |
else |
294 |
vec_quantities(num_file,6,bi,bj)=.FALSE. |
295 |
endif |
296 |
|
297 |
|
298 |
C=========================================================== |
299 |
c create files for model counterparts to observations |
300 |
C=========================================================== |
301 |
|
302 |
if (ProfNo(num_file,bi,bj).GT.0) then |
303 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
304 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
305 |
|
306 |
if (profilesfile_equi_type.EQ.1) then |
307 |
|
308 |
write(fnameequinc(1:80),'(2a,i3.3,a,i3.3,a)') |
309 |
& profilesfile(1:IL),'.',iG,'.',jG,'.equi.nc' |
310 |
write(adfnameequinc(1:80),'(3a,i3.3,a,i3.3,a)') 'ad', |
311 |
& profilesfile(1:IL),'.',iG,'.',jG,'.equi.nc' |
312 |
|
313 |
inquire( file=fnameequinc, exist=exst ) |
314 |
if (.NOT.exst) then |
315 |
call profiles_init_ncfile(num_file, |
316 |
& fiddata(num_file,bi,bj),fnameequinc, |
317 |
& fidforward(num_file,bi,bj),ProfNo(num_file,bi,bj), |
318 |
& ProfDepthNo(num_file,bi,bj), |
319 |
& bi,bj,myThid) |
320 |
call profiles_init_ncfile(num_file,fiddata(num_file,bi,bj), |
321 |
& adfnameequinc, fidadjoint(num_file,bi,bj),ProfNo(num_file,bi,bj), |
322 |
& ProfDepthNo(num_file,bi,bj),bi,bj, myThid) |
323 |
else |
324 |
err = NF_OPEN(fnameequinc,NF_WRITE,fidforward(num_file,bi,bj)) |
325 |
err = NF_OPEN(adfnameequinc,NF_WRITE,fidadjoint(num_file,bi,bj)) |
326 |
endif |
327 |
|
328 |
else |
329 |
|
330 |
write(fnameequinc(1:80),'(2a,i3.3,a,i3.3,a)') |
331 |
& profilesfile(1:IL),'.',iG,'.',jG,'.equi.bin' |
332 |
write(adfnameequinc(1:80),'(3a,i3.3,a,i3.3,a)') 'ad', |
333 |
& profilesfile(1:IL),'.',iG,'.',jG,'.equi.bin' |
334 |
|
335 |
inquire( file=fnameequinc, exist=exst ) |
336 |
if (.NOT.exst) then |
337 |
call profiles_init_ncfile(num_file,fiddata(num_file,bi,bj), |
338 |
& fnameequinc,fidforward(num_file,bi,bj), |
339 |
& ProfNo(num_file,bi,bj),ProfDepthNo(num_file,bi,bj), |
340 |
& bi,bj,myThid) |
341 |
call profiles_init_ncfile(num_file,fiddata(num_file,bi,bj), |
342 |
& adfnameequinc, fidadjoint(num_file,bi,bj),ProfNo(num_file,bi,bj), |
343 |
& ProfDepthNo(num_file,bi,bj),bi,bj, myThid) |
344 |
else |
345 |
call MDSFINDUNIT( fidforward(num_file,bi,bj) , mythid ) |
346 |
open( fidforward(num_file,bi,bj),file=fnameequinc, |
347 |
& form ='unformatted',status='unknown', access='direct', |
348 |
& recl= (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 ) |
349 |
call MDSFINDUNIT( fidadjoint(num_file,bi,bj) , mythid ) |
350 |
open( fidadjoint(num_file,bi,bj),file=adfnameequinc, |
351 |
& form ='unformatted',status='unknown', access='direct', |
352 |
& recl= (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 ) |
353 |
endif |
354 |
|
355 |
endif |
356 |
|
357 |
endif |
358 |
|
359 |
c ENDDO |
360 |
c ENDDO |
361 |
|
362 |
|
363 |
C=========================================================== |
364 |
else |
365 |
ProfNo(num_file,bi,bj)=0 |
366 |
do k=1,NVARMAX |
367 |
prof_num_var_cur(num_file,k,bi,bj)=0 |
368 |
vec_quantities(num_file,k,bi,bj)=.FALSE. |
369 |
enddo |
370 |
prof_num_var_tot(num_file,bi,bj)=0 |
371 |
do k=1,NOBSGLOB |
372 |
prof_time(num_file,k,bi,bj)=-999 |
373 |
prof_lon(num_file,k,bi,bj)=-999 |
374 |
prof_lat(num_file,k,bi,bj)=-999 |
375 |
prof_ind_glob(num_file,k,bi,bj)=-999 |
376 |
enddo |
377 |
|
378 |
endif !if (IL.NE.0) then |
379 |
enddo ! do num_file=1,NFILESPROFMAX |
380 |
C=========================================================== |
381 |
|
382 |
ENDDO |
383 |
ENDDO |
384 |
|
385 |
#endif |
386 |
|
387 |
END |
388 |
|