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