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(:)
38 INTEGER,
PUBLIC :: idim_gaus
40 INTEGER,
PUBLIC :: jdim_gaus
42 INTEGER,
ALLOCATABLE,
PUBLIC :: slmsk_gaus(:,:)
45 INTEGER,
ALLOCATABLE,
PUBLIC :: soilsnow_gaus(:,:)
49 REAL,
ALLOCATABLE,
PUBLIC :: dtref_gaus(:,:)
52 REAL,
ALLOCATABLE,
PUBLIC :: stc_inc_gaus(:,:,:)
55 REAL,
ALLOCATABLE,
PUBLIC :: slc_inc_gaus(:,:,:)
59 PUBLIC :: read_gsi_data
60 PUBLIC :: read_lat_lon_orog
62 public :: read_tf_clim_grb,get_tf_clm_dim
63 public :: read_salclm_gfs_nc,get_dim_nc
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)
161 type(nsst_data),
intent(in) :: nsst
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
496 call netcdf_err(error,
'defining xaxis_1')
498 error = nf90_def_dim(ncid,
"yaxis_1", jdim, dim_y)
499 call netcdf_err(error,
'defining yaxis_1')
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)
507 call netcdf_err(error,
'defining slc_inc');
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');
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
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)
656 error=nf90_inq_varid(ncid,
"tfinc", id_var)
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' )
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)
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)
684 end subroutine write_data
692 subroutine remove_checksum(ncid, id_var)
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')
708 call netcdf_err(error,
'deleting checksum' )
710 error= nf90_enddef(ncid)
711 call netcdf_err(error,
'ending define mode' )
715 end subroutine remove_checksum
731 SUBROUTINE read_lat_lon_orog(RLA,RLO,OROG,OROG_UF,&
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)
773 CALL netcdf_err(error,
'ERROR READING 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)
779 CALL netcdf_err(error,
'ERROR READING 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)
791 CALL netcdf_err(error,
'ERROR READING X ID' )
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)
796 CALL netcdf_err(error,
'ERROR READING Y ID' )
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)
867 END SUBROUTINE read_lat_lon_orog
875 SUBROUTINE netcdf_err( ERR, STRING )
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)
894 END SUBROUTINE netcdf_err
908 SUBROUTINE read_gsi_data(GSI_FILE, FILE_TYPE, LSOIL)
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' )
938 jdim_gaus = jdim_gaus - 2
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 946 ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
947 ALLOCATE(dtref_gaus(idim_gaus,jdim_gaus))
949 error=nf90_inq_varid(ncid,
"dtf", id_var)
950 CALL netcdf_err(error,
'READING dtf ID' )
951 error=nf90_get_var(ncid, id_var, dummy)
952 CALL netcdf_err(error,
'READING dtf' )
954 ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
955 ALLOCATE(slmsk_gaus(idim_gaus,jdim_gaus))
957 error=nf90_inq_varid(ncid,
"msk", id_var)
958 CALL netcdf_err(error,
'READING msk ID' )
959 error=nf90_get_var(ncid, id_var, idummy)
960 CALL netcdf_err(error,
'READING msk' )
965 slmsk_gaus(:,j) = idummy(:,j+1)
966 dtref_gaus(:,j) = dummy(:,j+1)
969 ELSEIF (file_type==
'LND')
then 971 ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
972 ALLOCATE(stc_inc_gaus(lsoil,idim_gaus,jdim_gaus))
973 ALLOCATE(slc_inc_gaus(lsoil,idim_gaus,jdim_gaus))
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" 982 CALL netcdf_err(error, trim(err_msg))
983 error=nf90_get_var(ncid, id_var, dummy)
984 err_msg =
"reading "//incvar//
" data" 985 CALL netcdf_err(error, err_msg)
988 stc_inc_gaus(k,:,j) = dummy(:,j+1)
991 incvar =
"slc"//k_ch//
"_inc" 992 error=nf90_inq_varid(ncid, incvar, id_var)
993 err_msg =
"reading "//incvar//
" ID" 994 CALL netcdf_err(error, trim(err_msg))
995 error=nf90_get_var(ncid, id_var, dummy)
996 err_msg =
"reading "//incvar//
" data" 997 CALL netcdf_err(error, err_msg)
1000 slc_inc_gaus(k,:,j) = dummy(:,j+1)
1005 ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
1006 ALLOCATE(soilsnow_gaus(idim_gaus,jdim_gaus))
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' )
1016 soilsnow_gaus(:,j) = idummy(:,j+1)
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)
1030 END SUBROUTINE read_gsi_data
1084 SUBROUTINE read_data(LSOIL,LENSFC,DO_NSST,IS_NOAHMP, &
1086 TSFFCS,SMCFCS,SWEFCS,STCFCS, &
1088 CVFCS,CVBFCS,CVTFCS,ALBFCS, &
1089 VEGFCS,SLIFCS,CNPFCS,F10M, &
1090 VETFCS,SOTFCS,ALFFCS, &
1092 SIHFCS,SICFCS,SITFCS, &
1093 TPRCP,SRFLAG,SNDFCS, &
1094 VMNFCS,VMXFCS,SLCFCS, &
1095 STCINC,SLCINC,LSOIL_INCR, &
1096 SLPFCS,ABSFCS,T2M,Q2M,SLMASK, &
1102 INTEGER,
INTENT(IN) :: lsoil, lensfc
1103 LOGICAL,
INTENT(IN) :: do_nsst
1105 CHARACTER(LEN=50),
OPTIONAL,
INTENT(IN) :: fname_inc
1106 INTEGER,
OPTIONAL,
INTENT(IN) :: lsoil_incr
1108 LOGICAL,
OPTIONAL,
INTENT(OUT) :: is_noahmp
1110 REAL,
OPTIONAL,
INTENT(OUT) :: cvfcs(lensfc), cvbfcs(lensfc)
1111 REAL,
OPTIONAL,
INTENT(OUT) :: cvtfcs(lensfc), albfcs(lensfc,4)
1112 REAL,
OPTIONAL,
INTENT(OUT) :: slifcs(lensfc), cnpfcs(lensfc)
1113 REAL,
OPTIONAL,
INTENT(OUT) :: vegfcs(lensfc), f10m(lensfc)
1114 REAL,
OPTIONAL,
INTENT(OUT) :: vetfcs(lensfc), sotfcs(lensfc)
1115 REAL,
OPTIONAL,
INTENT(OUT) :: tsffcs(lensfc), swefcs(lensfc)
1116 REAL,
OPTIONAL,
INTENT(OUT) :: tg3fcs(lensfc), zorfcs(lensfc)
1117 REAL,
OPTIONAL,
INTENT(OUT) :: alffcs(lensfc,2), ustar(lensfc)
1118 REAL,
OPTIONAL,
INTENT(OUT) :: fmm(lensfc), fhh(lensfc)
1119 REAL,
OPTIONAL,
INTENT(OUT) :: sihfcs(lensfc), sicfcs(lensfc)
1120 REAL,
OPTIONAL,
INTENT(OUT) :: sitfcs(lensfc), tprcp(lensfc)
1121 REAL,
OPTIONAL,
INTENT(OUT) :: srflag(lensfc), sndfcs(lensfc)
1122 REAL,
OPTIONAL,
INTENT(OUT) :: vmnfcs(lensfc), vmxfcs(lensfc)
1123 REAL,
OPTIONAL,
INTENT(OUT) :: slpfcs(lensfc), absfcs(lensfc)
1124 REAL,
OPTIONAL,
INTENT(OUT) :: t2m(lensfc), q2m(lensfc), slmask(lensfc)
1125 REAL,
OPTIONAL,
INTENT(OUT) :: slcfcs(lensfc,lsoil)
1126 REAL,
OPTIONAL,
INTENT(OUT) :: smcfcs(lensfc,lsoil)
1127 REAL,
OPTIONAL,
INTENT(OUT) :: stcfcs(lensfc,lsoil)
1128 REAL,
OPTIONAL,
INTENT(OUT) :: stcinc(lensfc,lsoil)
1129 REAL,
OPTIONAL,
INTENT(OUT) :: slcinc(lensfc,lsoil)
1130 REAL(KIND=4),
OPTIONAL,
INTENT(OUT) :: zsoil(lsoil)
1132 TYPE(nsst_data),
OPTIONAL :: nsst
1134 CHARACTER(LEN=3) :: rankch
1135 CHARACTER(LEN=50) :: fname
1136 CHARACTER(LEN=1) :: k_ch
1137 CHARACTER(LEN=10) :: incvar
1139 INTEGER :: error, error2, ncid, myrank
1140 INTEGER :: idim, jdim, id_dim
1141 INTEGER :: id_var, ierr, test, k
1143 LOGICAL :: jedi_incr_file
1145 REAL(KIND=8),
ALLOCATABLE :: dummy(:,:), dummy3d(:,:,:)
1147 IF (
PRESENT(fname_inc))
THEN 1150 CALL mpi_comm_rank(mpi_comm_world, myrank, error)
1152 WRITE(rankch,
'(I3.3)') (myrank+1)
1154 fname =
"./fnbgsi." // rankch
1158 print*,
"READ INPUT SFC DATA FROM: "//trim(fname)
1160 error=nf90_open(trim(fname),nf90_nowrite,ncid)
1161 CALL netcdf_err(error,
'OPENING FILE: '//trim(fname) )
1166 test=nf90_inq_dimid(ncid,
'xaxis_1', id_dim)
1168 IF ( test == nf90_noerr )
THEN 1169 jedi_incr_file=.false.
1171 error=nf90_inq_dimid(ncid,
'xaxis_1', id_dim)
1172 CALL netcdf_err(error,
'READING xaxis_1' )
1173 error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1174 CALL netcdf_err(error,
'READING xaxis_1' )
1176 error=nf90_inq_dimid(ncid,
'yaxis_1', id_dim)
1177 CALL netcdf_err(error,
'READING yaxis_1' )
1178 error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1179 CALL netcdf_err(error,
'READING yaxis_1' )
1182 jedi_incr_file=.true.
1184 error=nf90_inq_dimid(ncid,
'grid_xt', id_dim)
1185 CALL netcdf_err(error,
'READING grid_xt' )
1186 error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1187 CALL netcdf_err(error,
'READING grid_xt' )
1189 error=nf90_inq_dimid(ncid,
'grid_yt', id_dim)
1190 CALL netcdf_err(error,
'READING grid_yt' )
1191 error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1192 CALL netcdf_err(error,
'READING grid_yt' )
1196 IF ((idim*jdim) /= lensfc)
THEN 1197 print*,
'FATAL ERROR: DIMENSIONS WRONG.' 1198 CALL mpi_abort(mpi_comm_world, 88, ierr)
1204 IF(
PRESENT(is_noahmp))
THEN 1205 error=nf90_inq_varid(ncid,
"canliqxy", id_var)
1206 error2=nf90_inq_varid(ncid,
"tsnoxy", id_var)
1208 IF(error == 0 .AND. error2 == 0)
THEN 1210 print*,
"- WILL PROCESS FOR NOAH-MP LSM." 1214 ALLOCATE(dummy(idim,jdim))
1216 IF (
PRESENT(tsffcs))
THEN 1217 error=nf90_inq_varid(ncid,
"tsea", id_var)
1218 CALL netcdf_err(error,
'READING tsea ID' )
1219 error=nf90_get_var(ncid, id_var, dummy)
1220 CALL netcdf_err(error,
'READING tsea' )
1221 tsffcs = reshape(dummy, (/lensfc/))
1224 IF (
PRESENT(swefcs))
THEN 1225 error=nf90_inq_varid(ncid,
"sheleg", id_var)
1226 CALL netcdf_err(error,
'READING sheleg ID' )
1227 error=nf90_get_var(ncid, id_var, dummy)
1228 CALL netcdf_err(error,
'READING sheleg' )
1229 swefcs = reshape(dummy, (/lensfc/))
1232 IF (
PRESENT(tg3fcs))
THEN 1233 error=nf90_inq_varid(ncid,
"tg3", id_var)
1234 CALL netcdf_err(error,
'READING tg3 ID' )
1235 error=nf90_get_var(ncid, id_var, dummy)
1236 CALL netcdf_err(error,
'READING tg3' )
1237 tg3fcs = reshape(dummy, (/lensfc/))
1240 IF (
PRESENT(zorfcs))
THEN 1241 error=nf90_inq_varid(ncid,
"zorl", id_var)
1242 CALL netcdf_err(error,
'READING zorl ID' )
1243 error=nf90_get_var(ncid, id_var, dummy)
1244 CALL netcdf_err(error,
'READING zorl' )
1245 zorfcs = reshape(dummy, (/lensfc/))
1248 IF (
PRESENT(albfcs))
THEN 1250 error=nf90_inq_varid(ncid,
"alvsf", id_var)
1251 CALL netcdf_err(error,
'READING alvsf ID' )
1252 error=nf90_get_var(ncid, id_var, dummy)
1253 CALL netcdf_err(error,
'READING alvsf' )
1254 albfcs(:,1) = reshape(dummy, (/lensfc/))
1256 error=nf90_inq_varid(ncid,
"alvwf", id_var)
1257 CALL netcdf_err(error,
'READING alvwf ID' )
1258 error=nf90_get_var(ncid, id_var, dummy)
1259 CALL netcdf_err(error,
'READING alvwf' )
1260 albfcs(:,2) = reshape(dummy, (/lensfc/))
1262 error=nf90_inq_varid(ncid,
"alnsf", id_var)
1263 CALL netcdf_err(error,
'READING alnsf ID' )
1264 error=nf90_get_var(ncid, id_var, dummy)
1265 CALL netcdf_err(error,
'READING alnsf' )
1266 albfcs(:,3) = reshape(dummy, (/lensfc/))
1268 error=nf90_inq_varid(ncid,
"alnwf", id_var)
1269 CALL netcdf_err(error,
'READING alnwf ID' )
1270 error=nf90_get_var(ncid, id_var, dummy)
1271 CALL netcdf_err(error,
'READING alnwf' )
1272 albfcs(:,4) = reshape(dummy, (/lensfc/))
1276 IF (
PRESENT(slifcs))
THEN 1277 error=nf90_inq_varid(ncid,
"slmsk", id_var)
1278 CALL netcdf_err(error,
'READING slmsk ID' )
1279 error=nf90_get_var(ncid, id_var, dummy)
1280 CALL netcdf_err(error,
'READING slmsk' )
1281 slifcs = reshape(dummy, (/lensfc/))
1283 WHERE (slmask > 1.5) slmask=0.0
1286 IF (
PRESENT(cnpfcs))
THEN 1287 error=nf90_inq_varid(ncid,
"canopy", id_var)
1288 CALL netcdf_err(error,
'READING canopy ID' )
1289 error=nf90_get_var(ncid, id_var, dummy)
1290 CALL netcdf_err(error,
'READING canopy' )
1291 cnpfcs = reshape(dummy, (/lensfc/))
1294 IF (
PRESENT(vegfcs))
THEN 1295 error=nf90_inq_varid(ncid,
"vfrac", id_var)
1296 CALL netcdf_err(error,
'READING vfrac ID' )
1297 error=nf90_get_var(ncid, id_var, dummy)
1298 CALL netcdf_err(error,
'READING vfrac' )
1299 vegfcs = reshape(dummy, (/lensfc/))
1302 IF (
PRESENT(f10m))
THEN 1303 error=nf90_inq_varid(ncid,
"f10m", id_var)
1304 CALL netcdf_err(error,
'READING f10m ID' )
1305 error=nf90_get_var(ncid, id_var, dummy)
1306 CALL netcdf_err(error,
'READING f10m' )
1307 f10m = reshape(dummy, (/lensfc/))
1310 IF (
PRESENT(vetfcs))
THEN 1311 error=nf90_inq_varid(ncid,
"vtype", id_var)
1312 CALL netcdf_err(error,
'READING vtype ID' )
1313 error=nf90_get_var(ncid, id_var, dummy)
1314 CALL netcdf_err(error,
'READING vtype' )
1315 vetfcs = reshape(dummy, (/lensfc/))
1318 IF (
PRESENT(sotfcs))
THEN 1319 error=nf90_inq_varid(ncid,
"stype", id_var)
1320 CALL netcdf_err(error,
'READING stype ID' )
1321 error=nf90_get_var(ncid, id_var, dummy)
1322 CALL netcdf_err(error,
'READING stype' )
1323 sotfcs = reshape(dummy, (/lensfc/))
1326 IF (
PRESENT(alffcs))
THEN 1327 error=nf90_inq_varid(ncid,
"facsf", id_var)
1328 CALL netcdf_err(error,
'READING facsf ID' )
1329 error=nf90_get_var(ncid, id_var, dummy)
1330 CALL netcdf_err(error,
'READING facsf' )
1331 alffcs(:,1) = reshape(dummy, (/lensfc/))
1333 error=nf90_inq_varid(ncid,
"facwf", id_var)
1334 CALL netcdf_err(error,
'READING facwf ID' )
1335 error=nf90_get_var(ncid, id_var, dummy)
1336 CALL netcdf_err(error,
'READING facwf' )
1337 alffcs(:,2) = reshape(dummy, (/lensfc/))
1340 IF (
PRESENT(ustar))
THEN 1341 error=nf90_inq_varid(ncid,
"uustar", id_var)
1342 CALL netcdf_err(error,
'READING uustar ID' )
1343 error=nf90_get_var(ncid, id_var, dummy)
1344 CALL netcdf_err(error,
'READING uustar' )
1345 ustar = reshape(dummy, (/lensfc/))
1348 IF (
PRESENT(fmm))
THEN 1349 error=nf90_inq_varid(ncid,
"ffmm", id_var)
1350 CALL netcdf_err(error,
'READING ffmm ID' )
1351 error=nf90_get_var(ncid, id_var, dummy)
1352 CALL netcdf_err(error,
'READING ffmm' )
1353 fmm = reshape(dummy, (/lensfc/))
1356 IF (
PRESENT(fhh))
THEN 1357 error=nf90_inq_varid(ncid,
"ffhh", id_var)
1358 CALL netcdf_err(error,
'READING ffhh ID' )
1359 error=nf90_get_var(ncid, id_var, dummy)
1360 CALL netcdf_err(error,
'READING ffhh' )
1361 fhh = reshape(dummy, (/lensfc/))
1364 IF (
PRESENT(sihfcs))
THEN 1365 error=nf90_inq_varid(ncid,
"hice", id_var)
1366 CALL netcdf_err(error,
'READING hice ID' )
1367 error=nf90_get_var(ncid, id_var, dummy)
1368 CALL netcdf_err(error,
'READING hice' )
1369 sihfcs = reshape(dummy, (/lensfc/))
1372 IF (
PRESENT(sicfcs))
THEN 1373 error=nf90_inq_varid(ncid,
"fice", id_var)
1374 CALL netcdf_err(error,
'READING fice ID' )
1375 error=nf90_get_var(ncid, id_var, dummy)
1376 CALL netcdf_err(error,
'READING fice' )
1377 sicfcs = reshape(dummy, (/lensfc/))
1380 IF (
PRESENT(sitfcs))
THEN 1381 error=nf90_inq_varid(ncid,
"tisfc", id_var)
1382 CALL netcdf_err(error,
'READING tisfc ID' )
1383 error=nf90_get_var(ncid, id_var, dummy)
1384 CALL netcdf_err(error,
'READING tisfc' )
1385 sitfcs = reshape(dummy, (/lensfc/))
1388 IF (
PRESENT(tprcp))
THEN 1389 error=nf90_inq_varid(ncid,
"tprcp", id_var)
1390 CALL netcdf_err(error,
'READING tprcp ID' )
1391 error=nf90_get_var(ncid, id_var, dummy)
1392 CALL netcdf_err(error,
'READING tprcp' )
1393 tprcp = reshape(dummy, (/lensfc/))
1396 IF (
PRESENT(srflag))
THEN 1397 error=nf90_inq_varid(ncid,
"srflag", id_var)
1398 CALL netcdf_err(error,
'READING srflag ID' )
1399 error=nf90_get_var(ncid, id_var, dummy)
1400 CALL netcdf_err(error,
'READING srflag' )
1401 srflag = reshape(dummy, (/lensfc/))
1404 IF (
PRESENT(sndfcs))
THEN 1405 error=nf90_inq_varid(ncid,
"snwdph", id_var)
1406 CALL netcdf_err(error,
'READING snwdph ID' )
1407 error=nf90_get_var(ncid, id_var, dummy)
1408 CALL netcdf_err(error,
'READING snwdph' )
1409 sndfcs = reshape(dummy, (/lensfc/))
1412 IF (
PRESENT(vmnfcs))
THEN 1413 error=nf90_inq_varid(ncid,
"shdmin", id_var)
1414 CALL netcdf_err(error,
'READING shdmin ID' )
1415 error=nf90_get_var(ncid, id_var, dummy)
1416 CALL netcdf_err(error,
'READING shdmin' )
1417 vmnfcs = reshape(dummy, (/lensfc/))
1420 IF (
PRESENT(vmxfcs))
THEN 1421 error=nf90_inq_varid(ncid,
"shdmax", id_var)
1422 CALL netcdf_err(error,
'READING shdmax ID' )
1423 error=nf90_get_var(ncid, id_var, dummy)
1424 CALL netcdf_err(error,
'READING shdmax' )
1425 vmxfcs = reshape(dummy, (/lensfc/))
1428 IF (
PRESENT(slpfcs))
THEN 1429 error=nf90_inq_varid(ncid,
"slope", id_var)
1430 CALL netcdf_err(error,
'READING slope ID' )
1431 error=nf90_get_var(ncid, id_var, dummy)
1432 CALL netcdf_err(error,
'READING slope' )
1433 slpfcs = reshape(dummy, (/lensfc/))
1436 IF (
PRESENT(absfcs))
THEN 1437 error=nf90_inq_varid(ncid,
"snoalb", id_var)
1438 CALL netcdf_err(error,
'READING snoalb ID' )
1439 error=nf90_get_var(ncid, id_var, dummy)
1440 CALL netcdf_err(error,
'READING snoalb' )
1441 absfcs = reshape(dummy, (/lensfc/))
1444 IF (
PRESENT(t2m))
THEN 1445 error=nf90_inq_varid(ncid,
"t2m", id_var)
1446 CALL netcdf_err(error,
'READING t2m ID' )
1447 error=nf90_get_var(ncid, id_var, dummy)
1448 CALL netcdf_err(error,
'READING t2m' )
1449 t2m = reshape(dummy, (/lensfc/))
1452 IF (
PRESENT(q2m))
THEN 1453 error=nf90_inq_varid(ncid,
"q2m", id_var)
1454 CALL netcdf_err(error,
'READING q2m ID' )
1455 error=nf90_get_var(ncid, id_var, dummy)
1456 CALL netcdf_err(error,
'READING q2m' )
1457 q2m = reshape(dummy, (/lensfc/))
1460 nsst_read :
IF(do_nsst)
THEN 1463 print*,
"WILL READ NSST RECORDS." 1465 error=nf90_inq_varid(ncid,
"c_0", id_var)
1466 CALL netcdf_err(error,
'READING c_0 ID' )
1467 error=nf90_get_var(ncid, id_var, dummy)
1468 CALL netcdf_err(error,
'READING c_0' )
1469 nsst%C_0 = reshape(dummy, (/lensfc/))
1471 error=nf90_inq_varid(ncid,
"c_d", id_var)
1472 CALL netcdf_err(error,
'READING c_d ID' )
1473 error=nf90_get_var(ncid, id_var, dummy)
1474 CALL netcdf_err(error,
'READING c_d' )
1475 nsst%C_D = reshape(dummy, (/lensfc/))
1477 error=nf90_inq_varid(ncid,
"d_conv", id_var)
1478 CALL netcdf_err(error,
'READING d_conv ID' )
1479 error=nf90_get_var(ncid, id_var, dummy)
1480 CALL netcdf_err(error,
'READING d_conv' )
1481 nsst%D_CONV = reshape(dummy, (/lensfc/))
1483 error=nf90_inq_varid(ncid,
"dt_cool", id_var)
1484 CALL netcdf_err(error,
'READING dt_cool ID' )
1485 error=nf90_get_var(ncid, id_var, dummy)
1486 CALL netcdf_err(error,
'READING dt_cool' )
1487 nsst%DT_COOL = reshape(dummy, (/lensfc/))
1489 error=nf90_inq_varid(ncid,
"ifd", id_var)
1490 CALL netcdf_err(error,
'READING ifd ID' )
1491 error=nf90_get_var(ncid, id_var, dummy)
1492 CALL netcdf_err(error,
'READING ifd' )
1493 nsst%IFD = reshape(dummy, (/lensfc/))
1495 error=nf90_inq_varid(ncid,
"qrain", id_var)
1496 CALL netcdf_err(error,
'READING qrain ID' )
1497 error=nf90_get_var(ncid, id_var, dummy)
1498 CALL netcdf_err(error,
'READING qrain' )
1499 nsst%QRAIN = reshape(dummy, (/lensfc/))
1501 error=nf90_inq_varid(ncid,
"tref", id_var)
1502 CALL netcdf_err(error,
'READING tref ID' )
1503 error=nf90_get_var(ncid, id_var, dummy)
1504 CALL netcdf_err(error,
'READING tref' )
1505 nsst%TREF = reshape(dummy, (/lensfc/))
1507 error=nf90_inq_varid(ncid,
"w_0", id_var)
1508 CALL netcdf_err(error,
'READING w_0 ID' )
1509 error=nf90_get_var(ncid, id_var, dummy)
1510 CALL netcdf_err(error,
'READING w_0' )
1511 nsst%W_0 = reshape(dummy, (/lensfc/))
1513 error=nf90_inq_varid(ncid,
"w_d", id_var)
1514 CALL netcdf_err(error,
'READING w_d ID' )
1515 error=nf90_get_var(ncid, id_var, dummy)
1516 CALL netcdf_err(error,
'READING w_d' )
1517 nsst%W_D = reshape(dummy, (/lensfc/))
1519 error=nf90_inq_varid(ncid,
"xs", id_var)
1520 CALL netcdf_err(error,
'READING xs ID' )
1521 error=nf90_get_var(ncid, id_var, dummy)
1522 CALL netcdf_err(error,
'READING xs' )
1523 nsst%XS = reshape(dummy, (/lensfc/))
1525 error=nf90_inq_varid(ncid,
"xt", id_var)
1526 CALL netcdf_err(error,
'READING xt ID' )
1527 error=nf90_get_var(ncid, id_var, dummy)
1528 CALL netcdf_err(error,
'READING xt' )
1529 nsst%XT = reshape(dummy, (/lensfc/))
1531 error=nf90_inq_varid(ncid,
"xtts", id_var)
1532 CALL netcdf_err(error,
'READING xtts ID' )
1533 error=nf90_get_var(ncid, id_var, dummy)
1534 CALL netcdf_err(error,
'READING xtts' )
1535 nsst%XTTS = reshape(dummy, (/lensfc/))
1537 error=nf90_inq_varid(ncid,
"xu", id_var)
1538 CALL netcdf_err(error,
'READING xu ID' )
1539 error=nf90_get_var(ncid, id_var, dummy)
1540 CALL netcdf_err(error,
'READING xu' )
1541 nsst%XU = reshape(dummy, (/lensfc/))
1543 error=nf90_inq_varid(ncid,
"xv", id_var)
1544 CALL netcdf_err(error,
'READING xv ID' )
1545 error=nf90_get_var(ncid, id_var, dummy)
1546 CALL netcdf_err(error,
'READING xv' )
1547 nsst%XV = reshape(dummy, (/lensfc/))
1549 error=nf90_inq_varid(ncid,
"xz", id_var)
1550 CALL netcdf_err(error,
'READING xz ID' )
1551 error=nf90_get_var(ncid, id_var, dummy)
1552 CALL netcdf_err(error,
'READING xz' )
1553 nsst%XZ = reshape(dummy, (/lensfc/))
1555 error=nf90_inq_varid(ncid,
"xzts", id_var)
1556 CALL netcdf_err(error,
'READING xzts ID' )
1557 error=nf90_get_var(ncid, id_var, dummy)
1558 CALL netcdf_err(error,
'READING xzts' )
1559 nsst%XZTS = reshape(dummy, (/lensfc/))
1561 error=nf90_inq_varid(ncid,
"z_c", id_var)
1562 CALL netcdf_err(error,
'READING z_c ID' )
1563 error=nf90_get_var(ncid, id_var, dummy)
1564 CALL netcdf_err(error,
'READING z_c' )
1565 nsst%Z_C = reshape(dummy, (/lensfc/))
1567 error=nf90_inq_varid(ncid,
"zm", id_var)
1568 CALL netcdf_err(error,
'READING zm ID' )
1569 error=nf90_get_var(ncid, id_var, dummy)
1570 CALL netcdf_err(error,
'READING zm' )
1571 nsst%ZM = reshape(dummy, (/lensfc/))
1576 ALLOCATE(dummy3d(idim,jdim,lsoil))
1578 IF (
PRESENT(smcfcs))
THEN 1579 error=nf90_inq_varid(ncid,
"smc", id_var)
1580 CALL netcdf_err(error,
'READING smc ID' )
1581 error=nf90_get_var(ncid, id_var, dummy3d)
1582 CALL netcdf_err(error,
'READING smc' )
1583 smcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1586 IF (
PRESENT(slcfcs))
THEN 1587 error=nf90_inq_varid(ncid,
"slc", id_var)
1588 CALL netcdf_err(error,
'READING slc ID' )
1589 error=nf90_get_var(ncid, id_var, dummy3d)
1590 CALL netcdf_err(error,
'READING slc' )
1591 slcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1594 IF (
PRESENT(stcfcs))
THEN 1595 error=nf90_inq_varid(ncid,
"stc", id_var)
1596 CALL netcdf_err(error,
'READING stc ID' )
1597 error=nf90_get_var(ncid, id_var, dummy3d)
1598 CALL netcdf_err(error,
'READING stc' )
1599 stcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1603 IF (jedi_incr_file)
THEN 1604 IF (
PRESENT(slcinc))
THEN 1605 error=nf90_inq_varid(ncid,
"soill", id_var)
1606 CALL netcdf_err(error,
'READING soill ID' )
1607 error=nf90_get_var(ncid, id_var, dummy3d)
1608 CALL netcdf_err(error,
'READING slc increments' )
1609 slcinc = reshape(dummy3d, (/lensfc,lsoil/))
1612 IF (
PRESENT(stcinc))
THEN 1613 error=nf90_inq_varid(ncid,
"soilt", id_var)
1614 CALL netcdf_err(error,
'READING soilt ID' )
1615 error=nf90_get_var(ncid, id_var, dummy3d)
1616 CALL netcdf_err(error,
'READING stc increments' )
1617 stcinc = reshape(dummy3d, (/lensfc,lsoil/))
1621 IF (
PRESENT(stcinc))
THEN 1622 IF (.NOT.
PRESENT(lsoil_incr))
THEN 1623 write(6,*)
'FATAL ERROR variable lsoil_incr not declared.' 1624 CALL mpi_abort(mpi_comm_world, 134, error)
1626 DO k = 1, lsoil_incr
1627 WRITE(k_ch,
'(I1)') k
1629 incvar =
"soilt"//k_ch//
"_inc" 1630 print *,
"reading", incvar
1631 error=nf90_inq_varid(ncid,trim(incvar), id_var)
1632 CALL netcdf_err(error,
'READING soilt*_inc ID')
1633 error=nf90_get_var(ncid, id_var, dummy)
1634 CALL netcdf_err(error,
'READING soilt*_inc increments')
1636 stcinc(:,k) = reshape(dummy, (/lensfc/))
1639 IF (
PRESENT(slcinc))
THEN 1640 IF (.NOT.
PRESENT(lsoil_incr))
THEN 1641 write(6,*)
'FATAL ERROR variable lsoil_incr not declared.' 1642 CALL mpi_abort(mpi_comm_world, 136, error)
1644 DO k = 1, lsoil_incr
1645 WRITE(k_ch,
'(I1)') k
1647 incvar =
"slc"//k_ch//
"_inc" 1648 print *,
"reading", trim(incvar)
1649 error=nf90_inq_varid(ncid, trim(incvar), id_var)
1650 CALL netcdf_err(error,
'READING slc*_inc ID')
1651 error=nf90_get_var(ncid, id_var, dummy)
1652 CALL netcdf_err(error,
'READING slc*_inc increments')
1654 slcinc(:,k) = reshape(dummy, (/lensfc/))
1664 IF (
PRESENT(cvfcs)) cvfcs = 0.0
1665 IF (
PRESENT(cvtfcs)) cvtfcs = 0.0
1666 IF (
PRESENT(cvbfcs)) cvbfcs = 0.0
1671 IF (
PRESENT(zsoil))
THEN 1678 error = nf90_close(ncid)
1680 END SUBROUTINE read_data
1698 subroutine read_tf_clim_grb(file_sst,sst,rlats_sst,rlons_sst,mlat_sst,mlon_sst,mon)
1705 character(*) ,
intent(in ) :: file_sst
1706 integer ,
intent(in ) :: mon,mlat_sst,mlon_sst
1707 real,
dimension(mlat_sst) ,
intent( out) :: rlats_sst
1708 real,
dimension(mlon_sst) ,
intent( out) :: rlons_sst
1709 real,
dimension(mlon_sst,mlat_sst) ,
intent( out) :: sst
1712 integer,
parameter:: lu_sst = 21
1713 real,
parameter :: deg2rad = 3.141593/180.0
1716 logical(1),
allocatable,
dimension(:) :: lb
1720 integer :: iret,ni,nj
1721 integer :: mscan,kb1,ierr
1722 integer :: jincdir,i,iincdir,kb2,kb3,kf,kg,k,j,jf
1723 integer,
dimension(22):: jgds,kgds
1724 integer,
dimension(25):: jpds,kpds
1729 real,
allocatable,
dimension(:) :: f
1734 write(*,*)
' sstclm : ',file_sst
1735 call baopenr(lu_sst,trim(file_sst),iret)
1736 if (iret /= 0 )
then 1737 write(6,*)
'FATAL ERROR in read_tf_clm_grb: error opening sst file.' 1738 CALL mpi_abort(mpi_comm_world, 111, ierr)
1748 call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1757 allocate(lb(nlat_sst*nlon_sst))
1758 allocate(f(nlat_sst*nlon_sst))
1759 jf=nlat_sst*nlon_sst
1762 call getgb(lu_sst,0,jf,j,jpds,jgds,kf,k,kpds,kgds,lb,f,iret)
1764 write(6,*)
'FATAL ERROR in read_tf_clm_grb: error reading sst analysis data record.' 1766 CALL mpi_abort(mpi_comm_world, 111, ierr)
1769 if ( (nlat_sst /= mlat_sst) .or. (nlon_sst /= mlon_sst) )
then 1770 write(6,*)
'FATAL ERROR in read_rtg_org: inconsistent dimensions. mlat_sst,mlon_sst=',&
1771 mlat_sst,mlon_sst,
' -versus- nlat_sst,nlon_sst=',nlat_sst,nlon_sst
1773 CALL mpi_abort(mpi_comm_world, 111, ierr)
1779 dres = 180.0/
real(nlat_sst)
1780 ysst0 = 0.5*dres-90.0
1785 rlats_sst(j) = ysst0 +
real(j-1)*dres
1789 rlons_sst(i) = (xsst0 +
real(i-1)*dres)
1797 kb1=ibits(mscan,7,1)
1798 kb2=ibits(mscan,6,1)
1799 kb3=ibits(mscan,5,1)
1814 i=(ni+1)*kb1+(mod(k-1,ni)+1)*iincdir
1815 j=(nj+1)*(1-kb2)+jincdir*((k-1)/ni+1)
1817 j=(nj+1)*(1-kb2)+(mod(k-1,nj)+1)*jincdir
1818 i=(ni+1)*kb1+iincdir*((k-1)/nj+1)
1825 call baclose(lu_sst,iret)
1826 if (iret /= 0 )
then 1827 write(6,*)
'FATAL ERROR in read_tf_clm_grb: error closing sst file.' 1828 CALL mpi_abort(mpi_comm_world, 121, ierr)
1831 end subroutine read_tf_clim_grb
1840 subroutine get_tf_clm_dim(file_sst,mlat_sst,mlon_sst)
1846 character(*) ,
intent(in ) :: file_sst
1847 integer ,
intent(out) :: mlat_sst,mlon_sst
1850 integer,
parameter:: lu_sst = 21
1853 integer :: kf,kg,k,j,ierr
1854 integer,
dimension(22):: jgds,kgds
1855 integer,
dimension(25):: jpds,kpds
1860 call baopenr(lu_sst,trim(file_sst),iret)
1861 if (iret /= 0 )
then 1862 write(6,*)
'FATAL ERROR in get_tf_clm_dim: error opening sst file.' 1863 CALL mpi_abort(mpi_comm_world, 111, ierr)
1873 call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1878 write(*,*)
'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1880 call baclose(lu_sst,iret)
1881 if (iret /= 0 )
then 1882 write(6,*)
'FATAL ERROR in get_tf_clm_dim: error closing sst file.' 1883 CALL mpi_abort(mpi_comm_world, 121, ierr)
1885 end subroutine get_tf_clm_dim
1898 subroutine read_salclm_gfs_nc(filename,sal,xlats,xlons,nlat,nlon,itime)
1903 character (len=*),
intent(in) :: filename
1904 integer,
intent(in) :: nlat,nlon
1905 integer,
intent(in) :: itime
1906 real,
dimension(nlat),
intent(out) :: xlats
1907 real,
dimension(nlon),
intent(out) :: xlons
1908 real,
dimension(nlon,nlat),
intent(out) :: sal
1912 integer,
parameter :: ndims = 3
1913 character (len = *),
parameter :: lat_name =
"latitude" 1914 character (len = *),
parameter :: lon_name =
"longitude" 1915 character (len = *),
parameter :: t_name =
"time" 1916 character (len = *),
parameter :: sal_name=
"sal" 1917 integer :: time_varid,lon_varid, lat_varid, sal_varid
1920 integer,
dimension(ndims) :: start, count
1922 character (len = *),
parameter :: units =
"units" 1923 character (len = *),
parameter :: sal_units =
"psu" 1925 character (len = *),
parameter :: time_units =
"months" 1926 character (len = *),
parameter :: lat_units =
"degrees_north" 1927 character (len = *),
parameter :: lon_units =
"degrees_east" 1930 call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1933 call nc_check( nf90_inq_varid(ncid, t_name, time_varid) )
1934 call nc_check( nf90_inq_varid(ncid, lat_name, lat_varid) )
1935 call nc_check( nf90_inq_varid(ncid, lon_name, lon_varid) )
1939 call nc_check( nf90_get_var(ncid, lat_varid, xlats) )
1940 call nc_check( nf90_get_var(ncid, lon_varid, xlons) )
1943 call nc_check( nf90_inq_varid(ncid, sal_name,sal_varid) )
1947 start = (/ 1, 1, itime /)
1948 count = (/ nlon, nlat, 1 /)
1952 call nc_check( nf90_get_var(ncid, sal_varid, sal, start, count) )
1956 call nc_check( nf90_close(ncid) )
1961 end subroutine read_salclm_gfs_nc
1969 subroutine get_dim_nc(filename,nlat,nlon)
1973 character (len=*),
intent(in) :: filename
1974 integer,
intent(out) :: nlat,nlon
1976 character (len = *),
parameter :: lat_name =
"latitude" 1977 character (len = *),
parameter :: lon_name =
"longitude" 1979 integer :: latdimid,londimid
1982 call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1985 call nc_check( nf90_inq_dimid(ncid,lat_name,latdimid) )
1986 call nc_check( nf90_inq_dimid(ncid,lon_name,londimid) )
1987 call nc_check( nf90_inquire_dimension(ncid,latdimid,len=nlat) )
1988 call nc_check( nf90_inquire_dimension(ncid,londimid,len=nlon) )
1994 call nc_check( nf90_close(ncid) )
1999 end subroutine get_dim_nc
2006 subroutine nc_check(status)
2011 integer,
intent ( in) :: status
2014 if(status /= nf90_noerr)
then 2015 print *,
'FATAL ERROR:' 2016 print *, trim(nf90_strerror(status))
2017 CALL mpi_abort(mpi_comm_world, 122, ierr)
2019 end subroutine nc_check
2021 END MODULE read_write_data