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(:)
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, &
137 integer,
intent(in) :: lensfc, lsoil
138 integer,
intent(in) :: idim, jdim
140 logical,
intent(in) :: do_nsst
141 logical,
intent(in) :: inc_file
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)
163 integer :: dim_x, dim_y, dim_soil, dim_time, dims_3d(3)
165 real :: dum2d(idim,jdim), dum3d(idim,jdim,lsoil)
167 character(len=50) :: fnbgso
168 character(len=3) :: rankch
170 integer :: myrank, error, ncid, id_var
171 integer :: varid_stc, varid_slc
173 call mpi_comm_rank(mpi_comm_world, myrank, error)
175 write(rankch,
'(i3.3)') (myrank+1)
177 if (.NOT.(inc_file))
then 179 fnbgso =
"./fnbgso." // rankch
182 print*,
"update OUTPUT SFC DATA TO: ",trim(fnbgso)
184 error=nf90_open(trim(fnbgso),nf90_write,ncid)
185 CALL netcdf_err(error,
'OPENING FILE: '//trim(fnbgso) )
187 if(
present(slifcs))
then 188 error=nf90_inq_varid(ncid,
"slmsk", id_var)
190 dum2d = reshape(slifcs, (/idim,jdim/))
191 error = nf90_put_var( ncid, id_var, dum2d)
192 call netcdf_err(error,
'writing slmsk record' )
196 if(
present(tsffcs))
then 197 error=nf90_inq_varid(ncid,
"tsea", id_var)
199 dum2d = reshape(tsffcs, (/idim,jdim/))
200 error = nf90_put_var( ncid, id_var, dum2d)
201 call netcdf_err(error,
'writing tsea record' )
205 if(
present(swefcs))
then 206 error=nf90_inq_varid(ncid,
"sheleg", id_var)
208 dum2d = reshape(swefcs, (/idim,jdim/))
209 error = nf90_put_var( ncid, id_var, dum2d)
210 call netcdf_err(error,
'writing sheleg record' )
214 if(
present(tg3fcs))
then 215 error=nf90_inq_varid(ncid,
"tg3", id_var)
217 dum2d = reshape(tg3fcs, (/idim,jdim/))
218 error = nf90_put_var( ncid, id_var, dum2d)
223 if(
present(zorfcs))
then 224 error=nf90_inq_varid(ncid,
"zorl", id_var)
226 dum2d = reshape(zorfcs, (/idim,jdim/))
227 error = nf90_put_var( ncid, id_var, dum2d)
228 call netcdf_err(error,
'writing zorl record' )
232 if(
present(albfcs))
then 233 error=nf90_inq_varid(ncid,
"alvsf", id_var)
235 dum2d = reshape(albfcs(:,1), (/idim,jdim/))
236 error = nf90_put_var( ncid, id_var, dum2d)
237 call netcdf_err(error,
'writing alvsf record' )
240 error=nf90_inq_varid(ncid,
"alvwf", id_var)
242 dum2d = reshape(albfcs(:,2), (/idim,jdim/))
243 error = nf90_put_var( ncid, id_var, dum2d)
244 call netcdf_err(error,
'writing alvwf record' )
247 error=nf90_inq_varid(ncid,
"alnsf", id_var)
249 dum2d = reshape(albfcs(:,3), (/idim,jdim/))
250 error = nf90_put_var( ncid, id_var, dum2d)
251 call netcdf_err(error,
'writing alnsf record' )
254 error=nf90_inq_varid(ncid,
"alnwf", id_var)
256 dum2d = reshape(albfcs(:,4), (/idim,jdim/))
257 error = nf90_put_var( ncid, id_var, dum2d)
258 call netcdf_err(error,
'writing alnwf record' )
262 if(
present(alffcs))
then 263 error=nf90_inq_varid(ncid,
"facsf", id_var)
265 dum2d = reshape(alffcs(:,1), (/idim,jdim/))
266 error = nf90_put_var( ncid, id_var, dum2d)
267 call netcdf_err(error,
'writing facsf record' )
270 error=nf90_inq_varid(ncid,
"facwf", id_var)
272 dum2d = reshape(alffcs(:,2), (/idim,jdim/))
273 error = nf90_put_var( ncid, id_var, dum2d)
274 call netcdf_err(error,
'writing facwf record' )
278 if(
present(vegfcs))
then 279 error=nf90_inq_varid(ncid,
"vfrac", id_var)
281 dum2d = reshape(vegfcs, (/idim,jdim/))
282 error = nf90_put_var( ncid, id_var, dum2d)
283 call netcdf_err(error,
'writing vegfcs record' )
287 if(
present(cnpfcs))
then 288 error=nf90_inq_varid(ncid,
"canopy", id_var)
290 dum2d = reshape(cnpfcs, (/idim,jdim/))
291 error = nf90_put_var( ncid, id_var, dum2d)
292 call netcdf_err(error,
'writing canopy record' )
296 if(
present(f10m))
then 297 error=nf90_inq_varid(ncid,
"f10m", id_var)
299 dum2d = reshape(f10m, (/idim,jdim/))
300 error = nf90_put_var( ncid, id_var, dum2d)
301 call netcdf_err(error,
'writing f10m record' )
305 if(
present(t2m))
then 306 error=nf90_inq_varid(ncid,
"t2m", id_var)
308 dum2d = reshape(t2m, (/idim,jdim/))
309 error = nf90_put_var( ncid, id_var, dum2d)
314 if(
present(q2m))
then 315 error=nf90_inq_varid(ncid,
"q2m", id_var)
317 dum2d = reshape(q2m, (/idim,jdim/))
318 error = nf90_put_var( ncid, id_var, dum2d)
323 if(
present(vetfcs))
then 324 error=nf90_inq_varid(ncid,
"vtype", id_var)
326 dum2d = reshape(vetfcs, (/idim,jdim/))
327 error = nf90_put_var( ncid, id_var, dum2d)
328 call netcdf_err(error,
'writing vtype record' )
332 if(
present(sotfcs))
then 333 error=nf90_inq_varid(ncid,
"stype", id_var)
335 dum2d = reshape(sotfcs, (/idim,jdim/))
336 error = nf90_put_var( ncid, id_var, dum2d)
337 call netcdf_err(error,
'writing stype record' )
341 if(
present(ustar))
then 342 error=nf90_inq_varid(ncid,
"uustar", id_var)
344 dum2d = reshape(ustar, (/idim,jdim/))
345 error = nf90_put_var( ncid, id_var, dum2d)
346 call netcdf_err(error,
'writing uustar record' )
350 if(
present(fmm))
then 351 error=nf90_inq_varid(ncid,
"ffmm", id_var)
353 dum2d = reshape(fmm, (/idim,jdim/))
354 error = nf90_put_var( ncid, id_var, dum2d)
355 call netcdf_err(error,
'writing ffmm record' )
359 if(
present(fhh))
then 360 error=nf90_inq_varid(ncid,
"ffhh", id_var)
362 dum2d = reshape(fhh, (/idim,jdim/))
363 error = nf90_put_var( ncid, id_var, dum2d)
364 call netcdf_err(error,
'writing ffhh record' )
368 if(
present(sicfcs))
then 369 error=nf90_inq_varid(ncid,
"fice", id_var)
371 dum2d = reshape(sicfcs, (/idim,jdim/))
372 error = nf90_put_var( ncid, id_var, dum2d)
373 call netcdf_err(error,
'writing fice record' )
377 if(
present(sihfcs))
then 378 error=nf90_inq_varid(ncid,
"hice", id_var)
380 dum2d = reshape(sihfcs, (/idim,jdim/))
381 error = nf90_put_var( ncid, id_var, dum2d)
382 call netcdf_err(error,
'writing hice record' )
386 if(
present(sitfcs))
then 387 error=nf90_inq_varid(ncid,
"tisfc", id_var)
389 dum2d = reshape(sitfcs, (/idim,jdim/))
390 error = nf90_put_var( ncid, id_var, dum2d)
391 call netcdf_err(error,
'writing tisfc record' )
395 if(
present(tprcp))
then 396 error=nf90_inq_varid(ncid,
"tprcp", id_var)
398 dum2d = reshape(tprcp, (/idim,jdim/))
399 error = nf90_put_var( ncid, id_var, dum2d)
400 call netcdf_err(error,
'writing tprcp record' )
404 if(
present(srflag))
then 405 error=nf90_inq_varid(ncid,
"srflag", id_var)
407 dum2d = reshape(srflag, (/idim,jdim/))
408 error = nf90_put_var( ncid, id_var, dum2d)
409 call netcdf_err(error,
'writing srflag record' )
413 if(
present(swdfcs))
then 414 error=nf90_inq_varid(ncid,
"snwdph", id_var)
416 dum2d = reshape(swdfcs, (/idim,jdim/))
417 error = nf90_put_var( ncid, id_var, dum2d)
418 call netcdf_err(error,
'writing snwdph record' )
422 if(
present(vmnfcs))
then 423 error=nf90_inq_varid(ncid,
"shdmin", id_var)
425 dum2d = reshape(vmnfcs, (/idim,jdim/))
426 error = nf90_put_var( ncid, id_var, dum2d)
427 call netcdf_err(error,
'writing shdmin record' )
431 if(
present(vmxfcs))
then 432 error=nf90_inq_varid(ncid,
"shdmax", id_var)
434 dum2d = reshape(vmxfcs, (/idim,jdim/))
435 error = nf90_put_var( ncid, id_var, dum2d)
436 call netcdf_err(error,
'writing shdmax record' )
440 if(
present(slpfcs))
then 441 error=nf90_inq_varid(ncid,
"slope", id_var)
443 dum2d = reshape(slpfcs, (/idim,jdim/))
444 error = nf90_put_var( ncid, id_var, dum2d)
445 call netcdf_err(error,
'writing slope record' )
449 if(
present(absfcs))
then 450 error=nf90_inq_varid(ncid,
"snoalb", id_var)
452 dum2d = reshape(absfcs, (/idim,jdim/))
453 error = nf90_put_var( ncid, id_var, dum2d)
454 call netcdf_err(error,
'writing snoalb record' )
458 if(
present(slcfcs))
then 459 error=nf90_inq_varid(ncid,
"slc", id_var)
461 dum3d = reshape(slcfcs, (/idim,jdim,lsoil/))
462 error = nf90_put_var( ncid, id_var, dum3d)
467 if(
present(smcfcs))
then 468 error=nf90_inq_varid(ncid,
"smc", id_var)
470 dum3d = reshape(smcfcs, (/idim,jdim,lsoil/))
471 error = nf90_put_var( ncid, id_var, dum3d)
476 if(
present(stcfcs))
then 477 error=nf90_inq_varid(ncid,
"stc", id_var)
479 dum3d = reshape(stcfcs, (/idim,jdim,lsoil/))
480 error = nf90_put_var( ncid, id_var, dum3d)
487 fnbgso =
"./gaussian_interp." // rankch
489 print*,
"Write increments onto cubed sphere tiles to: ", trim(fnbgso)
491 error=nf90_create(trim(fnbgso),nf90_64bit_offset,ncid)
492 CALL netcdf_err(error,
'OPENING FILE: '//trim(fnbgso) )
495 error = nf90_def_dim(ncid,
"xaxis_1", idim, dim_x)
498 error = nf90_def_dim(ncid,
"yaxis_1", jdim, dim_y)
501 error = nf90_def_dim(ncid,
"soil_levels",lsoil, dim_soil)
502 call netcdf_err(error,
'defining soil_levels')
505 error=nf90_def_var(ncid,
"slc_inc", nf90_double, &
506 (/dim_x,dim_y,dim_soil/),varid_slc)
509 error=nf90_def_var(ncid,
"stc_inc", nf90_double, &
510 (/dim_x,dim_y,dim_soil/),varid_stc)
513 error = nf90_enddef(ncid)
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' )
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' )
528 error=nf90_inq_varid(ncid,
"tref", id_var)
530 dum2d = reshape(nsst%tref, (/idim,jdim/))
531 error = nf90_put_var( ncid, id_var, dum2d)
532 call netcdf_err(error,
'WRITING TREF RECORD' )
535 error=nf90_inq_varid(ncid,
"z_c", id_var)
537 dum2d = reshape(nsst%z_c, (/idim,jdim/))
538 error = nf90_put_var( ncid, id_var, dum2d)
542 error=nf90_inq_varid(ncid,
"c_0", id_var)
544 dum2d = reshape(nsst%c_0, (/idim,jdim/))
545 error = nf90_put_var( ncid, id_var, dum2d)
549 error=nf90_inq_varid(ncid,
"c_d", id_var)
551 dum2d = reshape(nsst%c_d, (/idim,jdim/))
552 error = nf90_put_var( ncid, id_var, dum2d)
556 error=nf90_inq_varid(ncid,
"w_0", id_var)
558 dum2d = reshape(nsst%w_0, (/idim,jdim/))
559 error = nf90_put_var( ncid, id_var, dum2d)
563 error=nf90_inq_varid(ncid,
"w_d", id_var)
565 dum2d = reshape(nsst%w_d, (/idim,jdim/))
566 error = nf90_put_var( ncid, id_var, dum2d)
570 error=nf90_inq_varid(ncid,
"xt", id_var)
572 dum2d = reshape(nsst%xt, (/idim,jdim/))
573 error = nf90_put_var( ncid, id_var, dum2d)
577 error=nf90_inq_varid(ncid,
"xs", id_var)
579 dum2d = reshape(nsst%xs, (/idim,jdim/))
580 error = nf90_put_var( ncid, id_var, dum2d)
584 error=nf90_inq_varid(ncid,
"xu", id_var)
586 dum2d = reshape(nsst%xu, (/idim,jdim/))
587 error = nf90_put_var( ncid, id_var, dum2d)
591 error=nf90_inq_varid(ncid,
"xv", id_var)
593 dum2d = reshape(nsst%xv, (/idim,jdim/))
594 error = nf90_put_var( ncid, id_var, dum2d)
598 error=nf90_inq_varid(ncid,
"xz", id_var)
600 dum2d = reshape(nsst%xz, (/idim,jdim/))
601 error = nf90_put_var( ncid, id_var, dum2d)
605 error=nf90_inq_varid(ncid,
"zm", id_var)
607 dum2d = reshape(nsst%zm, (/idim,jdim/))
608 error = nf90_put_var( ncid, id_var, dum2d)
612 error=nf90_inq_varid(ncid,
"xtts", id_var)
614 dum2d = reshape(nsst%xtts, (/idim,jdim/))
615 error = nf90_put_var( ncid, id_var, dum2d)
616 call netcdf_err(error,
'WRITING XTTS RECORD' )
619 error=nf90_inq_varid(ncid,
"xzts", id_var)
621 dum2d = reshape(nsst%xzts, (/idim,jdim/))
622 error = nf90_put_var( ncid, id_var, dum2d)
623 call netcdf_err(error,
'WRITING XZTS RECORD' )
626 error=nf90_inq_varid(ncid,
"d_conv", id_var)
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' )
633 error=nf90_inq_varid(ncid,
"ifd", id_var)
635 dum2d = reshape(nsst%ifd, (/idim,jdim/))
636 error = nf90_put_var( ncid, id_var, dum2d)
640 error=nf90_inq_varid(ncid,
"dt_cool", id_var)
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' )
647 error=nf90_inq_varid(ncid,
"qrain", id_var)
649 dum2d = reshape(nsst%qrain, (/idim,jdim/))
650 error = nf90_put_var( ncid, id_var, dum2d)
651 call netcdf_err(error,
'WRITING QRAIN RECORD' )
656 error=nf90_inq_varid(ncid,
"tfinc", id_var)
658 error=nf90_inq_dimid(ncid,
"xaxis_1", dim_x)
660 error=nf90_inq_dimid(ncid,
"yaxis_1", dim_y)
662 error=nf90_inq_dimid(ncid,
"Time", dim_time)
666 dims_3d(3) = dim_time
667 error=nf90_redef(ncid)
668 error = nf90_def_var(ncid,
'tfinc', nf90_double, dims_3d, id_var)
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)
676 dum2d = reshape(nsst%tfinc, (/idim,jdim/))
677 error = nf90_put_var( ncid, id_var, dum2d)
678 call netcdf_err(error,
'WRITING TFINC RECORD' )
682 error = nf90_close(ncid)
696 integer,
intent(in) :: ncid, id_var
700 error=nf90_inquire_attribute(ncid, id_var,
'checksum')
704 error = nf90_redef(ncid)
705 call netcdf_err(error,
'entering define mode' )
707 error=nf90_del_att(ncid, id_var,
'checksum')
710 error= nf90_enddef(ncid)
732 TILE_NUM,IDIM,JDIM,IJDIM,LANDFRAC)
739 INTEGER,
INTENT(IN) :: idim, jdim, ijdim
741 CHARACTER(LEN=5),
INTENT(OUT) :: tile_num
743 REAL,
INTENT(OUT) :: rla(ijdim),rlo(ijdim)
744 REAL,
INTENT(OUT) :: orog(ijdim),orog_uf(ijdim)
745 REAL(KIND=KIND_IO8),
INTENT(OUT),
OPTIONAL :: landfrac(ijdim)
747 CHARACTER(LEN=50) :: fnorog, fngrid
748 CHARACTER(LEN=3) :: rankch
750 INTEGER :: error, ncid, ncid_orog
751 INTEGER :: i, ii, j, jj, myrank
752 INTEGER :: id_dim, id_var, nx, ny
754 REAL,
ALLOCATABLE :: dummy(:,:), geolat(:,:), geolon(:,:)
755 REAL(KIND=4),
ALLOCATABLE :: dummy4(:,:)
757 CALL mpi_comm_rank(mpi_comm_world, myrank, error)
759 WRITE(rankch,
'(I3.3)') (myrank+1)
761 fngrid =
"./fngrid." // rankch
764 print*,
"READ FV3 GRID INFO FROM: "//trim(fngrid)
766 error=nf90_open(trim(fngrid),nf90_nowrite,ncid)
767 CALL netcdf_err(error,
'OPENING FILE: '//trim(fngrid) )
769 error=nf90_inq_dimid(ncid,
'nx', id_dim)
770 CALL netcdf_err(error,
'ERROR READING NX ID' )
772 error=nf90_inquire_dimension(ncid,id_dim,len=nx)
775 error=nf90_inq_dimid(ncid,
'ny', id_dim)
776 CALL netcdf_err(error,
'ERROR READING NY ID' )
778 error=nf90_inquire_dimension(ncid,id_dim,len=ny)
781 IF ((nx/2) /= idim .OR. (ny/2) /= jdim)
THEN 782 print*,
'FATAL ERROR: DIMENSIONS IN FILE: ',(nx/2),(ny/2)
783 print*,
'DO NOT MATCH GRID DIMENSIONS: ',idim,jdim
784 CALL mpi_abort(mpi_comm_world, 130, error)
787 ALLOCATE(geolon(nx+1,ny+1))
788 ALLOCATE(geolat(nx+1,ny+1))
790 error=nf90_inq_varid(ncid,
'x', id_var)
792 error=nf90_get_var(ncid, id_var, geolon)
793 CALL netcdf_err(error,
'ERROR READING X RECORD' )
795 error=nf90_inq_varid(ncid,
'y', id_var)
797 error=nf90_get_var(ncid, id_var, geolat)
798 CALL netcdf_err(error,
'ERROR READING Y RECORD' )
800 ALLOCATE(dummy(idim,jdim))
806 dummy(i,j) = geolon(ii,jj)
810 rlo = reshape(dummy, (/ijdim/))
818 dummy(i,j) = geolat(ii,jj)
822 rla = reshape(dummy, (/ijdim/))
824 DEALLOCATE(geolat, dummy)
826 error=nf90_inq_varid(ncid,
'tile', id_var)
827 CALL netcdf_err(error,
'ERROR READING TILE ID' )
828 error=nf90_get_var(ncid, id_var, tile_num)
829 CALL netcdf_err(error,
'ERROR READING TILE RECORD' )
831 error = nf90_close(ncid)
833 fnorog =
"./fnorog." // rankch
836 print*,
"READ FV3 OROG INFO FROM: "//trim(fnorog)
838 error=nf90_open(trim(fnorog),nf90_nowrite,ncid_orog)
839 CALL netcdf_err(error,
'OPENING FILE: '//trim(fnorog) )
841 ALLOCATE(dummy4(idim,jdim))
843 error=nf90_inq_varid(ncid_orog,
'orog_raw', id_var)
844 CALL netcdf_err(error,
'ERROR READING orog_raw ID' )
845 error=nf90_get_var(ncid_orog, id_var, dummy4)
846 CALL netcdf_err(error,
'ERROR READING orog_raw RECORD' )
847 orog_uf = reshape(dummy4, (/ijdim/))
849 error=nf90_inq_varid(ncid_orog,
'orog_filt', id_var)
850 CALL netcdf_err(error,
'ERROR READING orog_filt ID' )
851 error=nf90_get_var(ncid_orog, id_var, dummy4)
852 CALL netcdf_err(error,
'ERROR READING orog_filt RECORD' )
853 orog = reshape(dummy4, (/ijdim/))
855 IF(
PRESENT(landfrac))
THEN 856 error=nf90_inq_varid(ncid_orog,
'land_frac', id_var)
857 CALL netcdf_err(error,
'ERROR READING land_frac ID' )
858 error=nf90_get_var(ncid_orog, id_var, dummy4)
859 CALL netcdf_err(error,
'ERROR READING land_frac RECORD' )
860 landfrac = reshape(dummy4, (/ijdim/))
865 error = nf90_close(ncid_orog)
881 INTEGER,
INTENT(IN) :: ERR
882 CHARACTER(LEN=*),
INTENT(IN) :: STRING
883 CHARACTER(LEN=80) :: ERRMSG
886 IF( err == nf90_noerr )
RETURN 887 errmsg = nf90_strerror(err)
889 print*,
'FATAL ERROR: ', trim(string),
': ', trim(errmsg)
891 CALL mpi_abort(mpi_comm_world, 999, iret)
912 CHARACTER(LEN=*),
INTENT(IN) :: gsi_file
913 CHARACTER(LEN=3),
INTENT(IN) :: file_type
914 INTEGER,
INTENT(IN),
OPTIONAL :: lsoil
916 INTEGER :: error, id_dim, ncid
919 INTEGER(KIND=1),
ALLOCATABLE :: idummy(:,:)
921 REAL(KIND=8),
ALLOCATABLE :: dummy(:,:)
923 CHARACTER(LEN=1) :: k_ch
924 CHARACTER(LEN=10) :: incvar
925 CHARACTER(LEN=80) :: err_msg
929 print*,
"READ INPUT GSI DATA FROM: "//trim(gsi_file)
931 error=nf90_open(trim(gsi_file),nf90_nowrite,ncid)
932 CALL netcdf_err(error,
'OPENING FILE: '//trim(gsi_file) )
934 error=nf90_inq_dimid(ncid,
'latitude', id_dim)
935 CALL netcdf_err(error,
'READING latitude ID' )
936 error=nf90_inquire_dimension(ncid,id_dim,len=
jdim_gaus)
937 CALL netcdf_err(error,
'READING latitude length' )
940 error=nf90_inq_dimid(ncid,
'longitude', id_dim)
941 CALL netcdf_err(error,
'READING longitude ID' )
942 error=nf90_inquire_dimension(ncid,id_dim,len=
idim_gaus)
943 CALL netcdf_err(error,
'READING longitude length' )
945 IF (file_type==
'NST')
then 949 error=nf90_inq_varid(ncid,
"dtf", id_var)
951 error=nf90_get_var(ncid, id_var, dummy)
957 error=nf90_inq_varid(ncid,
"msk", id_var)
959 error=nf90_get_var(ncid, id_var, idummy)
969 ELSEIF (file_type==
'LND')
then 977 WRITE(k_ch,
'(I1)') k
979 incvar =
"soilt"//k_ch//
"_inc" 980 error=nf90_inq_varid(ncid, incvar, id_var)
981 err_msg =
"reading "//incvar//
" ID" 983 error=nf90_get_var(ncid, id_var, dummy)
984 err_msg =
"reading "//incvar//
" data" 991 incvar =
"slc"//k_ch//
"_inc" 992 error=nf90_inq_varid(ncid, incvar, id_var)
993 err_msg =
"reading "//incvar//
" ID" 995 error=nf90_get_var(ncid, id_var, dummy)
996 err_msg =
"reading "//incvar//
" data" 1008 error=nf90_inq_varid(ncid,
"soilsnow_mask", id_var)
1009 CALL netcdf_err(error,
'READING soilsnow_mask ID' )
1010 error=nf90_get_var(ncid, id_var, idummy)
1011 CALL netcdf_err(error,
'READING soilsnow_mask' )
1021 print *,
'WARNING: FILE_TYPE', file_type,
'not recognised.', &
1022 ', no increments read in' 1025 IF(
ALLOCATED(dummy))
DEALLOCATE(dummy)
1026 IF(
ALLOCATED(idummy))
DEALLOCATE(idummy)
1028 error = nf90_close(ncid)
1087 SUBROUTINE read_data(LSOIL,LENSFC,DO_NSST,DO_SNO_INC_JEDI,&
1088 DO_SOI_INC_JEDI,INC_FILE,IS_NOAHMP, &
1089 TSFFCS,SMCFCS,SWEFCS,STCFCS, &
1091 CVFCS,CVBFCS,CVTFCS,ALBFCS, &
1092 VEGFCS,SLIFCS,CNPFCS,F10M, &
1093 VETFCS,SOTFCS,ALFFCS, &
1095 SIHFCS,SICFCS,SITFCS, &
1096 TPRCP,SRFLAG,SNDFCS, &
1097 VMNFCS,VMXFCS,SLCFCS, &
1099 SLPFCS,ABSFCS,T2M,Q2M,SLMASK, &
1105 INTEGER,
INTENT(IN) :: lsoil, lensfc
1106 LOGICAL,
INTENT(IN) :: do_nsst, inc_file
1107 LOGICAL,
INTENT(IN) :: do_sno_inc_jedi, do_soi_inc_jedi
1109 LOGICAL,
OPTIONAL,
INTENT(OUT) :: is_noahmp
1111 REAL,
OPTIONAL,
INTENT(OUT) :: cvfcs(lensfc), cvbfcs(lensfc)
1112 REAL,
OPTIONAL,
INTENT(OUT) :: cvtfcs(lensfc), albfcs(lensfc,4)
1113 REAL,
OPTIONAL,
INTENT(OUT) :: slifcs(lensfc), cnpfcs(lensfc)
1114 REAL,
OPTIONAL,
INTENT(OUT) :: vegfcs(lensfc), f10m(lensfc)
1115 REAL,
OPTIONAL,
INTENT(OUT) :: vetfcs(lensfc), sotfcs(lensfc)
1116 REAL,
OPTIONAL,
INTENT(OUT) :: tsffcs(lensfc), swefcs(lensfc)
1117 REAL,
OPTIONAL,
INTENT(OUT) :: tg3fcs(lensfc), zorfcs(lensfc)
1118 REAL,
OPTIONAL,
INTENT(OUT) :: alffcs(lensfc,2), ustar(lensfc)
1119 REAL,
OPTIONAL,
INTENT(OUT) :: fmm(lensfc), fhh(lensfc)
1120 REAL,
OPTIONAL,
INTENT(OUT) :: sihfcs(lensfc), sicfcs(lensfc)
1121 REAL,
OPTIONAL,
INTENT(OUT) :: sitfcs(lensfc), tprcp(lensfc)
1122 REAL,
OPTIONAL,
INTENT(OUT) :: srflag(lensfc), sndfcs(lensfc)
1123 REAL,
OPTIONAL,
INTENT(OUT) :: vmnfcs(lensfc), vmxfcs(lensfc)
1124 REAL,
OPTIONAL,
INTENT(OUT) :: slpfcs(lensfc), absfcs(lensfc)
1125 REAL,
OPTIONAL,
INTENT(OUT) :: t2m(lensfc), q2m(lensfc), slmask(lensfc)
1126 REAL,
OPTIONAL,
INTENT(OUT) :: slcfcs(lensfc,lsoil)
1127 REAL,
OPTIONAL,
INTENT(OUT) :: smcfcs(lensfc,lsoil)
1128 REAL,
OPTIONAL,
INTENT(OUT) :: stcfcs(lensfc,lsoil)
1129 REAL,
OPTIONAL,
INTENT(OUT) :: stcinc(lensfc,lsoil)
1130 REAL,
OPTIONAL,
INTENT(OUT) :: slcinc(lensfc,lsoil)
1131 REAL(KIND=4),
OPTIONAL,
INTENT(OUT) :: zsoil(lsoil)
1136 CHARACTER(LEN=50) :: fnbgsi
1137 CHARACTER(LEN=3) :: rankch
1139 INTEGER :: error, error2, ncid, myrank
1140 INTEGER :: idim, jdim, id_dim
1141 INTEGER :: id_var, ierr
1143 REAL(KIND=8),
ALLOCATABLE :: dummy(:,:), dummy3d(:,:,:)
1145 CALL mpi_comm_rank(mpi_comm_world, myrank, error)
1147 WRITE(rankch,
'(I3.3)') (myrank+1)
1149 IF ((inc_file) .and. (do_sno_inc_jedi))
THEN 1150 fnbgsi =
"./snow_xainc." // rankch
1151 ELSEIF ((inc_file) .and. (do_soi_inc_jedi))
THEN 1152 fnbgsi =
"./soil_xainc." // rankch
1154 fnbgsi =
"./fnbgsi." // rankch
1158 print*,
"READ INPUT SFC DATA FROM: "//trim(fnbgsi)
1160 error=nf90_open(trim(fnbgsi),nf90_nowrite,ncid)
1161 CALL netcdf_err(error,
'OPENING FILE: '//trim(fnbgsi) )
1163 IF ((inc_file) .and. (do_soi_inc_jedi))
THEN 1165 error=nf90_inq_dimid(ncid,
'grid_xt', id_dim)
1167 error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1170 error=nf90_inq_dimid(ncid,
'grid_yt', id_dim)
1172 error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1177 error=nf90_inq_dimid(ncid,
'xaxis_1', id_dim)
1179 error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1182 error=nf90_inq_dimid(ncid,
'yaxis_1', id_dim)
1184 error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1189 IF ((idim*jdim) /= lensfc)
THEN 1190 print*,
'FATAL ERROR: DIMENSIONS WRONG.' 1191 CALL mpi_abort(mpi_comm_world, 88, ierr)
1197 IF(
PRESENT(is_noahmp))
THEN 1198 error=nf90_inq_varid(ncid,
"canliqxy", id_var)
1199 error2=nf90_inq_varid(ncid,
"tsnoxy", id_var)
1201 IF(error == 0 .AND. error2 == 0)
THEN 1203 print*,
"- WILL PROCESS FOR NOAH-MP LSM." 1207 ALLOCATE(dummy(idim,jdim))
1209 IF (
PRESENT(tsffcs))
THEN 1210 error=nf90_inq_varid(ncid,
"tsea", id_var)
1212 error=nf90_get_var(ncid, id_var, dummy)
1214 tsffcs = reshape(dummy, (/lensfc/))
1217 IF (
PRESENT(swefcs))
THEN 1218 error=nf90_inq_varid(ncid,
"sheleg", id_var)
1220 error=nf90_get_var(ncid, id_var, dummy)
1222 swefcs = reshape(dummy, (/lensfc/))
1225 IF (
PRESENT(tg3fcs))
THEN 1226 error=nf90_inq_varid(ncid,
"tg3", id_var)
1228 error=nf90_get_var(ncid, id_var, dummy)
1230 tg3fcs = reshape(dummy, (/lensfc/))
1233 IF (
PRESENT(zorfcs))
THEN 1234 error=nf90_inq_varid(ncid,
"zorl", id_var)
1236 error=nf90_get_var(ncid, id_var, dummy)
1238 zorfcs = reshape(dummy, (/lensfc/))
1241 IF (
PRESENT(albfcs))
THEN 1243 error=nf90_inq_varid(ncid,
"alvsf", id_var)
1245 error=nf90_get_var(ncid, id_var, dummy)
1247 albfcs(:,1) = reshape(dummy, (/lensfc/))
1249 error=nf90_inq_varid(ncid,
"alvwf", id_var)
1251 error=nf90_get_var(ncid, id_var, dummy)
1253 albfcs(:,2) = reshape(dummy, (/lensfc/))
1255 error=nf90_inq_varid(ncid,
"alnsf", id_var)
1257 error=nf90_get_var(ncid, id_var, dummy)
1259 albfcs(:,3) = reshape(dummy, (/lensfc/))
1261 error=nf90_inq_varid(ncid,
"alnwf", id_var)
1263 error=nf90_get_var(ncid, id_var, dummy)
1265 albfcs(:,4) = reshape(dummy, (/lensfc/))
1269 IF (
PRESENT(slifcs))
THEN 1270 error=nf90_inq_varid(ncid,
"slmsk", id_var)
1272 error=nf90_get_var(ncid, id_var, dummy)
1274 slifcs = reshape(dummy, (/lensfc/))
1276 WHERE (slmask > 1.5) slmask=0.0
1279 IF (
PRESENT(cnpfcs))
THEN 1280 error=nf90_inq_varid(ncid,
"canopy", id_var)
1282 error=nf90_get_var(ncid, id_var, dummy)
1284 cnpfcs = reshape(dummy, (/lensfc/))
1287 IF (
PRESENT(vegfcs))
THEN 1288 error=nf90_inq_varid(ncid,
"vfrac", id_var)
1290 error=nf90_get_var(ncid, id_var, dummy)
1292 vegfcs = reshape(dummy, (/lensfc/))
1295 IF (
PRESENT(f10m))
THEN 1296 error=nf90_inq_varid(ncid,
"f10m", id_var)
1298 error=nf90_get_var(ncid, id_var, dummy)
1300 f10m = reshape(dummy, (/lensfc/))
1303 IF (
PRESENT(vetfcs))
THEN 1304 error=nf90_inq_varid(ncid,
"vtype", id_var)
1306 error=nf90_get_var(ncid, id_var, dummy)
1308 vetfcs = reshape(dummy, (/lensfc/))
1311 IF (
PRESENT(sotfcs))
THEN 1312 error=nf90_inq_varid(ncid,
"stype", id_var)
1314 error=nf90_get_var(ncid, id_var, dummy)
1316 sotfcs = reshape(dummy, (/lensfc/))
1319 IF (
PRESENT(alffcs))
THEN 1320 error=nf90_inq_varid(ncid,
"facsf", id_var)
1322 error=nf90_get_var(ncid, id_var, dummy)
1324 alffcs(:,1) = reshape(dummy, (/lensfc/))
1326 error=nf90_inq_varid(ncid,
"facwf", id_var)
1328 error=nf90_get_var(ncid, id_var, dummy)
1330 alffcs(:,2) = reshape(dummy, (/lensfc/))
1333 IF (
PRESENT(ustar))
THEN 1334 error=nf90_inq_varid(ncid,
"uustar", id_var)
1336 error=nf90_get_var(ncid, id_var, dummy)
1338 ustar = reshape(dummy, (/lensfc/))
1341 IF (
PRESENT(fmm))
THEN 1342 error=nf90_inq_varid(ncid,
"ffmm", id_var)
1344 error=nf90_get_var(ncid, id_var, dummy)
1346 fmm = reshape(dummy, (/lensfc/))
1349 IF (
PRESENT(fhh))
THEN 1350 error=nf90_inq_varid(ncid,
"ffhh", id_var)
1352 error=nf90_get_var(ncid, id_var, dummy)
1354 fhh = reshape(dummy, (/lensfc/))
1357 IF (
PRESENT(sihfcs))
THEN 1358 error=nf90_inq_varid(ncid,
"hice", id_var)
1360 error=nf90_get_var(ncid, id_var, dummy)
1362 sihfcs = reshape(dummy, (/lensfc/))
1365 IF (
PRESENT(sicfcs))
THEN 1366 error=nf90_inq_varid(ncid,
"fice", id_var)
1368 error=nf90_get_var(ncid, id_var, dummy)
1370 sicfcs = reshape(dummy, (/lensfc/))
1373 IF (
PRESENT(sitfcs))
THEN 1374 error=nf90_inq_varid(ncid,
"tisfc", id_var)
1376 error=nf90_get_var(ncid, id_var, dummy)
1378 sitfcs = reshape(dummy, (/lensfc/))
1381 IF (
PRESENT(tprcp))
THEN 1382 error=nf90_inq_varid(ncid,
"tprcp", id_var)
1384 error=nf90_get_var(ncid, id_var, dummy)
1386 tprcp = reshape(dummy, (/lensfc/))
1389 IF (
PRESENT(srflag))
THEN 1390 error=nf90_inq_varid(ncid,
"srflag", id_var)
1392 error=nf90_get_var(ncid, id_var, dummy)
1394 srflag = reshape(dummy, (/lensfc/))
1397 IF (
PRESENT(sndfcs))
THEN 1398 error=nf90_inq_varid(ncid,
"snwdph", id_var)
1400 error=nf90_get_var(ncid, id_var, dummy)
1402 sndfcs = reshape(dummy, (/lensfc/))
1405 IF (
PRESENT(vmnfcs))
THEN 1406 error=nf90_inq_varid(ncid,
"shdmin", id_var)
1408 error=nf90_get_var(ncid, id_var, dummy)
1410 vmnfcs = reshape(dummy, (/lensfc/))
1413 IF (
PRESENT(vmxfcs))
THEN 1414 error=nf90_inq_varid(ncid,
"shdmax", id_var)
1416 error=nf90_get_var(ncid, id_var, dummy)
1418 vmxfcs = reshape(dummy, (/lensfc/))
1421 IF (
PRESENT(slpfcs))
THEN 1422 error=nf90_inq_varid(ncid,
"slope", id_var)
1424 error=nf90_get_var(ncid, id_var, dummy)
1426 slpfcs = reshape(dummy, (/lensfc/))
1429 IF (
PRESENT(absfcs))
THEN 1430 error=nf90_inq_varid(ncid,
"snoalb", id_var)
1432 error=nf90_get_var(ncid, id_var, dummy)
1434 absfcs = reshape(dummy, (/lensfc/))
1437 IF (
PRESENT(t2m))
THEN 1438 error=nf90_inq_varid(ncid,
"t2m", id_var)
1440 error=nf90_get_var(ncid, id_var, dummy)
1442 t2m = reshape(dummy, (/lensfc/))
1445 IF (
PRESENT(q2m))
THEN 1446 error=nf90_inq_varid(ncid,
"q2m", id_var)
1448 error=nf90_get_var(ncid, id_var, dummy)
1450 q2m = reshape(dummy, (/lensfc/))
1453 nsst_read :
IF(do_nsst)
THEN 1456 print*,
"WILL READ NSST RECORDS." 1458 error=nf90_inq_varid(ncid,
"c_0", id_var)
1460 error=nf90_get_var(ncid, id_var, dummy)
1462 nsst%C_0 = reshape(dummy, (/lensfc/))
1464 error=nf90_inq_varid(ncid,
"c_d", id_var)
1466 error=nf90_get_var(ncid, id_var, dummy)
1468 nsst%C_D = reshape(dummy, (/lensfc/))
1470 error=nf90_inq_varid(ncid,
"d_conv", id_var)
1472 error=nf90_get_var(ncid, id_var, dummy)
1474 nsst%D_CONV = reshape(dummy, (/lensfc/))
1476 error=nf90_inq_varid(ncid,
"dt_cool", id_var)
1477 CALL netcdf_err(error,
'READING dt_cool ID' )
1478 error=nf90_get_var(ncid, id_var, dummy)
1480 nsst%DT_COOL = reshape(dummy, (/lensfc/))
1482 error=nf90_inq_varid(ncid,
"ifd", id_var)
1484 error=nf90_get_var(ncid, id_var, dummy)
1486 nsst%IFD = reshape(dummy, (/lensfc/))
1488 error=nf90_inq_varid(ncid,
"qrain", id_var)
1490 error=nf90_get_var(ncid, id_var, dummy)
1492 nsst%QRAIN = reshape(dummy, (/lensfc/))
1494 error=nf90_inq_varid(ncid,
"tref", id_var)
1496 error=nf90_get_var(ncid, id_var, dummy)
1498 nsst%TREF = reshape(dummy, (/lensfc/))
1500 error=nf90_inq_varid(ncid,
"w_0", id_var)
1502 error=nf90_get_var(ncid, id_var, dummy)
1504 nsst%W_0 = reshape(dummy, (/lensfc/))
1506 error=nf90_inq_varid(ncid,
"w_d", id_var)
1508 error=nf90_get_var(ncid, id_var, dummy)
1510 nsst%W_D = reshape(dummy, (/lensfc/))
1512 error=nf90_inq_varid(ncid,
"xs", id_var)
1514 error=nf90_get_var(ncid, id_var, dummy)
1516 nsst%XS = reshape(dummy, (/lensfc/))
1518 error=nf90_inq_varid(ncid,
"xt", id_var)
1520 error=nf90_get_var(ncid, id_var, dummy)
1522 nsst%XT = reshape(dummy, (/lensfc/))
1524 error=nf90_inq_varid(ncid,
"xtts", id_var)
1526 error=nf90_get_var(ncid, id_var, dummy)
1528 nsst%XTTS = reshape(dummy, (/lensfc/))
1530 error=nf90_inq_varid(ncid,
"xu", id_var)
1532 error=nf90_get_var(ncid, id_var, dummy)
1534 nsst%XU = reshape(dummy, (/lensfc/))
1536 error=nf90_inq_varid(ncid,
"xv", id_var)
1538 error=nf90_get_var(ncid, id_var, dummy)
1540 nsst%XV = reshape(dummy, (/lensfc/))
1542 error=nf90_inq_varid(ncid,
"xz", id_var)
1544 error=nf90_get_var(ncid, id_var, dummy)
1546 nsst%XZ = reshape(dummy, (/lensfc/))
1548 error=nf90_inq_varid(ncid,
"xzts", id_var)
1550 error=nf90_get_var(ncid, id_var, dummy)
1552 nsst%XZTS = reshape(dummy, (/lensfc/))
1554 error=nf90_inq_varid(ncid,
"z_c", id_var)
1556 error=nf90_get_var(ncid, id_var, dummy)
1558 nsst%Z_C = reshape(dummy, (/lensfc/))
1560 error=nf90_inq_varid(ncid,
"zm", id_var)
1562 error=nf90_get_var(ncid, id_var, dummy)
1564 nsst%ZM = reshape(dummy, (/lensfc/))
1570 ALLOCATE(dummy3d(idim,jdim,lsoil))
1572 IF (
PRESENT(smcfcs))
THEN 1573 error=nf90_inq_varid(ncid,
"smc", id_var)
1575 error=nf90_get_var(ncid, id_var, dummy3d)
1577 smcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1580 IF (
PRESENT(slcfcs))
THEN 1581 error=nf90_inq_varid(ncid,
"slc", id_var)
1583 error=nf90_get_var(ncid, id_var, dummy3d)
1585 slcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1588 IF (
PRESENT(stcfcs))
THEN 1589 error=nf90_inq_varid(ncid,
"stc", id_var)
1591 error=nf90_get_var(ncid, id_var, dummy3d)
1593 stcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1596 IF (
PRESENT(slcinc))
THEN 1597 error=nf90_inq_varid(ncid,
"soill", id_var)
1599 error=nf90_get_var(ncid, id_var, dummy3d)
1600 CALL netcdf_err(error,
'READING slc increments' )
1601 slcinc = reshape(dummy3d, (/lensfc,lsoil/))
1604 IF (
PRESENT(stcinc))
THEN 1605 error=nf90_inq_varid(ncid,
"soilt", id_var)
1607 error=nf90_get_var(ncid, id_var, dummy3d)
1608 CALL netcdf_err(error,
'READING stc increments' )
1609 stcinc = reshape(dummy3d, (/lensfc,lsoil/))
1616 IF (
PRESENT(cvfcs)) cvfcs = 0.0
1617 IF (
PRESENT(cvtfcs)) cvtfcs = 0.0
1618 IF (
PRESENT(cvbfcs)) cvbfcs = 0.0
1623 IF (
PRESENT(zsoil))
THEN 1630 error = nf90_close(ncid)
1650 subroutine read_tf_clim_grb(file_sst,sst,rlats_sst,rlons_sst,mlat_sst,mlon_sst,mon)
1657 character(*) ,
intent(in ) :: file_sst
1658 integer ,
intent(in ) :: mon,mlat_sst,mlon_sst
1659 real,
dimension(mlat_sst) ,
intent( out) :: rlats_sst
1660 real,
dimension(mlon_sst) ,
intent( out) :: rlons_sst
1661 real,
dimension(mlon_sst,mlat_sst) ,
intent( out) :: sst
1664 integer,
parameter:: lu_sst = 21
1665 real,
parameter :: deg2rad = 3.141593/180.0
1668 logical(1),
allocatable,
dimension(:) :: lb
1672 integer :: iret,ni,nj
1673 integer :: mscan,kb1,ierr
1674 integer :: jincdir,i,iincdir,kb2,kb3,kf,kg,k,j,jf
1675 integer,
dimension(22):: jgds,kgds
1676 integer,
dimension(25):: jpds,kpds
1681 real,
allocatable,
dimension(:) :: f
1686 write(*,*)
' sstclm : ',file_sst
1687 call baopenr(lu_sst,trim(file_sst),iret)
1688 if (iret /= 0 )
then 1689 write(6,*)
'FATAL ERROR in read_tf_clm_grb: error opening sst file.' 1690 CALL mpi_abort(mpi_comm_world, 111, ierr)
1700 call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1709 allocate(lb(nlat_sst*nlon_sst))
1710 allocate(f(nlat_sst*nlon_sst))
1711 jf=nlat_sst*nlon_sst
1714 call getgb(lu_sst,0,jf,j,jpds,jgds,kf,k,kpds,kgds,lb,f,iret)
1716 write(6,*)
'FATAL ERROR in read_tf_clm_grb: error reading sst analysis data record.' 1718 CALL mpi_abort(mpi_comm_world, 111, ierr)
1721 if ( (nlat_sst /= mlat_sst) .or. (nlon_sst /= mlon_sst) )
then 1722 write(6,*)
'FATAL ERROR in read_rtg_org: inconsistent dimensions. mlat_sst,mlon_sst=',&
1723 mlat_sst,mlon_sst,
' -versus- nlat_sst,nlon_sst=',nlat_sst,nlon_sst
1725 CALL mpi_abort(mpi_comm_world, 111, ierr)
1731 dres = 180.0/
real(nlat_sst)
1732 ysst0 = 0.5*dres-90.0
1737 rlats_sst(j) = ysst0 +
real(j-1)*dres
1741 rlons_sst(i) = (xsst0 +
real(i-1)*dres)
1749 kb1=ibits(mscan,7,1)
1750 kb2=ibits(mscan,6,1)
1751 kb3=ibits(mscan,5,1)
1766 i=(ni+1)*kb1+(mod(k-1,ni)+1)*iincdir
1767 j=(nj+1)*(1-kb2)+jincdir*((k-1)/ni+1)
1769 j=(nj+1)*(1-kb2)+(mod(k-1,nj)+1)*jincdir
1770 i=(ni+1)*kb1+iincdir*((k-1)/nj+1)
1777 call baclose(lu_sst,iret)
1778 if (iret /= 0 )
then 1779 write(6,*)
'FATAL ERROR in read_tf_clm_grb: error closing sst file.' 1780 CALL mpi_abort(mpi_comm_world, 121, ierr)
1798 character(*) ,
intent(in ) :: file_sst
1799 integer ,
intent(out) :: mlat_sst,mlon_sst
1802 integer,
parameter:: lu_sst = 21
1805 integer :: kf,kg,k,j,ierr
1806 integer,
dimension(22):: jgds,kgds
1807 integer,
dimension(25):: jpds,kpds
1812 call baopenr(lu_sst,trim(file_sst),iret)
1813 if (iret /= 0 )
then 1814 write(6,*)
'FATAL ERROR in get_tf_clm_dim: error opening sst file.' 1815 CALL mpi_abort(mpi_comm_world, 111, ierr)
1825 call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1830 write(*,*)
'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1832 call baclose(lu_sst,iret)
1833 if (iret /= 0 )
then 1834 write(6,*)
'FATAL ERROR in get_tf_clm_dim: error closing sst file.' 1835 CALL mpi_abort(mpi_comm_world, 121, ierr)
1855 character (len=*),
intent(in) :: filename
1856 integer,
intent(in) :: nlat,nlon
1857 integer,
intent(in) :: itime
1858 real,
dimension(nlat),
intent(out) :: xlats
1859 real,
dimension(nlon),
intent(out) :: xlons
1860 real,
dimension(nlon,nlat),
intent(out) :: sal
1864 integer,
parameter :: ndims = 3
1865 character (len = *),
parameter :: lat_name =
"latitude" 1866 character (len = *),
parameter :: lon_name =
"longitude" 1867 character (len = *),
parameter :: t_name =
"time" 1868 character (len = *),
parameter :: sal_name=
"sal" 1869 integer :: time_varid,lon_varid, lat_varid, sal_varid
1872 integer,
dimension(ndims) :: start, count
1874 character (len = *),
parameter :: units =
"units" 1875 character (len = *),
parameter :: sal_units =
"psu" 1877 character (len = *),
parameter :: time_units =
"months" 1878 character (len = *),
parameter :: lat_units =
"degrees_north" 1879 character (len = *),
parameter :: lon_units =
"degrees_east" 1882 call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1885 call nc_check( nf90_inq_varid(ncid, t_name, time_varid) )
1886 call nc_check( nf90_inq_varid(ncid, lat_name, lat_varid) )
1887 call nc_check( nf90_inq_varid(ncid, lon_name, lon_varid) )
1891 call nc_check( nf90_get_var(ncid, lat_varid, xlats) )
1892 call nc_check( nf90_get_var(ncid, lon_varid, xlons) )
1895 call nc_check( nf90_inq_varid(ncid, sal_name,sal_varid) )
1899 start = (/ 1, 1, itime /)
1900 count = (/ nlon, nlat, 1 /)
1904 call nc_check( nf90_get_var(ncid, sal_varid, sal, start, count) )
1925 character (len=*),
intent(in) :: filename
1926 integer,
intent(out) :: nlat,nlon
1928 character (len = *),
parameter :: lat_name =
"latitude" 1929 character (len = *),
parameter :: lon_name =
"longitude" 1931 integer :: latdimid,londimid
1934 call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1937 call nc_check( nf90_inq_dimid(ncid,lat_name,latdimid) )
1938 call nc_check( nf90_inq_dimid(ncid,lon_name,londimid) )
1939 call nc_check( nf90_inquire_dimension(ncid,latdimid,len=nlat) )
1940 call nc_check( nf90_inquire_dimension(ncid,londimid,len=nlon) )
1963 integer,
intent ( in) :: status
1966 if(status /= nf90_noerr)
then 1967 print *,
'FATAL ERROR:' 1968 print *, trim(nf90_strerror(status))
1969 CALL mpi_abort(mpi_comm_world, 122, ierr)
subroutine, public read_data(LSOIL, LENSFC, DO_NSST, DO_SNO_INC_JEDI, DO_SOI_INC_JEDI, INC_FILE, IS_NOAHMP, 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, SLPFCS, ABSFCS, T2M, Q2M, SLMASK, ZSOIL, NSST)
Read the first guess surface records and nsst records (if selected) for a single cubed-sphere tile...
subroutine, public get_tf_clm_dim(file_sst, mlat_sst, mlon_sst)
Get the i/j dimensions of RTG SST climatology file.
integer, dimension(:,:), allocatable, public soilsnow_gaus
GSI soil / snow mask for land 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.
integer, public idim_gaus
'i' dimension of GSI gaussian grid.
subroutine, public read_salclm_gfs_nc(filename, sal, xlats, xlons, nlat, nlon, itime)
Read the woa05 salinity monthly climatology file.
integer, dimension(:,:), allocatable, public slmsk_gaus
GSI land mask on the gaussian grid.
Holds machine dependent constants for global_cycle.
subroutine, public get_dim_nc(filename, nlat, nlon)
Get the i/j dimensions of the data from a NetCDF file.
real, dimension(:,:,:), allocatable, public slc_inc_gaus
GSI soil moisture increments on the gaussian grid.
subroutine, public read_gsi_data(GSI_FILE, FILE_TYPE, LSOIL)
Read increment file from the GSI containing either the foundation temperature increments and mask...
real, dimension(:,:,:), allocatable, public stc_inc_gaus
GSI soil temperature increments on the gaussian grid.
subroutine remove_checksum(ncid, id_var)
Remove the checksum attribute from a netcdf record.
real, dimension(:,:), allocatable, public dtref_gaus
GSI foundation temperature increment on the gaussian grid.
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...
This module contains routines that read and write data.
integer, public jdim_gaus
'j' dimension of GSI gaussian grid.
subroutine nc_check(status)
Check the NetCDF status code.
subroutine netcdf_err(ERR, STRING)
If a NetCDF call returns an error, print out a user-supplied message and the NetCDF library message...
subroutine, public read_lat_lon_orog(RLA, RLO, OROG, OROG_UF, TILE_NUM, IDIM, JDIM, IJDIM, LANDFRAC)
Read latitude and longitude for the cubed-sphere tile from the 'grid' file.