global_cycle 1.14.0
Loading...
Searching...
No Matches
read_write_data.f90
Go to the documentation of this file.
1
5
9
10 USE netcdf
11
12 PRIVATE
13
14! DATA STRUCTURE TO HOLD ALL NSST RECORDS.
15
16 TYPE, PUBLIC :: nsst_data
17 REAL, ALLOCATABLE :: C_0(:)
18 REAL, ALLOCATABLE :: C_D(:)
19 REAL, ALLOCATABLE :: D_CONV(:)
20 REAL, ALLOCATABLE :: DT_COOL(:)
21 REAL, ALLOCATABLE :: IFD(:)
22 REAL, ALLOCATABLE :: QRAIN(:)
23 REAL, ALLOCATABLE :: TREF(:)
24 REAL, ALLOCATABLE :: TFINC(:)
25 REAL, ALLOCATABLE :: W_0(:)
26 REAL, ALLOCATABLE :: W_D(:)
27 REAL, ALLOCATABLE :: XS(:)
28 REAL, ALLOCATABLE :: XT(:)
29 REAL, ALLOCATABLE :: XTTS(:)
30 REAL, ALLOCATABLE :: XU(:)
31 REAL, ALLOCATABLE :: XV(:)
32 REAL, ALLOCATABLE :: XZ(:)
33 REAL, ALLOCATABLE :: XZTS(:)
34 REAL, ALLOCATABLE :: Z_C(:)
35 REAL, ALLOCATABLE :: ZM(:)
36 END TYPE nsst_data
37
38 INTEGER, PUBLIC :: idim_gaus
40 INTEGER, PUBLIC :: jdim_gaus
42 INTEGER, ALLOCATABLE, PUBLIC :: slmsk_gaus(:,:)
44
45 INTEGER, ALLOCATABLE, PUBLIC :: soilsnow_gaus(:,:)
48
49 REAL, ALLOCATABLE, PUBLIC :: dtref_gaus(:,:)
51
52 REAL, ALLOCATABLE, PUBLIC :: stc_inc_gaus(:,:,:)
54
55 REAL, ALLOCATABLE, PUBLIC :: slc_inc_gaus(:,:,:)
57
58 PUBLIC :: read_data
59 PUBLIC :: read_gsi_data
60 PUBLIC :: read_lat_lon_orog
61 PUBLIC :: write_data
64
65 CONTAINS
66
122
123 subroutine write_data(lensfc,idim,jdim,lsoil, &
124 do_nsst,inc_file,nsst,slifcs,tsffcs,vegfcs,swefcs, &
125 tg3fcs,zorfcs,albfcs,alffcs, &
126 cnpfcs,f10m,t2m,q2m,vetfcs, &
127 sotfcs,ustar,fmm,fhh,sicfcs, &
128 sihfcs,sitfcs,tprcp,srflag, &
129 swdfcs,vmnfcs,vmxfcs,slpfcs, &
130 absfcs,slcfcs,smcfcs,stcfcs, &
131 stcinc, slcinc)
132
133 use mpi
134
135 implicit none
136
137 integer, intent(in) :: lensfc, lsoil
138 integer, intent(in) :: idim, jdim
139
140 logical, intent(in) :: do_nsst
141 logical, intent(in) :: inc_file
142
143 real, intent(in), optional :: slifcs(lensfc),tsffcs(lensfc)
144 real, intent(in), optional :: swefcs(lensfc),tg3fcs(lensfc)
145 real, intent(in), optional :: zorfcs(lensfc),albfcs(lensfc,4)
146 real, intent(in), optional :: alffcs(lensfc,2),cnpfcs(lensfc)
147 real, intent(in), optional :: f10m(lensfc),t2m(lensfc)
148 real, intent(in), optional :: q2m(lensfc),vegfcs(lensfc)
149 real, intent(in), optional :: vetfcs(lensfc),sotfcs(lensfc)
150 real, intent(in), optional :: ustar(lensfc),fmm(lensfc)
151 real, intent(in), optional :: fhh(lensfc), sicfcs(lensfc)
152 real, intent(in), optional :: sihfcs(lensfc), sitfcs(lensfc)
153 real, intent(in), optional :: tprcp(lensfc), srflag(lensfc)
154 real, intent(in), optional :: swdfcs(lensfc), vmnfcs(lensfc)
155 real, intent(in), optional :: vmxfcs(lensfc), slpfcs(lensfc)
156 real, intent(in), optional :: absfcs(lensfc), slcfcs(lensfc,lsoil)
157 real, intent(in), optional :: smcfcs(lensfc,lsoil), stcfcs(lensfc,lsoil)
158 real, intent(in), optional :: stcinc(lensfc,lsoil)
159 real, intent(in), optional :: slcinc(lensfc,lsoil)
160
161 type(nsst_data), intent(in) :: nsst
162
163 integer :: dim_x, dim_y, dim_soil, dim_time, dims_3d(3)
164
165 real :: dum2d(idim,jdim), dum3d(idim,jdim,lsoil)
166
167 character(len=50) :: fnbgso
168 character(len=3) :: rankch
169
170 integer :: myrank, error, ncid, id_var
171 integer :: varid_stc, varid_slc
172
173 call mpi_comm_rank(mpi_comm_world, myrank, error)
174
175 write(rankch, '(i3.3)') (myrank+1)
176
177 if (.NOT.(inc_file)) then
178
179 fnbgso = "./fnbgso." // rankch
180
181 print*
182 print*,"update OUTPUT SFC DATA TO: ",trim(fnbgso)
183
184 error=nf90_open(trim(fnbgso),nf90_write,ncid)
185 CALL netcdf_err(error, 'OPENING FILE: '//trim(fnbgso) )
186
187 if(present(slifcs)) then
188 error=nf90_inq_varid(ncid, "slmsk", id_var)
189 call netcdf_err(error, 'reading slmsk id' )
190 dum2d = reshape(slifcs, (/idim,jdim/))
191 error = nf90_put_var( ncid, id_var, dum2d)
192 call netcdf_err(error, 'writing slmsk record' )
193 call remove_checksum(ncid, id_var)
194 endif
195
196 if(present(tsffcs)) then
197 error=nf90_inq_varid(ncid, "tsea", id_var)
198 call netcdf_err(error, 'reading tsea id' )
199 dum2d = reshape(tsffcs, (/idim,jdim/))
200 error = nf90_put_var( ncid, id_var, dum2d)
201 call netcdf_err(error, 'writing tsea record' )
202 call remove_checksum(ncid, id_var)
203 endif
204
205 if(present(swefcs)) then
206 error=nf90_inq_varid(ncid, "sheleg", id_var)
207 call netcdf_err(error, 'reading sheleg id' )
208 dum2d = reshape(swefcs, (/idim,jdim/))
209 error = nf90_put_var( ncid, id_var, dum2d)
210 call netcdf_err(error, 'writing sheleg record' )
211 call remove_checksum(ncid, id_var)
212 endif
213
214 if(present(tg3fcs)) then
215 error=nf90_inq_varid(ncid, "tg3", id_var)
216 call netcdf_err(error, 'reading tg3 id' )
217 dum2d = reshape(tg3fcs, (/idim,jdim/))
218 error = nf90_put_var( ncid, id_var, dum2d)
219 call netcdf_err(error, 'writing tg3 record' )
220 call remove_checksum(ncid, id_var)
221 endif
222
223 if(present(zorfcs)) then
224 error=nf90_inq_varid(ncid, "zorl", id_var)
225 call netcdf_err(error, 'reading zorl id' )
226 dum2d = reshape(zorfcs, (/idim,jdim/))
227 error = nf90_put_var( ncid, id_var, dum2d)
228 call netcdf_err(error, 'writing zorl record' )
229 call remove_checksum(ncid, id_var)
230 endif
231
232 if(present(albfcs)) then
233 error=nf90_inq_varid(ncid, "alvsf", id_var)
234 call netcdf_err(error, 'reading alvsf id' )
235 dum2d = reshape(albfcs(:,1), (/idim,jdim/))
236 error = nf90_put_var( ncid, id_var, dum2d)
237 call netcdf_err(error, 'writing alvsf record' )
238 call remove_checksum(ncid, id_var)
239
240 error=nf90_inq_varid(ncid, "alvwf", id_var)
241 call netcdf_err(error, 'reading alvwf id' )
242 dum2d = reshape(albfcs(:,2), (/idim,jdim/))
243 error = nf90_put_var( ncid, id_var, dum2d)
244 call netcdf_err(error, 'writing alvwf record' )
245 call remove_checksum(ncid, id_var)
246
247 error=nf90_inq_varid(ncid, "alnsf", id_var)
248 call netcdf_err(error, 'reading alnsf id' )
249 dum2d = reshape(albfcs(:,3), (/idim,jdim/))
250 error = nf90_put_var( ncid, id_var, dum2d)
251 call netcdf_err(error, 'writing alnsf record' )
252 call remove_checksum(ncid, id_var)
253
254 error=nf90_inq_varid(ncid, "alnwf", id_var)
255 call netcdf_err(error, 'reading alnwf id' )
256 dum2d = reshape(albfcs(:,4), (/idim,jdim/))
257 error = nf90_put_var( ncid, id_var, dum2d)
258 call netcdf_err(error, 'writing alnwf record' )
259 call remove_checksum(ncid, id_var)
260 endif
261
262 if(present(alffcs)) then
263 error=nf90_inq_varid(ncid, "facsf", id_var)
264 call netcdf_err(error, 'reading facsf id' )
265 dum2d = reshape(alffcs(:,1), (/idim,jdim/))
266 error = nf90_put_var( ncid, id_var, dum2d)
267 call netcdf_err(error, 'writing facsf record' )
268 call remove_checksum(ncid, id_var)
269
270 error=nf90_inq_varid(ncid, "facwf", id_var)
271 call netcdf_err(error, 'reading facwf id' )
272 dum2d = reshape(alffcs(:,2), (/idim,jdim/))
273 error = nf90_put_var( ncid, id_var, dum2d)
274 call netcdf_err(error, 'writing facwf record' )
275 call remove_checksum(ncid, id_var)
276 endif
277
278 if(present(vegfcs)) then
279 error=nf90_inq_varid(ncid, "vfrac", id_var)
280 call netcdf_err(error, 'reading vfrac id' )
281 dum2d = reshape(vegfcs, (/idim,jdim/))
282 error = nf90_put_var( ncid, id_var, dum2d)
283 call netcdf_err(error, 'writing vegfcs record' )
284 call remove_checksum(ncid, id_var)
285 endif
286
287 if(present(cnpfcs)) then
288 error=nf90_inq_varid(ncid, "canopy", id_var)
289 call netcdf_err(error, 'reading canopy id' )
290 dum2d = reshape(cnpfcs, (/idim,jdim/))
291 error = nf90_put_var( ncid, id_var, dum2d)
292 call netcdf_err(error, 'writing canopy record' )
293 call remove_checksum(ncid, id_var)
294 endif
295
296 if(present(f10m)) then
297 error=nf90_inq_varid(ncid, "f10m", id_var)
298 call netcdf_err(error, 'reading f10m id' )
299 dum2d = reshape(f10m, (/idim,jdim/))
300 error = nf90_put_var( ncid, id_var, dum2d)
301 call netcdf_err(error, 'writing f10m record' )
302 call remove_checksum(ncid, id_var)
303 endif
304
305 if(present(t2m)) then
306 error=nf90_inq_varid(ncid, "t2m", id_var)
307 call netcdf_err(error, 'reading t2m id' )
308 dum2d = reshape(t2m, (/idim,jdim/))
309 error = nf90_put_var( ncid, id_var, dum2d)
310 call netcdf_err(error, 'writing t2m record' )
311 call remove_checksum(ncid, id_var)
312 endif
313
314 if(present(q2m)) then
315 error=nf90_inq_varid(ncid, "q2m", id_var)
316 call netcdf_err(error, 'reading q2m id' )
317 dum2d = reshape(q2m, (/idim,jdim/))
318 error = nf90_put_var( ncid, id_var, dum2d)
319 call netcdf_err(error, 'writing q2m record' )
320 call remove_checksum(ncid, id_var)
321 endif
322
323 if(present(vetfcs)) then
324 error=nf90_inq_varid(ncid, "vtype", id_var)
325 call netcdf_err(error, 'reading vtype id' )
326 dum2d = reshape(vetfcs, (/idim,jdim/))
327 error = nf90_put_var( ncid, id_var, dum2d)
328 call netcdf_err(error, 'writing vtype record' )
329 call remove_checksum(ncid, id_var)
330 endif
331
332 if(present(sotfcs)) then
333 error=nf90_inq_varid(ncid, "stype", id_var)
334 call netcdf_err(error, 'reading stype id' )
335 dum2d = reshape(sotfcs, (/idim,jdim/))
336 error = nf90_put_var( ncid, id_var, dum2d)
337 call netcdf_err(error, 'writing stype record' )
338 call remove_checksum(ncid, id_var)
339 endif
340
341 if(present(ustar)) then
342 error=nf90_inq_varid(ncid, "uustar", id_var)
343 call netcdf_err(error, 'reading uustar id' )
344 dum2d = reshape(ustar, (/idim,jdim/))
345 error = nf90_put_var( ncid, id_var, dum2d)
346 call netcdf_err(error, 'writing uustar record' )
347 call remove_checksum(ncid, id_var)
348 endif
349
350 if(present(fmm)) then
351 error=nf90_inq_varid(ncid, "ffmm", id_var)
352 call netcdf_err(error, 'reading ffmm id' )
353 dum2d = reshape(fmm, (/idim,jdim/))
354 error = nf90_put_var( ncid, id_var, dum2d)
355 call netcdf_err(error, 'writing ffmm record' )
356 call remove_checksum(ncid, id_var)
357 endif
358
359 if(present(fhh)) then
360 error=nf90_inq_varid(ncid, "ffhh", id_var)
361 call netcdf_err(error, 'reading ffhh id' )
362 dum2d = reshape(fhh, (/idim,jdim/))
363 error = nf90_put_var( ncid, id_var, dum2d)
364 call netcdf_err(error, 'writing ffhh record' )
365 call remove_checksum(ncid, id_var)
366 endif
367
368 if(present(sicfcs)) then
369 error=nf90_inq_varid(ncid, "fice", id_var)
370 call netcdf_err(error, 'reading fice id' )
371 dum2d = reshape(sicfcs, (/idim,jdim/))
372 error = nf90_put_var( ncid, id_var, dum2d)
373 call netcdf_err(error, 'writing fice record' )
374 call remove_checksum(ncid, id_var)
375 endif
376
377 if(present(sihfcs)) then
378 error=nf90_inq_varid(ncid, "hice", id_var)
379 call netcdf_err(error, 'reading hice id' )
380 dum2d = reshape(sihfcs, (/idim,jdim/))
381 error = nf90_put_var( ncid, id_var, dum2d)
382 call netcdf_err(error, 'writing hice record' )
383 call remove_checksum(ncid, id_var)
384 endif
385
386 if(present(sitfcs)) then
387 error=nf90_inq_varid(ncid, "tisfc", id_var)
388 call netcdf_err(error, 'reading tisfc id' )
389 dum2d = reshape(sitfcs, (/idim,jdim/))
390 error = nf90_put_var( ncid, id_var, dum2d)
391 call netcdf_err(error, 'writing tisfc record' )
392 call remove_checksum(ncid, id_var)
393 endif
394
395 if(present(tprcp)) then
396 error=nf90_inq_varid(ncid, "tprcp", id_var)
397 call netcdf_err(error, 'reading tprcp id' )
398 dum2d = reshape(tprcp, (/idim,jdim/))
399 error = nf90_put_var( ncid, id_var, dum2d)
400 call netcdf_err(error, 'writing tprcp record' )
401 call remove_checksum(ncid, id_var)
402 endif
403
404 if(present(srflag)) then
405 error=nf90_inq_varid(ncid, "srflag", id_var)
406 call netcdf_err(error, 'reading srflag id' )
407 dum2d = reshape(srflag, (/idim,jdim/))
408 error = nf90_put_var( ncid, id_var, dum2d)
409 call netcdf_err(error, 'writing srflag record' )
410 call remove_checksum(ncid, id_var)
411 endif
412
413 if(present(swdfcs)) then
414 error=nf90_inq_varid(ncid, "snwdph", id_var)
415 call netcdf_err(error, 'reading snwdph id' )
416 dum2d = reshape(swdfcs, (/idim,jdim/))
417 error = nf90_put_var( ncid, id_var, dum2d)
418 call netcdf_err(error, 'writing snwdph record' )
419 call remove_checksum(ncid, id_var)
420 endif
421
422 if(present(vmnfcs)) then
423 error=nf90_inq_varid(ncid, "shdmin", id_var)
424 call netcdf_err(error, 'reading shdmin id' )
425 dum2d = reshape(vmnfcs, (/idim,jdim/))
426 error = nf90_put_var( ncid, id_var, dum2d)
427 call netcdf_err(error, 'writing shdmin record' )
428 call remove_checksum(ncid, id_var)
429 endif
430
431 if(present(vmxfcs)) then
432 error=nf90_inq_varid(ncid, "shdmax", id_var)
433 call netcdf_err(error, 'reading shdmax id' )
434 dum2d = reshape(vmxfcs, (/idim,jdim/))
435 error = nf90_put_var( ncid, id_var, dum2d)
436 call netcdf_err(error, 'writing shdmax record' )
437 call remove_checksum(ncid, id_var)
438 endif
439
440 if(present(slpfcs)) then
441 error=nf90_inq_varid(ncid, "slope", id_var)
442 call netcdf_err(error, 'reading slope id' )
443 dum2d = reshape(slpfcs, (/idim,jdim/))
444 error = nf90_put_var( ncid, id_var, dum2d)
445 call netcdf_err(error, 'writing slope record' )
446 call remove_checksum(ncid, id_var)
447 endif
448
449 if(present(absfcs)) then
450 error=nf90_inq_varid(ncid, "snoalb", id_var)
451 call netcdf_err(error, 'reading snoalb id' )
452 dum2d = reshape(absfcs, (/idim,jdim/))
453 error = nf90_put_var( ncid, id_var, dum2d)
454 call netcdf_err(error, 'writing snoalb record' )
455 call remove_checksum(ncid, id_var)
456 endif
457
458 if(present(slcfcs)) then
459 error=nf90_inq_varid(ncid, "slc", id_var)
460 call netcdf_err(error, 'reading slc id' )
461 dum3d = reshape(slcfcs, (/idim,jdim,lsoil/))
462 error = nf90_put_var( ncid, id_var, dum3d)
463 call netcdf_err(error, 'writing slc record' )
464 call remove_checksum(ncid, id_var)
465 endif
466
467 if(present(smcfcs)) then
468 error=nf90_inq_varid(ncid, "smc", id_var)
469 call netcdf_err(error, 'reading smc id' )
470 dum3d = reshape(smcfcs, (/idim,jdim,lsoil/))
471 error = nf90_put_var( ncid, id_var, dum3d)
472 call netcdf_err(error, 'writing smc record' )
473 call remove_checksum(ncid, id_var)
474 endif
475
476 if(present(stcfcs)) then
477 error=nf90_inq_varid(ncid, "stc", id_var)
478 call netcdf_err(error, 'reading stc id' )
479 dum3d = reshape(stcfcs, (/idim,jdim,lsoil/))
480 error = nf90_put_var( ncid, id_var, dum3d)
481 call netcdf_err(error, 'writing stc record' )
482 call remove_checksum(ncid, id_var)
483 endif
484
485 else
486
487 fnbgso = "./gaussian_interp." // rankch
488 print*
489 print*,"Write increments onto cubed sphere tiles to: ", trim(fnbgso)
490
491 error=nf90_create(trim(fnbgso),nf90_64bit_offset,ncid)
492 CALL netcdf_err(error, 'OPENING FILE: '//trim(fnbgso) )
493
494 ! Define dimensions in the file.
495 error = nf90_def_dim(ncid, "xaxis_1", idim, dim_x)
496 call netcdf_err(error, 'defining xaxis_1')
497
498 error = nf90_def_dim(ncid, "yaxis_1", jdim, dim_y)
499 call netcdf_err(error, 'defining yaxis_1')
500
501 error = nf90_def_dim(ncid, "soil_levels",lsoil, dim_soil)
502 call netcdf_err(error, 'defining soil_levels')
503
504 ! Define variables in the file.
505 error=nf90_def_var(ncid, "slc_inc", nf90_double, &
506 (/dim_x,dim_y,dim_soil/),varid_slc)
507 call netcdf_err(error, 'defining slc_inc');
508
509 error=nf90_def_var(ncid, "stc_inc", nf90_double, &
510 (/dim_x,dim_y,dim_soil/),varid_stc)
511 call netcdf_err(error, 'defining stc_inc');
512
513 error = nf90_enddef(ncid)
514
515 ! Put variables in the file.
516 dum3d = reshape(stcinc, (/idim,jdim,lsoil/))
517 error = nf90_put_var( ncid, varid_stc, dum3d)
518 call netcdf_err(error, 'writing stc_inc record' )
519
520 dum3d = reshape(slcinc, (/idim,jdim,lsoil/))
521 error = nf90_put_var( ncid, varid_slc, dum3d)
522 call netcdf_err(error, 'writing slc_inc record' )
523
524 endif
525
526 if(do_nsst) then
527
528 error=nf90_inq_varid(ncid, "tref", id_var)
529 call netcdf_err(error, 'reading tref id' )
530 dum2d = reshape(nsst%tref, (/idim,jdim/))
531 error = nf90_put_var( ncid, id_var, dum2d)
532 call netcdf_err(error, 'WRITING TREF RECORD' )
533 call remove_checksum(ncid, id_var)
534
535 error=nf90_inq_varid(ncid, "z_c", id_var)
536 call netcdf_err(error, 'reading z_c id' )
537 dum2d = reshape(nsst%z_c, (/idim,jdim/))
538 error = nf90_put_var( ncid, id_var, dum2d)
539 call netcdf_err(error, 'WRITING Z_C RECORD' )
540 call remove_checksum(ncid, id_var)
541
542 error=nf90_inq_varid(ncid, "c_0", id_var)
543 call netcdf_err(error, 'reading c_0 id' )
544 dum2d = reshape(nsst%c_0, (/idim,jdim/))
545 error = nf90_put_var( ncid, id_var, dum2d)
546 call netcdf_err(error, 'WRITING C_0 RECORD' )
547 call remove_checksum(ncid, id_var)
548
549 error=nf90_inq_varid(ncid, "c_d", id_var)
550 call netcdf_err(error, 'reading c_d id' )
551 dum2d = reshape(nsst%c_d, (/idim,jdim/))
552 error = nf90_put_var( ncid, id_var, dum2d)
553 call netcdf_err(error, 'WRITING C_D RECORD' )
554 call remove_checksum(ncid, id_var)
555
556 error=nf90_inq_varid(ncid, "w_0", id_var)
557 call netcdf_err(error, 'reading w_0 id' )
558 dum2d = reshape(nsst%w_0, (/idim,jdim/))
559 error = nf90_put_var( ncid, id_var, dum2d)
560 call netcdf_err(error, 'WRITING W_0 RECORD' )
561 call remove_checksum(ncid, id_var)
562
563 error=nf90_inq_varid(ncid, "w_d", id_var)
564 call netcdf_err(error, 'reading w_d id' )
565 dum2d = reshape(nsst%w_d, (/idim,jdim/))
566 error = nf90_put_var( ncid, id_var, dum2d)
567 call netcdf_err(error, 'WRITING W_D RECORD' )
568 call remove_checksum(ncid, id_var)
569
570 error=nf90_inq_varid(ncid, "xt", id_var)
571 call netcdf_err(error, 'reading xt id' )
572 dum2d = reshape(nsst%xt, (/idim,jdim/))
573 error = nf90_put_var( ncid, id_var, dum2d)
574 call netcdf_err(error, 'WRITING XT RECORD' )
575 call remove_checksum(ncid, id_var)
576
577 error=nf90_inq_varid(ncid, "xs", id_var)
578 call netcdf_err(error, 'reading xs id' )
579 dum2d = reshape(nsst%xs, (/idim,jdim/))
580 error = nf90_put_var( ncid, id_var, dum2d)
581 call netcdf_err(error, 'WRITING XS RECORD' )
582 call remove_checksum(ncid, id_var)
583
584 error=nf90_inq_varid(ncid, "xu", id_var)
585 call netcdf_err(error, 'reading xu id' )
586 dum2d = reshape(nsst%xu, (/idim,jdim/))
587 error = nf90_put_var( ncid, id_var, dum2d)
588 call netcdf_err(error, 'WRITING XU RECORD' )
589 call remove_checksum(ncid, id_var)
590
591 error=nf90_inq_varid(ncid, "xv", id_var)
592 call netcdf_err(error, 'reading xv id' )
593 dum2d = reshape(nsst%xv, (/idim,jdim/))
594 error = nf90_put_var( ncid, id_var, dum2d)
595 call netcdf_err(error, 'WRITING XV RECORD' )
596 call remove_checksum(ncid, id_var)
597
598 error=nf90_inq_varid(ncid, "xz", id_var)
599 call netcdf_err(error, 'reading xz id' )
600 dum2d = reshape(nsst%xz, (/idim,jdim/))
601 error = nf90_put_var( ncid, id_var, dum2d)
602 call netcdf_err(error, 'WRITING XZ RECORD' )
603 call remove_checksum(ncid, id_var)
604
605 error=nf90_inq_varid(ncid, "zm", id_var)
606 call netcdf_err(error, 'reading zm id' )
607 dum2d = reshape(nsst%zm, (/idim,jdim/))
608 error = nf90_put_var( ncid, id_var, dum2d)
609 call netcdf_err(error, 'WRITING ZM RECORD' )
610 call remove_checksum(ncid, id_var)
611
612 error=nf90_inq_varid(ncid, "xtts", id_var)
613 call netcdf_err(error, 'reading xtts id' )
614 dum2d = reshape(nsst%xtts, (/idim,jdim/))
615 error = nf90_put_var( ncid, id_var, dum2d)
616 call netcdf_err(error, 'WRITING XTTS RECORD' )
617 call remove_checksum(ncid, id_var)
618
619 error=nf90_inq_varid(ncid, "xzts", id_var)
620 call netcdf_err(error, 'reading xzts id' )
621 dum2d = reshape(nsst%xzts, (/idim,jdim/))
622 error = nf90_put_var( ncid, id_var, dum2d)
623 call netcdf_err(error, 'WRITING XZTS RECORD' )
624 call remove_checksum(ncid, id_var)
625
626 error=nf90_inq_varid(ncid, "d_conv", id_var)
627 call netcdf_err(error, 'reading d_conv id' )
628 dum2d = reshape(nsst%d_conv, (/idim,jdim/))
629 error = nf90_put_var( ncid, id_var, dum2d)
630 call netcdf_err(error, 'WRITING D_CONV RECORD' )
631 call remove_checksum(ncid, id_var)
632
633 error=nf90_inq_varid(ncid, "ifd", id_var)
634 call netcdf_err(error, 'reading idf id' )
635 dum2d = reshape(nsst%ifd, (/idim,jdim/))
636 error = nf90_put_var( ncid, id_var, dum2d)
637 call netcdf_err(error, 'WRITING IFD RECORD' )
638 call remove_checksum(ncid, id_var)
639
640 error=nf90_inq_varid(ncid, "dt_cool", id_var)
641 call netcdf_err(error, 'reading dt_cool id' )
642 dum2d = reshape(nsst%dt_cool, (/idim,jdim/))
643 error = nf90_put_var( ncid, id_var, dum2d)
644 call netcdf_err(error, 'WRITING DT_COOL RECORD' )
645 call remove_checksum(ncid, id_var)
646
647 error=nf90_inq_varid(ncid, "qrain", id_var)
648 call netcdf_err(error, 'reading qrain id' )
649 dum2d = reshape(nsst%qrain, (/idim,jdim/))
650 error = nf90_put_var( ncid, id_var, dum2d)
651 call netcdf_err(error, 'WRITING QRAIN RECORD' )
652 call remove_checksum(ncid, id_var)
653
654! Some files don't include 'tfinc', which is just diagnostic.
655! If missing, then add it to the restart file.
656 error=nf90_inq_varid(ncid, "tfinc", id_var)
657 if (error /= 0) then
658 error=nf90_inq_dimid(ncid, "xaxis_1", dim_x)
659 call netcdf_err(error, 'finding xaxis_1' )
660 error=nf90_inq_dimid(ncid, "yaxis_1", dim_y)
661 call netcdf_err(error, 'finding yaxis_1' )
662 error=nf90_inq_dimid(ncid, "Time", dim_time)
663 call netcdf_err(error, 'finding Time' )
664 dims_3d(1) = dim_x
665 dims_3d(2) = dim_y
666 dims_3d(3) = dim_time
667 error=nf90_redef(ncid)
668 error = nf90_def_var(ncid, 'tfinc', nf90_double, dims_3d, id_var)
669 call netcdf_err(error, 'DEFINING tfinc' )
670 error = nf90_put_att(ncid, id_var, "long_name", "tfinc")
671 call netcdf_err(error, 'DEFINING tfinc LONG NAME' )
672 error = nf90_put_att(ncid, id_var, "units", "none")
673 call netcdf_err(error, 'DEFINING tfinc UNITS' )
674 error=nf90_enddef(ncid)
675 endif
676 dum2d = reshape(nsst%tfinc, (/idim,jdim/))
677 error = nf90_put_var( ncid, id_var, dum2d)
678 call netcdf_err(error, 'WRITING TFINC RECORD' )
679
680 endif
681
682 error = nf90_close(ncid)
683
684 end subroutine write_data
685
692 subroutine remove_checksum(ncid, id_var)
693
694 implicit none
695
696 integer, intent(in) :: ncid, id_var
697
698 integer :: error
699
700 error=nf90_inquire_attribute(ncid, id_var, 'checksum')
701
702 if (error == 0) then ! attribute was found
703
704 error = nf90_redef(ncid)
705 call netcdf_err(error, 'entering define mode' )
706
707 error=nf90_del_att(ncid, id_var, 'checksum')
708 call netcdf_err(error, 'deleting checksum' )
709
710 error= nf90_enddef(ncid)
711 call netcdf_err(error, 'ending define mode' )
712
713 endif
714
715 end subroutine remove_checksum
716
732 SUBROUTINE read_lat_lon_orog(RLA,RLO,OROG,OROG_UF,&
733 TILE_NUM,IDIM,JDIM,IJDIM,LANDFRAC,LAKEFRAC)
734
735 USE mpi
736 USE machine
737
738 IMPLICIT NONE
739
740 INTEGER, INTENT(IN) :: idim, jdim, ijdim
741
742 CHARACTER(LEN=5), INTENT(OUT) :: tile_num
743
744 REAL, INTENT(OUT) :: rla(ijdim),rlo(ijdim)
745 REAL, INTENT(OUT) :: orog(ijdim),orog_uf(ijdim)
746 REAL(kind=kind_io8), INTENT(OUT), OPTIONAL :: landfrac(ijdim)
747 REAL(kind=kind_io8), INTENT(OUT), OPTIONAL :: lakefrac(ijdim)
748
749 CHARACTER(LEN=50) :: fnorog, fngrid
750 CHARACTER(LEN=3) :: rankch
751
752 INTEGER :: error, ncid, ncid_orog
753 INTEGER :: i, ii, j, jj, myrank
754 INTEGER :: id_dim, id_var, nx, ny
755
756 REAL, ALLOCATABLE :: dummy(:,:), geolat(:,:), geolon(:,:)
757 REAL(kind=4), allocatable :: dummy4(:,:)
758
759 CALL mpi_comm_rank(mpi_comm_world, myrank, error)
760
761 WRITE(rankch, '(I3.3)') (myrank+1)
762
763 fngrid = "./fngrid." // rankch
764
765 print*
766 print*, "READ FV3 GRID INFO FROM: "//trim(fngrid)
767
768 error=nf90_open(trim(fngrid),nf90_nowrite,ncid)
769 CALL netcdf_err(error, 'OPENING FILE: '//trim(fngrid) )
770
771 error=nf90_inq_dimid(ncid, 'nx', id_dim)
772 CALL netcdf_err(error, 'ERROR READING NX ID' )
773
774 error=nf90_inquire_dimension(ncid,id_dim,len=nx)
775 CALL netcdf_err(error, 'ERROR READING NX' )
776
777 error=nf90_inq_dimid(ncid, 'ny', id_dim)
778 CALL netcdf_err(error, 'ERROR READING NY ID' )
779
780 error=nf90_inquire_dimension(ncid,id_dim,len=ny)
781 CALL netcdf_err(error, 'ERROR READING NY' )
782
783 IF ((nx/2) /= idim .OR. (ny/2) /= jdim) THEN
784 print*,'FATAL ERROR: DIMENSIONS IN FILE: ',(nx/2),(ny/2)
785 print*,'DO NOT MATCH GRID DIMENSIONS: ',idim,jdim
786 CALL mpi_abort(mpi_comm_world, 130, error)
787 ENDIF
788
789 ALLOCATE(geolon(nx+1,ny+1))
790 ALLOCATE(geolat(nx+1,ny+1))
791
792 error=nf90_inq_varid(ncid, 'x', id_var)
793 CALL netcdf_err(error, 'ERROR READING X ID' )
794 error=nf90_get_var(ncid, id_var, geolon)
795 CALL netcdf_err(error, 'ERROR READING X RECORD' )
796
797 error=nf90_inq_varid(ncid, 'y', id_var)
798 CALL netcdf_err(error, 'ERROR READING Y ID' )
799 error=nf90_get_var(ncid, id_var, geolat)
800 CALL netcdf_err(error, 'ERROR READING Y RECORD' )
801
802 ALLOCATE(dummy(idim,jdim))
803
804 DO j = 1, jdim
805 DO i = 1, idim
806 ii = 2*i
807 jj = 2*j
808 dummy(i,j) = geolon(ii,jj)
809 ENDDO
810 ENDDO
811
812 rlo = reshape(dummy, (/ijdim/))
813
814 DEALLOCATE(geolon)
815
816 DO j = 1, jdim
817 DO i = 1, idim
818 ii = 2*i
819 jj = 2*j
820 dummy(i,j) = geolat(ii,jj)
821 ENDDO
822 ENDDO
823
824 rla = reshape(dummy, (/ijdim/))
825
826 DEALLOCATE(geolat, dummy)
827
828 error=nf90_inq_varid(ncid, 'tile', id_var)
829 CALL netcdf_err(error, 'ERROR READING TILE ID' )
830 error=nf90_get_var(ncid, id_var, tile_num)
831 CALL netcdf_err(error, 'ERROR READING TILE RECORD' )
832
833 error = nf90_close(ncid)
834
835 fnorog = "./fnorog." // rankch
836
837 print*
838 print*, "READ FV3 OROG INFO FROM: "//trim(fnorog)
839
840 error=nf90_open(trim(fnorog),nf90_nowrite,ncid_orog)
841 CALL netcdf_err(error, 'OPENING FILE: '//trim(fnorog) )
842
843 ALLOCATE(dummy4(idim,jdim))
844
845 error=nf90_inq_varid(ncid_orog, 'orog_raw', id_var)
846 CALL netcdf_err(error, 'ERROR READING orog_raw ID' )
847 error=nf90_get_var(ncid_orog, id_var, dummy4)
848 CALL netcdf_err(error, 'ERROR READING orog_raw RECORD' )
849 orog_uf = reshape(dummy4, (/ijdim/))
850
851 error=nf90_inq_varid(ncid_orog, 'orog_filt', id_var)
852 CALL netcdf_err(error, 'ERROR READING orog_filt ID' )
853 error=nf90_get_var(ncid_orog, id_var, dummy4)
854 CALL netcdf_err(error, 'ERROR READING orog_filt RECORD' )
855 orog = reshape(dummy4, (/ijdim/))
856
857 IF(PRESENT(landfrac))THEN
858 error=nf90_inq_varid(ncid_orog, 'land_frac', id_var)
859 CALL netcdf_err(error, 'ERROR READING land_frac ID' )
860 error=nf90_get_var(ncid_orog, id_var, dummy4)
861 CALL netcdf_err(error, 'ERROR READING land_frac RECORD' )
862 landfrac = reshape(dummy4, (/ijdim/))
863 ENDIF
864
865 IF(PRESENT(lakefrac))THEN
866 error=nf90_inq_varid(ncid_orog, 'lake_frac', id_var)
867 CALL netcdf_err(error, 'ERROR READING lake_frac ID' )
868 error=nf90_get_var(ncid_orog, id_var, dummy4)
869 CALL netcdf_err(error, 'ERROR READING lake_frac RECORD' )
870 lakefrac = reshape(dummy4, (/ijdim/))
871 ENDIF
872
873 DEALLOCATE(dummy4)
874
875 error = nf90_close(ncid_orog)
876
877 END SUBROUTINE read_lat_lon_orog
878
885 SUBROUTINE netcdf_err( ERR, STRING )
886
887 USE mpi
888
889 IMPLICIT NONE
890
891 INTEGER, INTENT(IN) :: ERR
892 CHARACTER(LEN=*), INTENT(IN) :: STRING
893 CHARACTER(LEN=80) :: ERRMSG
894 INTEGER :: IRET
895
896 IF( err == nf90_noerr )RETURN
897 errmsg = nf90_strerror(err)
898 print*,''
899 print*,'FATAL ERROR: ', trim(string), ': ', trim(errmsg)
900 print*,'STOP.'
901 CALL mpi_abort(mpi_comm_world, 999, iret)
902
903 RETURN
904 END SUBROUTINE netcdf_err
905
918 SUBROUTINE read_gsi_data(GSI_FILE, FILE_TYPE, LSOIL)
919
920 IMPLICIT NONE
921
922 CHARACTER(LEN=*), INTENT(IN) :: gsi_file
923 CHARACTER(LEN=3), INTENT(IN) :: file_type
924 INTEGER, INTENT(IN), OPTIONAL :: lsoil
925
926 INTEGER :: error, id_dim, ncid
927 INTEGER :: id_var, j
928
929 INTEGER(KIND=1), ALLOCATABLE :: idummy(:,:)
930
931 REAL(kind=8), allocatable :: dummy(:,:)
932
933 CHARACTER(LEN=1) :: k_ch
934 CHARACTER(LEN=10) :: incvar
935 CHARACTER(LEN=80) :: err_msg
936 INTEGER :: k
937
938 print*
939 print*, "READ INPUT GSI DATA FROM: "//trim(gsi_file)
940
941 error=nf90_open(trim(gsi_file),nf90_nowrite,ncid)
942 CALL netcdf_err(error, 'OPENING FILE: '//trim(gsi_file) )
943
944 error=nf90_inq_dimid(ncid, 'latitude', id_dim)
945 CALL netcdf_err(error, 'READING latitude ID' )
946 error=nf90_inquire_dimension(ncid,id_dim,len=jdim_gaus)
947 CALL netcdf_err(error, 'READING latitude length' )
948 jdim_gaus = jdim_gaus - 2 ! WILL IGNORE POLE POINTS
949
950 error=nf90_inq_dimid(ncid, 'longitude', id_dim)
951 CALL netcdf_err(error, 'READING longitude ID' )
952 error=nf90_inquire_dimension(ncid,id_dim,len=idim_gaus)
953 CALL netcdf_err(error, 'READING longitude length' )
954
955 IF (file_type=='NST') then
956 ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
958
959 error=nf90_inq_varid(ncid, "dtf", id_var)
960 CALL netcdf_err(error, 'READING dtf ID' )
961 error=nf90_get_var(ncid, id_var, dummy)
962 CALL netcdf_err(error, 'READING dtf' )
963
964 ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
966
967 error=nf90_inq_varid(ncid, "msk", id_var)
968 CALL netcdf_err(error, 'READING msk ID' )
969 error=nf90_get_var(ncid, id_var, idummy)
970 CALL netcdf_err(error, 'READING msk' )
971
972 ! REMOVE POLE POINTS.
973
974 DO j = 1, jdim_gaus
975 slmsk_gaus(:,j) = idummy(:,j+1)
976 dtref_gaus(:,j) = dummy(:,j+1)
977 ENDDO
978
979 ELSEIF (file_type=='LND') then
980
981 ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
982 ALLOCATE(stc_inc_gaus(lsoil,idim_gaus,jdim_gaus))
983 ALLOCATE(slc_inc_gaus(lsoil,idim_gaus,jdim_gaus))
984
985 ! read in soil temperature increments in each layer
986 DO k = 1, lsoil
987 WRITE(k_ch, '(I1)') k
988
989 incvar = "soilt"//k_ch//"_inc"
990 error=nf90_inq_varid(ncid, incvar, id_var)
991 err_msg = "reading "//incvar//" ID"
992 CALL netcdf_err(error, trim(err_msg))
993 error=nf90_get_var(ncid, id_var, dummy)
994 err_msg = "reading "//incvar//" data"
995 CALL netcdf_err(error, err_msg)
996
997 DO j = 1, jdim_gaus
998 stc_inc_gaus(k,:,j) = dummy(:,j+1)
999 ENDDO
1000
1001 incvar = "slc"//k_ch//"_inc"
1002 error=nf90_inq_varid(ncid, incvar, id_var)
1003 err_msg = "reading "//incvar//" ID"
1004 CALL netcdf_err(error, trim(err_msg))
1005 error=nf90_get_var(ncid, id_var, dummy)
1006 err_msg = "reading "//incvar//" data"
1007 CALL netcdf_err(error, err_msg)
1008
1009 DO j = 1, jdim_gaus
1010 slc_inc_gaus(k,:,j) = dummy(:,j+1)
1011 ENDDO
1012
1013 ENDDO
1014
1015 ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
1017
1018 error=nf90_inq_varid(ncid, "soilsnow_mask", id_var)
1019 CALL netcdf_err(error, 'READING soilsnow_mask ID' )
1020 error=nf90_get_var(ncid, id_var, idummy)
1021 CALL netcdf_err(error, 'READING soilsnow_mask' )
1022
1023 ! REMOVE POLE POINTS.
1024
1025 DO j = 1, jdim_gaus
1026 soilsnow_gaus(:,j) = idummy(:,j+1)
1027 ENDDO
1028
1029
1030 ELSE
1031 print *, 'WARNING: FILE_TYPE', file_type, 'not recognised.', &
1032 ', no increments read in'
1033 ENDIF
1034
1035 IF(ALLOCATED(dummy)) DEALLOCATE(dummy)
1036 IF(ALLOCATED(idummy)) DEALLOCATE(idummy)
1037
1038 error = nf90_close(ncid)
1039
1040 END SUBROUTINE read_gsi_data
1041
1094 SUBROUTINE read_data(LSOIL,LENSFC,DO_NSST,IS_NOAHMP, &
1095 FNAME_INC, &
1096 TSFFCS,SMCFCS,SWEFCS,STCFCS, &
1097 TG3FCS,ZORFCS, &
1098 CVFCS,CVBFCS,CVTFCS,ALBFCS, &
1099 VEGFCS,SLIFCS,CNPFCS,F10M, &
1100 VETFCS,SOTFCS,ALFFCS, &
1101 USTAR,FMM,FHH, &
1102 SIHFCS,SICFCS,SITFCS, &
1103 TPRCP,SRFLAG,SNDFCS, &
1104 VMNFCS,VMXFCS,SLCFCS, &
1105 STCINC,SLCINC,LSOIL_INCR, &
1106 SLPFCS,ABSFCS,T2M,Q2M,SLMASK, &
1107 ZSOIL,NSST)
1108 USE mpi
1109
1110 IMPLICIT NONE
1111
1112 INTEGER, INTENT(IN) :: lsoil, lensfc
1113 LOGICAL, INTENT(IN) :: do_nsst
1114
1115 CHARACTER(LEN=50), OPTIONAL, INTENT(IN) :: fname_inc
1116 INTEGER, OPTIONAL, INTENT(IN) :: lsoil_incr
1117
1118 LOGICAL, OPTIONAL, INTENT(OUT) :: is_noahmp
1119
1120 REAL, OPTIONAL, INTENT(OUT) :: cvfcs(lensfc), cvbfcs(lensfc)
1121 REAL, OPTIONAL, INTENT(OUT) :: cvtfcs(lensfc), albfcs(lensfc,4)
1122 REAL, OPTIONAL, INTENT(OUT) :: slifcs(lensfc), cnpfcs(lensfc)
1123 REAL, OPTIONAL, INTENT(OUT) :: vegfcs(lensfc), f10m(lensfc)
1124 REAL, OPTIONAL, INTENT(OUT) :: vetfcs(lensfc), sotfcs(lensfc)
1125 REAL, OPTIONAL, INTENT(OUT) :: tsffcs(lensfc), swefcs(lensfc)
1126 REAL, OPTIONAL, INTENT(OUT) :: tg3fcs(lensfc), zorfcs(lensfc)
1127 REAL, OPTIONAL, INTENT(OUT) :: alffcs(lensfc,2), ustar(lensfc)
1128 REAL, OPTIONAL, INTENT(OUT) :: fmm(lensfc), fhh(lensfc)
1129 REAL, OPTIONAL, INTENT(OUT) :: sihfcs(lensfc), sicfcs(lensfc)
1130 REAL, OPTIONAL, INTENT(OUT) :: sitfcs(lensfc), tprcp(lensfc)
1131 REAL, OPTIONAL, INTENT(OUT) :: srflag(lensfc), sndfcs(lensfc)
1132 REAL, OPTIONAL, INTENT(OUT) :: vmnfcs(lensfc), vmxfcs(lensfc)
1133 REAL, OPTIONAL, INTENT(OUT) :: slpfcs(lensfc), absfcs(lensfc)
1134 REAL, OPTIONAL, INTENT(OUT) :: t2m(lensfc), q2m(lensfc), slmask(lensfc)
1135 REAL, OPTIONAL, INTENT(OUT) :: slcfcs(lensfc,lsoil)
1136 REAL, OPTIONAL, INTENT(OUT) :: smcfcs(lensfc,lsoil)
1137 REAL, OPTIONAL, INTENT(OUT) :: stcfcs(lensfc,lsoil)
1138 REAL, OPTIONAL, INTENT(OUT) :: stcinc(lensfc,lsoil)
1139 REAL, OPTIONAL, INTENT(OUT) :: slcinc(lensfc,lsoil)
1140 REAL(kind=4), optional, INTENT(OUT) :: zsoil(lsoil)
1141
1142 TYPE(nsst_data), OPTIONAL :: nsst ! intent(out) will crash
1143 ! because subtypes are allocated in main.
1144 CHARACTER(LEN=3) :: rankch
1145 CHARACTER(LEN=50) :: fname
1146 CHARACTER(LEN=1) :: k_ch
1147 CHARACTER(LEN=10) :: incvar
1148
1149 INTEGER :: error, error2, ncid, myrank
1150 INTEGER :: idim, jdim, id_dim
1151 INTEGER :: id_var, ierr, test, k
1152
1153 LOGICAL :: jedi_incr_file
1154
1155 REAL(kind=8), allocatable :: dummy(:,:), dummy3d(:,:,:)
1156
1157 IF (PRESENT(fname_inc)) THEN
1158 fname = fname_inc
1159 ELSE
1160 CALL mpi_comm_rank(mpi_comm_world, myrank, error)
1161
1162 WRITE(rankch, '(I3.3)') (myrank+1)
1163
1164 fname = "./fnbgsi." // rankch
1165 ENDIF
1166
1167 print*
1168 print*, "READ INPUT SFC DATA FROM: "//trim(fname)
1169
1170 error=nf90_open(trim(fname),nf90_nowrite,ncid)
1171 CALL netcdf_err(error, 'OPENING FILE: '//trim(fname) )
1172
1173! Use the coordinate names to test whether this is
1174! a JEDI increment file
1175
1176 test=nf90_inq_dimid(ncid, 'xaxis_1', id_dim)
1177
1178 IF ( test == nf90_noerr ) THEN
1179 jedi_incr_file=.false.
1180
1181 error=nf90_inq_dimid(ncid, 'xaxis_1', id_dim)
1182 CALL netcdf_err(error, 'READING xaxis_1' )
1183 error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1184 CALL netcdf_err(error, 'READING xaxis_1' )
1185
1186 error=nf90_inq_dimid(ncid, 'yaxis_1', id_dim)
1187 CALL netcdf_err(error, 'READING yaxis_1' )
1188 error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1189 CALL netcdf_err(error, 'READING yaxis_1' )
1190
1191 ELSE
1192 jedi_incr_file=.true.
1193
1194 error=nf90_inq_dimid(ncid, 'grid_xt', id_dim)
1195 CALL netcdf_err(error, 'READING grid_xt' )
1196 error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1197 CALL netcdf_err(error, 'READING grid_xt' )
1198
1199 error=nf90_inq_dimid(ncid, 'grid_yt', id_dim)
1200 CALL netcdf_err(error, 'READING grid_yt' )
1201 error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1202 CALL netcdf_err(error, 'READING grid_yt' )
1203
1204 ENDIF
1205
1206 IF ((idim*jdim) /= lensfc) THEN
1207 print*,'FATAL ERROR: DIMENSIONS WRONG.'
1208 CALL mpi_abort(mpi_comm_world, 88, ierr)
1209 ENDIF
1210
1211! Check for records that indicate the restart file is
1212! for the Noah-MP land surface model.
1213
1214 IF(PRESENT(is_noahmp))THEN
1215 error=nf90_inq_varid(ncid, "canliqxy", id_var)
1216 error2=nf90_inq_varid(ncid, "tsnoxy", id_var)
1217 is_noahmp=.false.
1218 IF(error == 0 .AND. error2 == 0) THEN
1219 is_noahmp=.true.
1220 print*,"- WILL PROCESS FOR NOAH-MP LSM."
1221 ENDIF
1222 ENDIF
1223
1224 ALLOCATE(dummy(idim,jdim))
1225
1226 IF (PRESENT(tsffcs)) THEN
1227 error=nf90_inq_varid(ncid, "tsea", id_var)
1228 CALL netcdf_err(error, 'READING tsea ID' )
1229 error=nf90_get_var(ncid, id_var, dummy)
1230 CALL netcdf_err(error, 'READING tsea' )
1231 tsffcs = reshape(dummy, (/lensfc/))
1232 ENDIF
1233
1234 IF (PRESENT(swefcs)) THEN
1235 error=nf90_inq_varid(ncid, "sheleg", id_var)
1236 CALL netcdf_err(error, 'READING sheleg ID' )
1237 error=nf90_get_var(ncid, id_var, dummy)
1238 CALL netcdf_err(error, 'READING sheleg' )
1239 swefcs = reshape(dummy, (/lensfc/))
1240 ENDIF
1241
1242 IF (PRESENT(tg3fcs)) THEN
1243 error=nf90_inq_varid(ncid, "tg3", id_var)
1244 CALL netcdf_err(error, 'READING tg3 ID' )
1245 error=nf90_get_var(ncid, id_var, dummy)
1246 CALL netcdf_err(error, 'READING tg3' )
1247 tg3fcs = reshape(dummy, (/lensfc/))
1248 ENDIF
1249
1250 IF (PRESENT(zorfcs)) THEN
1251 error=nf90_inq_varid(ncid, "zorl", id_var)
1252 CALL netcdf_err(error, 'READING zorl ID' )
1253 error=nf90_get_var(ncid, id_var, dummy)
1254 CALL netcdf_err(error, 'READING zorl' )
1255 zorfcs = reshape(dummy, (/lensfc/))
1256 ENDIF
1257
1258 IF (PRESENT(albfcs)) THEN
1259
1260 error=nf90_inq_varid(ncid, "alvsf", id_var)
1261 CALL netcdf_err(error, 'READING alvsf ID' )
1262 error=nf90_get_var(ncid, id_var, dummy)
1263 CALL netcdf_err(error, 'READING alvsf' )
1264 albfcs(:,1) = reshape(dummy, (/lensfc/))
1265
1266 error=nf90_inq_varid(ncid, "alvwf", id_var)
1267 CALL netcdf_err(error, 'READING alvwf ID' )
1268 error=nf90_get_var(ncid, id_var, dummy)
1269 CALL netcdf_err(error, 'READING alvwf' )
1270 albfcs(:,2) = reshape(dummy, (/lensfc/))
1271
1272 error=nf90_inq_varid(ncid, "alnsf", id_var)
1273 CALL netcdf_err(error, 'READING alnsf ID' )
1274 error=nf90_get_var(ncid, id_var, dummy)
1275 CALL netcdf_err(error, 'READING alnsf' )
1276 albfcs(:,3) = reshape(dummy, (/lensfc/))
1277
1278 error=nf90_inq_varid(ncid, "alnwf", id_var)
1279 CALL netcdf_err(error, 'READING alnwf ID' )
1280 error=nf90_get_var(ncid, id_var, dummy)
1281 CALL netcdf_err(error, 'READING alnwf' )
1282 albfcs(:,4) = reshape(dummy, (/lensfc/))
1283
1284 ENDIF
1285
1286 IF (PRESENT(slifcs)) THEN
1287 error=nf90_inq_varid(ncid, "slmsk", id_var)
1288 CALL netcdf_err(error, 'READING slmsk ID' )
1289 error=nf90_get_var(ncid, id_var, dummy)
1290 CALL netcdf_err(error, 'READING slmsk' )
1291 slifcs = reshape(dummy, (/lensfc/))
1292 slmask = slifcs
1293 WHERE (slmask > 1.5) slmask=0.0 ! remove sea ice
1294 ENDIF
1295
1296 IF (PRESENT(cnpfcs)) THEN
1297 error=nf90_inq_varid(ncid, "canopy", id_var)
1298 CALL netcdf_err(error, 'READING canopy ID' )
1299 error=nf90_get_var(ncid, id_var, dummy)
1300 CALL netcdf_err(error, 'READING canopy' )
1301 cnpfcs = reshape(dummy, (/lensfc/))
1302 ENDIF
1303
1304 IF (PRESENT(vegfcs)) THEN
1305 error=nf90_inq_varid(ncid, "vfrac", id_var)
1306 CALL netcdf_err(error, 'READING vfrac ID' )
1307 error=nf90_get_var(ncid, id_var, dummy)
1308 CALL netcdf_err(error, 'READING vfrac' )
1309 vegfcs = reshape(dummy, (/lensfc/))
1310 ENDIF
1311
1312 IF (PRESENT(f10m)) THEN
1313 error=nf90_inq_varid(ncid, "f10m", id_var)
1314 CALL netcdf_err(error, 'READING f10m ID' )
1315 error=nf90_get_var(ncid, id_var, dummy)
1316 CALL netcdf_err(error, 'READING f10m' )
1317 f10m = reshape(dummy, (/lensfc/))
1318 ENDIF
1319
1320 IF (PRESENT(vetfcs)) THEN
1321 error=nf90_inq_varid(ncid, "vtype", id_var)
1322 CALL netcdf_err(error, 'READING vtype ID' )
1323 error=nf90_get_var(ncid, id_var, dummy)
1324 CALL netcdf_err(error, 'READING vtype' )
1325 vetfcs = reshape(dummy, (/lensfc/))
1326 ENDIF
1327
1328 IF (PRESENT(sotfcs)) THEN
1329 error=nf90_inq_varid(ncid, "stype", id_var)
1330 CALL netcdf_err(error, 'READING stype ID' )
1331 error=nf90_get_var(ncid, id_var, dummy)
1332 CALL netcdf_err(error, 'READING stype' )
1333 sotfcs = reshape(dummy, (/lensfc/))
1334 ENDIF
1335
1336 IF (PRESENT(alffcs)) THEN
1337 error=nf90_inq_varid(ncid, "facsf", id_var)
1338 CALL netcdf_err(error, 'READING facsf ID' )
1339 error=nf90_get_var(ncid, id_var, dummy)
1340 CALL netcdf_err(error, 'READING facsf' )
1341 alffcs(:,1) = reshape(dummy, (/lensfc/))
1342
1343 error=nf90_inq_varid(ncid, "facwf", id_var)
1344 CALL netcdf_err(error, 'READING facwf ID' )
1345 error=nf90_get_var(ncid, id_var, dummy)
1346 CALL netcdf_err(error, 'READING facwf' )
1347 alffcs(:,2) = reshape(dummy, (/lensfc/))
1348 ENDIF
1349
1350 IF (PRESENT(ustar)) THEN
1351 error=nf90_inq_varid(ncid, "uustar", id_var)
1352 CALL netcdf_err(error, 'READING uustar ID' )
1353 error=nf90_get_var(ncid, id_var, dummy)
1354 CALL netcdf_err(error, 'READING uustar' )
1355 ustar = reshape(dummy, (/lensfc/))
1356 ENDIF
1357
1358 IF (PRESENT(fmm)) THEN
1359 error=nf90_inq_varid(ncid, "ffmm", id_var)
1360 CALL netcdf_err(error, 'READING ffmm ID' )
1361 error=nf90_get_var(ncid, id_var, dummy)
1362 CALL netcdf_err(error, 'READING ffmm' )
1363 fmm = reshape(dummy, (/lensfc/))
1364 ENDIF
1365
1366 IF (PRESENT(fhh)) THEN
1367 error=nf90_inq_varid(ncid, "ffhh", id_var)
1368 CALL netcdf_err(error, 'READING ffhh ID' )
1369 error=nf90_get_var(ncid, id_var, dummy)
1370 CALL netcdf_err(error, 'READING ffhh' )
1371 fhh = reshape(dummy, (/lensfc/))
1372 ENDIF
1373
1374 IF (PRESENT(sihfcs)) THEN
1375 error=nf90_inq_varid(ncid, "hice", id_var)
1376 CALL netcdf_err(error, 'READING hice ID' )
1377 error=nf90_get_var(ncid, id_var, dummy)
1378 CALL netcdf_err(error, 'READING hice' )
1379 sihfcs = reshape(dummy, (/lensfc/))
1380 ENDIF
1381
1382 IF (PRESENT(sicfcs)) THEN
1383 error=nf90_inq_varid(ncid, "fice", id_var)
1384 CALL netcdf_err(error, 'READING fice ID' )
1385 error=nf90_get_var(ncid, id_var, dummy)
1386 CALL netcdf_err(error, 'READING fice' )
1387 sicfcs = reshape(dummy, (/lensfc/))
1388 ENDIF
1389
1390 IF (PRESENT(sitfcs)) THEN
1391 error=nf90_inq_varid(ncid, "tisfc", id_var)
1392 CALL netcdf_err(error, 'READING tisfc ID' )
1393 error=nf90_get_var(ncid, id_var, dummy)
1394 CALL netcdf_err(error, 'READING tisfc' )
1395 sitfcs = reshape(dummy, (/lensfc/))
1396 ENDIF
1397
1398 IF (PRESENT(tprcp)) THEN
1399 error=nf90_inq_varid(ncid, "tprcp", id_var)
1400 CALL netcdf_err(error, 'READING tprcp ID' )
1401 error=nf90_get_var(ncid, id_var, dummy)
1402 CALL netcdf_err(error, 'READING tprcp' )
1403 tprcp = reshape(dummy, (/lensfc/))
1404 ENDIF
1405
1406 IF (PRESENT(srflag)) THEN
1407 error=nf90_inq_varid(ncid, "srflag", id_var)
1408 CALL netcdf_err(error, 'READING srflag ID' )
1409 error=nf90_get_var(ncid, id_var, dummy)
1410 CALL netcdf_err(error, 'READING srflag' )
1411 srflag = reshape(dummy, (/lensfc/))
1412 ENDIF
1413
1414 IF (PRESENT(sndfcs)) THEN
1415 error=nf90_inq_varid(ncid, "snwdph", id_var)
1416 CALL netcdf_err(error, 'READING snwdph ID' )
1417 error=nf90_get_var(ncid, id_var, dummy)
1418 CALL netcdf_err(error, 'READING snwdph' )
1419 sndfcs = reshape(dummy, (/lensfc/))
1420 ENDIF
1421
1422 IF (PRESENT(vmnfcs)) THEN
1423 error=nf90_inq_varid(ncid, "shdmin", id_var)
1424 CALL netcdf_err(error, 'READING shdmin ID' )
1425 error=nf90_get_var(ncid, id_var, dummy)
1426 CALL netcdf_err(error, 'READING shdmin' )
1427 vmnfcs = reshape(dummy, (/lensfc/))
1428 ENDIF
1429
1430 IF (PRESENT(vmxfcs)) THEN
1431 error=nf90_inq_varid(ncid, "shdmax", id_var)
1432 CALL netcdf_err(error, 'READING shdmax ID' )
1433 error=nf90_get_var(ncid, id_var, dummy)
1434 CALL netcdf_err(error, 'READING shdmax' )
1435 vmxfcs = reshape(dummy, (/lensfc/))
1436 ENDIF
1437
1438 IF (PRESENT(slpfcs)) THEN
1439 error=nf90_inq_varid(ncid, "slope", id_var)
1440 CALL netcdf_err(error, 'READING slope ID' )
1441 error=nf90_get_var(ncid, id_var, dummy)
1442 CALL netcdf_err(error, 'READING slope' )
1443 slpfcs = reshape(dummy, (/lensfc/))
1444 ENDIF
1445
1446 IF (PRESENT(absfcs)) THEN
1447 error=nf90_inq_varid(ncid, "snoalb", id_var)
1448 CALL netcdf_err(error, 'READING snoalb ID' )
1449 error=nf90_get_var(ncid, id_var, dummy)
1450 CALL netcdf_err(error, 'READING snoalb' )
1451 absfcs = reshape(dummy, (/lensfc/))
1452 ENDIF
1453
1454 IF (PRESENT(t2m)) THEN
1455 error=nf90_inq_varid(ncid, "t2m", id_var)
1456 CALL netcdf_err(error, 'READING t2m ID' )
1457 error=nf90_get_var(ncid, id_var, dummy)
1458 CALL netcdf_err(error, 'READING t2m' )
1459 t2m = reshape(dummy, (/lensfc/))
1460 ENDIF
1461
1462 IF (PRESENT(q2m)) THEN
1463 error=nf90_inq_varid(ncid, "q2m", id_var)
1464 CALL netcdf_err(error, 'READING q2m ID' )
1465 error=nf90_get_var(ncid, id_var, dummy)
1466 CALL netcdf_err(error, 'READING q2m' )
1467 q2m = reshape(dummy, (/lensfc/))
1468 ENDIF
1469
1470 nsst_read : IF(do_nsst) THEN
1471
1472 print*
1473 print*,"WILL READ NSST RECORDS."
1474
1475 error=nf90_inq_varid(ncid, "c_0", id_var)
1476 CALL netcdf_err(error, 'READING c_0 ID' )
1477 error=nf90_get_var(ncid, id_var, dummy)
1478 CALL netcdf_err(error, 'READING c_0' )
1479 nsst%C_0 = reshape(dummy, (/lensfc/))
1480
1481 error=nf90_inq_varid(ncid, "c_d", id_var)
1482 CALL netcdf_err(error, 'READING c_d ID' )
1483 error=nf90_get_var(ncid, id_var, dummy)
1484 CALL netcdf_err(error, 'READING c_d' )
1485 nsst%C_D = reshape(dummy, (/lensfc/))
1486
1487 error=nf90_inq_varid(ncid, "d_conv", id_var)
1488 CALL netcdf_err(error, 'READING d_conv ID' )
1489 error=nf90_get_var(ncid, id_var, dummy)
1490 CALL netcdf_err(error, 'READING d_conv' )
1491 nsst%D_CONV = reshape(dummy, (/lensfc/))
1492
1493 error=nf90_inq_varid(ncid, "dt_cool", id_var)
1494 CALL netcdf_err(error, 'READING dt_cool ID' )
1495 error=nf90_get_var(ncid, id_var, dummy)
1496 CALL netcdf_err(error, 'READING dt_cool' )
1497 nsst%DT_COOL = reshape(dummy, (/lensfc/))
1498
1499 error=nf90_inq_varid(ncid, "ifd", id_var)
1500 CALL netcdf_err(error, 'READING ifd ID' )
1501 error=nf90_get_var(ncid, id_var, dummy)
1502 CALL netcdf_err(error, 'READING ifd' )
1503 nsst%IFD = reshape(dummy, (/lensfc/))
1504
1505 error=nf90_inq_varid(ncid, "qrain", id_var)
1506 CALL netcdf_err(error, 'READING qrain ID' )
1507 error=nf90_get_var(ncid, id_var, dummy)
1508 CALL netcdf_err(error, 'READING qrain' )
1509 nsst%QRAIN = reshape(dummy, (/lensfc/))
1510
1511 error=nf90_inq_varid(ncid, "tref", id_var)
1512 CALL netcdf_err(error, 'READING tref ID' )
1513 error=nf90_get_var(ncid, id_var, dummy)
1514 CALL netcdf_err(error, 'READING tref' )
1515 nsst%TREF = reshape(dummy, (/lensfc/))
1516
1517 error=nf90_inq_varid(ncid, "w_0", id_var)
1518 CALL netcdf_err(error, 'READING w_0 ID' )
1519 error=nf90_get_var(ncid, id_var, dummy)
1520 CALL netcdf_err(error, 'READING w_0' )
1521 nsst%W_0 = reshape(dummy, (/lensfc/))
1522
1523 error=nf90_inq_varid(ncid, "w_d", id_var)
1524 CALL netcdf_err(error, 'READING w_d ID' )
1525 error=nf90_get_var(ncid, id_var, dummy)
1526 CALL netcdf_err(error, 'READING w_d' )
1527 nsst%W_D = reshape(dummy, (/lensfc/))
1528
1529 error=nf90_inq_varid(ncid, "xs", id_var)
1530 CALL netcdf_err(error, 'READING xs ID' )
1531 error=nf90_get_var(ncid, id_var, dummy)
1532 CALL netcdf_err(error, 'READING xs' )
1533 nsst%XS = reshape(dummy, (/lensfc/))
1534
1535 error=nf90_inq_varid(ncid, "xt", id_var)
1536 CALL netcdf_err(error, 'READING xt ID' )
1537 error=nf90_get_var(ncid, id_var, dummy)
1538 CALL netcdf_err(error, 'READING xt' )
1539 nsst%XT = reshape(dummy, (/lensfc/))
1540
1541 error=nf90_inq_varid(ncid, "xtts", id_var)
1542 CALL netcdf_err(error, 'READING xtts ID' )
1543 error=nf90_get_var(ncid, id_var, dummy)
1544 CALL netcdf_err(error, 'READING xtts' )
1545 nsst%XTTS = reshape(dummy, (/lensfc/))
1546
1547 error=nf90_inq_varid(ncid, "xu", id_var)
1548 CALL netcdf_err(error, 'READING xu ID' )
1549 error=nf90_get_var(ncid, id_var, dummy)
1550 CALL netcdf_err(error, 'READING xu' )
1551 nsst%XU = reshape(dummy, (/lensfc/))
1552
1553 error=nf90_inq_varid(ncid, "xv", id_var)
1554 CALL netcdf_err(error, 'READING xv ID' )
1555 error=nf90_get_var(ncid, id_var, dummy)
1556 CALL netcdf_err(error, 'READING xv' )
1557 nsst%XV = reshape(dummy, (/lensfc/))
1558
1559 error=nf90_inq_varid(ncid, "xz", id_var)
1560 CALL netcdf_err(error, 'READING xz ID' )
1561 error=nf90_get_var(ncid, id_var, dummy)
1562 CALL netcdf_err(error, 'READING xz' )
1563 nsst%XZ = reshape(dummy, (/lensfc/))
1564
1565 error=nf90_inq_varid(ncid, "xzts", id_var)
1566 CALL netcdf_err(error, 'READING xzts ID' )
1567 error=nf90_get_var(ncid, id_var, dummy)
1568 CALL netcdf_err(error, 'READING xzts' )
1569 nsst%XZTS = reshape(dummy, (/lensfc/))
1570
1571 error=nf90_inq_varid(ncid, "z_c", id_var)
1572 CALL netcdf_err(error, 'READING z_c ID' )
1573 error=nf90_get_var(ncid, id_var, dummy)
1574 CALL netcdf_err(error, 'READING z_c' )
1575 nsst%Z_C = reshape(dummy, (/lensfc/))
1576
1577 error=nf90_inq_varid(ncid, "zm", id_var)
1578 CALL netcdf_err(error, 'READING zm ID' )
1579 error=nf90_get_var(ncid, id_var, dummy)
1580 CALL netcdf_err(error, 'READING zm' )
1581 nsst%ZM = reshape(dummy, (/lensfc/))
1582
1583 END IF nsst_read
1584
1585
1586 ALLOCATE(dummy3d(idim,jdim,lsoil))
1587
1588 IF (PRESENT(smcfcs)) THEN
1589 error=nf90_inq_varid(ncid, "smc", id_var)
1590 CALL netcdf_err(error, 'READING smc ID' )
1591 error=nf90_get_var(ncid, id_var, dummy3d)
1592 CALL netcdf_err(error, 'READING smc' )
1593 smcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1594 ENDIF
1595
1596 IF (PRESENT(slcfcs)) THEN
1597 error=nf90_inq_varid(ncid, "slc", id_var)
1598 CALL netcdf_err(error, 'READING slc ID' )
1599 error=nf90_get_var(ncid, id_var, dummy3d)
1600 CALL netcdf_err(error, 'READING slc' )
1601 slcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1602 ENDIF
1603
1604 IF (PRESENT(stcfcs)) THEN
1605 error=nf90_inq_varid(ncid, "stc", id_var)
1606 CALL netcdf_err(error, 'READING stc ID' )
1607 error=nf90_get_var(ncid, id_var, dummy3d)
1608 CALL netcdf_err(error, 'READING stc' )
1609 stcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1610 ENDIF
1611
1612! use dimension to control this one
1613 IF (jedi_incr_file) THEN
1614 IF (PRESENT(slcinc)) THEN
1615 error=nf90_inq_varid(ncid, "soill", id_var)
1616 CALL netcdf_err(error, 'READING soill ID' )
1617 error=nf90_get_var(ncid, id_var, dummy3d)
1618 CALL netcdf_err(error, 'READING slc increments' )
1619 slcinc = reshape(dummy3d, (/lensfc,lsoil/))
1620 ENDIF
1621
1622 IF (PRESENT(stcinc)) THEN
1623 error=nf90_inq_varid(ncid, "soilt", id_var)
1624 CALL netcdf_err(error, 'READING soilt ID' )
1625 error=nf90_get_var(ncid, id_var, dummy3d)
1626 CALL netcdf_err(error, 'READING stc increments' )
1627 stcinc = reshape(dummy3d, (/lensfc,lsoil/))
1628 ENDIF
1629 ELSE
1630 ! THIS IS A REGRIDDED GSI FILE
1631 IF (PRESENT(stcinc)) THEN
1632 IF (.NOT.PRESENT(lsoil_incr)) THEN
1633 write(6,*)'FATAL ERROR variable lsoil_incr not declared.'
1634 CALL mpi_abort(mpi_comm_world, 134, error)
1635 END IF
1636 DO k = 1, lsoil_incr
1637 WRITE(k_ch, '(I1)') k
1638
1639 incvar = "soilt"//k_ch//"_inc"
1640 print *, "reading", incvar
1641 error=nf90_inq_varid(ncid,trim(incvar), id_var)
1642 CALL netcdf_err(error, 'READING soilt*_inc ID')
1643 error=nf90_get_var(ncid, id_var, dummy)
1644 CALL netcdf_err(error, 'READING soilt*_inc increments')
1645
1646 stcinc(:,k) = reshape(dummy, (/lensfc/))
1647 ENDDO
1648 ENDIF
1649 IF (PRESENT(slcinc)) THEN
1650 IF (.NOT.PRESENT(lsoil_incr)) THEN
1651 write(6,*)'FATAL ERROR variable lsoil_incr not declared.'
1652 CALL mpi_abort(mpi_comm_world, 136, error)
1653 END IF
1654 DO k = 1, lsoil_incr
1655 WRITE(k_ch, '(I1)') k
1656
1657 incvar = "slc"//k_ch//"_inc"
1658 print *, "reading", trim(incvar)
1659 error=nf90_inq_varid(ncid, trim(incvar), id_var)
1660 CALL netcdf_err(error, 'READING slc*_inc ID')
1661 error=nf90_get_var(ncid, id_var, dummy)
1662 CALL netcdf_err(error, 'READING slc*_inc increments')
1663
1664 slcinc(:,k) = reshape(dummy, (/lensfc/))
1665 ENDDO
1666 ENDIF
1667 ENDIF
1668
1669 DEALLOCATE(dummy)
1670 DEALLOCATE(dummy3d)
1671
1672! cloud fields not in warm restart files. set to zero?
1673
1674 IF (PRESENT(cvfcs)) cvfcs = 0.0
1675 IF (PRESENT(cvtfcs)) cvtfcs = 0.0
1676 IF (PRESENT(cvbfcs)) cvbfcs = 0.0
1677
1678! soil layer thicknesses not in warm restart files. hardwire
1679! for now.
1680
1681 IF (PRESENT(zsoil)) THEN
1682 zsoil(1) = -0.1
1683 zsoil(2) = -0.4
1684 zsoil(3) = -1.0
1685 zsoil(4) = -2.0
1686 ENDIF
1687
1688 error = nf90_close(ncid)
1689
1690 END SUBROUTINE read_data
1691
1708subroutine read_tf_clim_grb(file_sst,sst,rlats_sst,rlons_sst,mlat_sst,mlon_sst,mon)
1709
1710 use mpi
1711
1712 implicit none
1713
1714! declare passed variables and arrays
1715 character(*) , intent(in ) :: file_sst
1716 integer , intent(in ) :: mon,mlat_sst,mlon_sst
1717 real, dimension(mlat_sst) , intent( out) :: rlats_sst
1718 real, dimension(mlon_sst) , intent( out) :: rlons_sst
1719 real, dimension(mlon_sst,mlat_sst) , intent( out) :: sst
1720
1721! declare local parameters
1722 integer,parameter:: lu_sst = 21 ! fortran unit number of grib sst file
1723 real, parameter :: deg2rad = 3.141593/180.0
1724
1725! declare local variables and arrays
1726 logical(1), allocatable, dimension(:) :: lb
1727
1728 integer :: nlat_sst
1729 integer :: nlon_sst
1730 integer :: iret,ni,nj
1731 integer :: mscan,kb1,ierr
1732 integer :: jincdir,i,iincdir,kb2,kb3,kf,kg,k,j,jf
1733 integer, dimension(22):: jgds,kgds
1734 integer, dimension(25):: jpds,kpds
1735
1736 real :: xsst0
1737 real :: ysst0
1738 real :: dres
1739 real, allocatable, dimension(:) :: f
1740
1741!************+******************************************************************************
1742!
1743! open sst analysis file (grib)
1744 write(*,*) ' sstclm : ',file_sst
1745 call baopenr(lu_sst,trim(file_sst),iret)
1746 if (iret /= 0 ) then
1747 write(6,*)'FATAL ERROR in read_tf_clm_grb: error opening sst file.'
1748 CALL mpi_abort(mpi_comm_world, 111, ierr)
1749 endif
1750
1751! define sst variables for read
1752 j=-1
1753 jpds=-1
1754 jgds=-1
1755 jpds(5)=11 ! sst variable
1756 jpds(6)=1 ! surface
1757 jpds(9) = mon
1758 call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1759
1760 nlat_sst = kgds(3) ! number points on longitude circle (360)
1761 nlon_sst = kgds(2) ! number points on latitude circle (720)
1762
1763! write(*,*) 'nlat_sst, nlon_sst, mon : ',nlat_sst, nlon_sst, mon
1764! write(*,*) 'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1765
1766! allocate arrays
1767 allocate(lb(nlat_sst*nlon_sst))
1768 allocate(f(nlat_sst*nlon_sst))
1769 jf=nlat_sst*nlon_sst
1770
1771! read in the analysis
1772 call getgb(lu_sst,0,jf,j,jpds,jgds,kf,k,kpds,kgds,lb,f,iret)
1773 if (iret /= 0) then
1774 write(6,*)'FATAL ERROR in read_tf_clm_grb: error reading sst analysis data record.'
1775 deallocate(lb,f)
1776 CALL mpi_abort(mpi_comm_world, 111, ierr)
1777 endif
1778
1779 if ( (nlat_sst /= mlat_sst) .or. (nlon_sst /= mlon_sst) ) then
1780 write(6,*)'FATAL ERROR in read_rtg_org: inconsistent dimensions. mlat_sst,mlon_sst=',&
1781 mlat_sst,mlon_sst,' -versus- nlat_sst,nlon_sst=',nlat_sst,nlon_sst
1782 deallocate(lb,f)
1783 CALL mpi_abort(mpi_comm_world, 111, ierr)
1784 endif
1785
1786!
1787! get xlats and xlons
1788!
1789 dres = 180.0/real(nlat_sst)
1790 ysst0 = 0.5*dres-90.0
1791 xsst0 = 0.5*dres
1792
1793! get lat_sst & lon_sst
1794 do j = 1, nlat_sst
1795 rlats_sst(j) = ysst0 + real(j-1)*dres
1796 enddo
1797
1798 do i = 1, nlon_sst
1799 rlons_sst(i) = (xsst0 + real(i-1)*dres)
1800 enddo
1801
1802! load dimensions and grid specs. check for unusual values
1803 ni=kgds(2) ! 720 for 0.5 x 0.5
1804 nj=kgds(3) ! 360 for 0.5 x 0.5 resolution
1805
1806 mscan=kgds(11)
1807 kb1=ibits(mscan,7,1) ! i scan direction
1808 kb2=ibits(mscan,6,1) ! j scan direction
1809 kb3=ibits(mscan,5,1) ! (i,j) or (j,i)
1810
1811! get i and j scanning directions from kb1 and kb2.
1812! 0 yields +1, 1 yields -1. +1 is west to east, -1 is east to west.
1813 iincdir = 1-kb1*2
1814
1815! 0 yields -1, 1 yields +1. +1 is south to north, -1 is north to south.
1816 jincdir = kb2*2 - 1
1817
1818! write(*,*) 'read_tf_clim_grb,iincdir,jincdir : ',iincdir,jincdir
1819 do k=1,kf
1820
1821! kb3 from scan mode indicates if i points are consecutive
1822! or if j points are consecutive
1823 if(kb3==0)then ! (i,j)
1824 i=(ni+1)*kb1+(mod(k-1,ni)+1)*iincdir
1825 j=(nj+1)*(1-kb2)+jincdir*((k-1)/ni+1)
1826 else ! (j,i)
1827 j=(nj+1)*(1-kb2)+(mod(k-1,nj)+1)*jincdir
1828 i=(ni+1)*kb1+iincdir*((k-1)/nj+1)
1829 endif
1830 sst(i,j) = f(k)
1831 end do
1832
1833 deallocate(lb,f)
1834
1835 call baclose(lu_sst,iret)
1836 if (iret /= 0 ) then
1837 write(6,*)'FATAL ERROR in read_tf_clm_grb: error closing sst file.'
1838 CALL mpi_abort(mpi_comm_world, 121, ierr)
1839 endif
1840
1841end subroutine read_tf_clim_grb
1842
1850subroutine get_tf_clm_dim(file_sst,mlat_sst,mlon_sst)
1851 use mpi
1852
1853 implicit none
1854
1855! declare passed variables and arrays
1856 character(*) , intent(in ) :: file_sst
1857 integer , intent(out) :: mlat_sst,mlon_sst
1858
1859! declare local parameters
1860 integer,parameter:: lu_sst = 21 ! fortran unit number of grib sst file
1861
1862 integer :: iret
1863 integer :: kf,kg,k,j,ierr
1864 integer, dimension(22):: jgds,kgds
1865 integer, dimension(25):: jpds,kpds
1866
1867!************+******************************************************************************
1868!
1869! open sst analysis file (grib)
1870 call baopenr(lu_sst,trim(file_sst),iret)
1871 if (iret /= 0 ) then
1872 write(6,*)'FATAL ERROR in get_tf_clm_dim: error opening sst file.'
1873 CALL mpi_abort(mpi_comm_world, 111, ierr)
1874 endif
1875
1876! define sst variables for read
1877 j=-1
1878 jpds=-1
1879 jgds=-1
1880 jpds(5)=11 ! sst variable
1881 jpds(6)=1 ! surface
1882 jpds(9) = 1
1883 call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1884
1885 mlat_sst = kgds(3) ! number points on longitude circle (360)
1886 mlon_sst = kgds(2) ! number points on latitude circle (720)
1887
1888 write(*,*) 'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1889
1890 call baclose(lu_sst,iret)
1891 if (iret /= 0 ) then
1892 write(6,*)'FATAL ERROR in get_tf_clm_dim: error closing sst file.'
1893 CALL mpi_abort(mpi_comm_world, 121, ierr)
1894 endif
1895end subroutine get_tf_clm_dim
1896
1908subroutine read_salclm_gfs_nc(filename,sal,xlats,xlons,nlat,nlon,itime)
1909 use netcdf
1910 implicit none
1911
1912 ! This is the name of the data file we will read.
1913 character (len=*), intent(in) :: filename
1914 integer, intent(in) :: nlat,nlon
1915 integer, intent(in) :: itime
1916 real, dimension(nlat), intent(out) :: xlats
1917 real, dimension(nlon), intent(out) :: xlons
1918 real, dimension(nlon,nlat), intent(out) :: sal
1919! Local variables
1920 integer :: ncid
1921
1922 integer, parameter :: ndims = 3
1923 character (len = *), parameter :: lat_name = "latitude"
1924 character (len = *), parameter :: lon_name = "longitude"
1925 character (len = *), parameter :: t_name = "time"
1926 character (len = *), parameter :: sal_name="sal"
1927 integer :: time_varid,lon_varid, lat_varid, sal_varid
1928
1929 ! The start and count arrays will tell the netCDF library where to read our data.
1930 integer, dimension(ndims) :: start, count
1931
1932 character (len = *), parameter :: units = "units"
1933 character (len = *), parameter :: sal_units = "psu"
1934 ! PSU (Practical SalinitUnit). 1 PSU = 1g/kg
1935 character (len = *), parameter :: time_units = "months"
1936 character (len = *), parameter :: lat_units = "degrees_north"
1937 character (len = *), parameter :: lon_units = "degrees_east"
1938
1939! Open the file.
1940 call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1941
1942! Get the varids of time, latitude, longitude & depth coordinate variables.
1943 call nc_check( nf90_inq_varid(ncid, t_name, time_varid) )
1944 call nc_check( nf90_inq_varid(ncid, lat_name, lat_varid) )
1945 call nc_check( nf90_inq_varid(ncid, lon_name, lon_varid) )
1946
1947! Read the time, latitude and longitude data.
1948! call nc_check( nf90_get_var(ncid, time_varid, ntime) )
1949 call nc_check( nf90_get_var(ncid, lat_varid, xlats) )
1950 call nc_check( nf90_get_var(ncid, lon_varid, xlons) )
1951
1952! Get the varids of the sal netCDF variables.
1953 call nc_check( nf90_inq_varid(ncid, sal_name,sal_varid) )
1954
1955! Read 1 record of nlat*nlon values, starting at the beginning
1956! of the record (the (1, 1, 1, rec) element in the netCDF file).
1957 start = (/ 1, 1, itime /)
1958 count = (/ nlon, nlat, 1 /)
1959
1960! write(*,*) 'read_salclm_gfs_nc itime : ',itime
1961! Read the sal data from the file, one record at a time.
1962 call nc_check( nf90_get_var(ncid, sal_varid, sal, start, count) )
1963
1964! Close the file. This frees up any internal netCDF resources
1965! associated with the file.
1966 call nc_check( nf90_close(ncid) )
1967
1968! If we got this far, everything worked as expected. Yipee!
1969! print *,"*** SUCCESS reading file ", filename, "!"
1970
1971end subroutine read_salclm_gfs_nc
1972
1979subroutine get_dim_nc(filename,nlat,nlon)
1980 use netcdf
1981 implicit none
1982
1983 character (len=*), intent(in) :: filename
1984 integer, intent(out) :: nlat,nlon
1985! Local variables
1986 character (len = *), parameter :: lat_name = "latitude"
1987 character (len = *), parameter :: lon_name = "longitude"
1988 integer :: ncid
1989 integer :: latdimid,londimid
1990
1991! Open the file.
1992 call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1993
1994! Get dimensions
1995 call nc_check( nf90_inq_dimid(ncid,lat_name,latdimid) )
1996 call nc_check( nf90_inq_dimid(ncid,lon_name,londimid) )
1997 call nc_check( nf90_inquire_dimension(ncid,latdimid,len=nlat) )
1998 call nc_check( nf90_inquire_dimension(ncid,londimid,len=nlon) )
1999
2000! write(*,'(a,1x,a6,2I8)') 'get_dim_nc, file, nlat, nlon : ',filename,nlat,nlon
2001
2002! Close the file. This frees up any internal netCDF resources
2003! associated with the file.
2004 call nc_check( nf90_close(ncid) )
2005
2006! If we got this far, everything worked as expected. Yipee!
2007! print *,"*** SUCCESS get dimensions from nc file ", filename, "!"
2008
2009end subroutine get_dim_nc
2010
2016subroutine nc_check(status)
2017
2018 use mpi
2019 use netcdf
2020
2021 integer, intent ( in) :: status
2022 integer :: ierr
2023
2024 if(status /= nf90_noerr) then
2025 print *, 'FATAL ERROR:'
2026 print *, trim(nf90_strerror(status))
2027 CALL mpi_abort(mpi_comm_world, 122, ierr)
2028 end if
2029end subroutine nc_check
2030
2031 END MODULE read_write_data
Holds machine dependent constants for global_cycle.
Definition machine.f90:7
integer, parameter kind_io8
Kind type for 8-byte float point variables.
Definition machine.f90:17
This module contains routines that read and write data.
integer, public jdim_gaus
'j' dimension of GSI gaussian grid.
subroutine, public get_tf_clm_dim(file_sst, mlat_sst, mlon_sst)
Get the i/j dimensions of RTG SST climatology file.
subroutine, public read_data(lsoil, lensfc, do_nsst, is_noahmp, fname_inc, tsffcs, smcfcs, swefcs, stcfcs, tg3fcs, zorfcs, cvfcs, cvbfcs, cvtfcs, albfcs, vegfcs, slifcs, cnpfcs, f10m, vetfcs, sotfcs, alffcs, ustar, fmm, fhh, sihfcs, sicfcs, sitfcs, tprcp, srflag, sndfcs, vmnfcs, vmxfcs, slcfcs, stcinc, slcinc, lsoil_incr, slpfcs, absfcs, t2m, q2m, slmask, zsoil, nsst)
Read the first guess surface records and nsst records (if selected) for a single cubed-sphere tile.
integer, dimension(:,:), allocatable, public soilsnow_gaus
GSI soil / snow mask for land on the gaussian grid.
real, dimension(:,:,:), allocatable, public slc_inc_gaus
GSI soil moisture increments on the gaussian grid.
subroutine, public read_lat_lon_orog(rla, rlo, orog, orog_uf, tile_num, idim, jdim, ijdim, landfrac, lakefrac)
Read latitude and longitude for the cubed-sphere tile from the 'grid' file.
real, dimension(:,:,:), allocatable, public stc_inc_gaus
GSI soil temperature increments on the gaussian grid.
subroutine, public read_tf_clim_grb(file_sst, sst, rlats_sst, rlons_sst, mlat_sst, mlon_sst, mon)
Read a GRIB1 sst climatological analysis file.
subroutine, public write_data(lensfc, idim, jdim, lsoil, do_nsst, inc_file, nsst, slifcs, tsffcs, vegfcs, swefcs, tg3fcs, zorfcs, albfcs, alffcs, cnpfcs, f10m, t2m, q2m, vetfcs, sotfcs, ustar, fmm, fhh, sicfcs, sihfcs, sitfcs, tprcp, srflag, swdfcs, vmnfcs, vmxfcs, slpfcs, absfcs, slcfcs, smcfcs, stcfcs, stcinc, slcinc)
Update surface records - and nsst records if selected - on a single cubed-sphere tile to a pre-existi...
subroutine, public read_gsi_data(gsi_file, file_type, lsoil)
Read increment file from the GSI containing either the foundation temperature increments and mask,...
subroutine, public get_dim_nc(filename, nlat, nlon)
Get the i/j dimensions of the data from a NetCDF file.
real, dimension(:,:), allocatable, public dtref_gaus
GSI foundation temperature increment on the gaussian grid.
subroutine nc_check(status)
Check the NetCDF status code.
integer, dimension(:,:), allocatable, public slmsk_gaus
GSI land mask on the gaussian grid.
subroutine, public read_salclm_gfs_nc(filename, sal, xlats, xlons, nlat, nlon, itime)
Read the woa05 salinity monthly climatology file.
subroutine netcdf_err(err, string)
If a NetCDF call returns an error, print out a user-supplied message and the NetCDF library message.
subroutine remove_checksum(ncid, id_var)
Remove the checksum attribute from a netcdf record.
integer, public idim_gaus
'i' dimension of GSI gaussian grid.