global_cycle  1.10.0
read_write_data.f90
Go to the documentation of this file.
1 
5 
8 MODULE read_write_data
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  PUBLIC :: read_data
56  PUBLIC :: read_gsi_data
57  PUBLIC :: read_lat_lon_orog
58  PUBLIC :: write_data
59  public :: read_tf_clim_grb,get_tf_clm_dim
60  public :: read_salclm_gfs_nc,get_dim_nc
61 
62  CONTAINS
63 
116 
117  subroutine write_data(lensfc,idim,jdim,lsoil, &
118  do_nsst,nsst,slifcs,tsffcs,vegfcs,swefcs, &
119  tg3fcs,zorfcs,albfcs,alffcs, &
120  cnpfcs,f10m,t2m,q2m,vetfcs, &
121  sotfcs,ustar,fmm,fhh,sicfcs, &
122  sihfcs,sitfcs,tprcp,srflag, &
123  swdfcs,vmnfcs,vmxfcs,slpfcs, &
124  absfcs,slcfcs,smcfcs,stcfcs)
126  use mpi
127 
128  implicit none
129 
130  integer, intent(in) :: lensfc, lsoil
131  integer, intent(in) :: idim, jdim
132 
133  logical, intent(in) :: do_nsst
134 
135  real, intent(in), optional :: slifcs(lensfc),tsffcs(lensfc)
136  real, intent(in), optional :: swefcs(lensfc),tg3fcs(lensfc)
137  real, intent(in), optional :: zorfcs(lensfc),albfcs(lensfc,4)
138  real, intent(in), optional :: alffcs(lensfc,2),cnpfcs(lensfc)
139  real, intent(in), optional :: f10m(lensfc),t2m(lensfc)
140  real, intent(in), optional :: q2m(lensfc),vegfcs(lensfc)
141  real, intent(in), optional :: vetfcs(lensfc),sotfcs(lensfc)
142  real, intent(in), optional :: ustar(lensfc),fmm(lensfc)
143  real, intent(in), optional :: fhh(lensfc), sicfcs(lensfc)
144  real, intent(in), optional :: sihfcs(lensfc), sitfcs(lensfc)
145  real, intent(in), optional :: tprcp(lensfc), srflag(lensfc)
146  real, intent(in), optional :: swdfcs(lensfc), vmnfcs(lensfc)
147  real, intent(in), optional :: vmxfcs(lensfc), slpfcs(lensfc)
148  real, intent(in), optional :: absfcs(lensfc), slcfcs(lensfc,lsoil)
149  real, intent(in), optional :: smcfcs(lensfc,lsoil), stcfcs(lensfc,lsoil)
150 
151  type(nsst_data), intent(in) :: nsst
152 
153  integer :: dim_x, dim_y, dim_time, dims_3d(3)
154 
155  real :: dum2d(idim,jdim), dum3d(idim,jdim,lsoil)
156 
157  character(len=50) :: fnbgso
158  character(len=3) :: rankch
159 
160  integer :: myrank, error, ncid, id_var
161 
162  call mpi_comm_rank(mpi_comm_world, myrank, error)
163 
164  write(rankch, '(i3.3)') (myrank+1)
165 
166  fnbgso = "./fnbgso." // rankch
167 
168  print*
169  print*,"update OUTPUT SFC DATA TO: ",trim(fnbgso)
170 
171  error=nf90_open(trim(fnbgso),nf90_write,ncid)
172  CALL netcdf_err(error, 'OPENING FILE: '//trim(fnbgso) )
173 
174  if(present(slifcs)) then
175  error=nf90_inq_varid(ncid, "slmsk", id_var)
176  call netcdf_err(error, 'reading slmsk id' )
177  dum2d = reshape(slifcs, (/idim,jdim/))
178  error = nf90_put_var( ncid, id_var, dum2d)
179  call netcdf_err(error, 'writing slmsk record' )
180  call remove_checksum(ncid, id_var)
181  endif
182 
183  if(present(tsffcs)) then
184  error=nf90_inq_varid(ncid, "tsea", id_var)
185  call netcdf_err(error, 'reading tsea id' )
186  dum2d = reshape(tsffcs, (/idim,jdim/))
187  error = nf90_put_var( ncid, id_var, dum2d)
188  call netcdf_err(error, 'writing tsea record' )
189  call remove_checksum(ncid, id_var)
190  endif
191 
192  if(present(swefcs)) then
193  error=nf90_inq_varid(ncid, "sheleg", id_var)
194  call netcdf_err(error, 'reading sheleg id' )
195  dum2d = reshape(swefcs, (/idim,jdim/))
196  error = nf90_put_var( ncid, id_var, dum2d)
197  call netcdf_err(error, 'writing sheleg record' )
198  call remove_checksum(ncid, id_var)
199  endif
200 
201  if(present(tg3fcs)) then
202  error=nf90_inq_varid(ncid, "tg3", id_var)
203  call netcdf_err(error, 'reading tg3 id' )
204  dum2d = reshape(tg3fcs, (/idim,jdim/))
205  error = nf90_put_var( ncid, id_var, dum2d)
206  call netcdf_err(error, 'writing tg3 record' )
207  call remove_checksum(ncid, id_var)
208  endif
209 
210  if(present(zorfcs)) then
211  error=nf90_inq_varid(ncid, "zorl", id_var)
212  call netcdf_err(error, 'reading zorl id' )
213  dum2d = reshape(zorfcs, (/idim,jdim/))
214  error = nf90_put_var( ncid, id_var, dum2d)
215  call netcdf_err(error, 'writing zorl record' )
216  call remove_checksum(ncid, id_var)
217  endif
218 
219  if(present(albfcs)) then
220  error=nf90_inq_varid(ncid, "alvsf", id_var)
221  call netcdf_err(error, 'reading alvsf id' )
222  dum2d = reshape(albfcs(:,1), (/idim,jdim/))
223  error = nf90_put_var( ncid, id_var, dum2d)
224  call netcdf_err(error, 'writing alvsf record' )
225  call remove_checksum(ncid, id_var)
226 
227  error=nf90_inq_varid(ncid, "alvwf", id_var)
228  call netcdf_err(error, 'reading alvwf id' )
229  dum2d = reshape(albfcs(:,2), (/idim,jdim/))
230  error = nf90_put_var( ncid, id_var, dum2d)
231  call netcdf_err(error, 'writing alvwf record' )
232  call remove_checksum(ncid, id_var)
233 
234  error=nf90_inq_varid(ncid, "alnsf", id_var)
235  call netcdf_err(error, 'reading alnsf id' )
236  dum2d = reshape(albfcs(:,3), (/idim,jdim/))
237  error = nf90_put_var( ncid, id_var, dum2d)
238  call netcdf_err(error, 'writing alnsf record' )
239  call remove_checksum(ncid, id_var)
240 
241  error=nf90_inq_varid(ncid, "alnwf", id_var)
242  call netcdf_err(error, 'reading alnwf id' )
243  dum2d = reshape(albfcs(:,4), (/idim,jdim/))
244  error = nf90_put_var( ncid, id_var, dum2d)
245  call netcdf_err(error, 'writing alnwf record' )
246  call remove_checksum(ncid, id_var)
247  endif
248 
249  if(present(alffcs)) then
250  error=nf90_inq_varid(ncid, "facsf", id_var)
251  call netcdf_err(error, 'reading facsf id' )
252  dum2d = reshape(alffcs(:,1), (/idim,jdim/))
253  error = nf90_put_var( ncid, id_var, dum2d)
254  call netcdf_err(error, 'writing facsf record' )
255  call remove_checksum(ncid, id_var)
256 
257  error=nf90_inq_varid(ncid, "facwf", id_var)
258  call netcdf_err(error, 'reading facwf id' )
259  dum2d = reshape(alffcs(:,2), (/idim,jdim/))
260  error = nf90_put_var( ncid, id_var, dum2d)
261  call netcdf_err(error, 'writing facwf record' )
262  call remove_checksum(ncid, id_var)
263  endif
264 
265  if(present(vegfcs)) then
266  error=nf90_inq_varid(ncid, "vfrac", id_var)
267  call netcdf_err(error, 'reading vfrac id' )
268  dum2d = reshape(vegfcs, (/idim,jdim/))
269  error = nf90_put_var( ncid, id_var, dum2d)
270  call netcdf_err(error, 'writing vegfcs record' )
271  call remove_checksum(ncid, id_var)
272  endif
273 
274  if(present(cnpfcs)) then
275  error=nf90_inq_varid(ncid, "canopy", id_var)
276  call netcdf_err(error, 'reading canopy id' )
277  dum2d = reshape(cnpfcs, (/idim,jdim/))
278  error = nf90_put_var( ncid, id_var, dum2d)
279  call netcdf_err(error, 'writing canopy record' )
280  call remove_checksum(ncid, id_var)
281  endif
282 
283  if(present(f10m)) then
284  error=nf90_inq_varid(ncid, "f10m", id_var)
285  call netcdf_err(error, 'reading f10m id' )
286  dum2d = reshape(f10m, (/idim,jdim/))
287  error = nf90_put_var( ncid, id_var, dum2d)
288  call netcdf_err(error, 'writing f10m record' )
289  call remove_checksum(ncid, id_var)
290  endif
291 
292  if(present(t2m)) then
293  error=nf90_inq_varid(ncid, "t2m", id_var)
294  call netcdf_err(error, 'reading t2m id' )
295  dum2d = reshape(t2m, (/idim,jdim/))
296  error = nf90_put_var( ncid, id_var, dum2d)
297  call netcdf_err(error, 'writing t2m record' )
298  call remove_checksum(ncid, id_var)
299  endif
300 
301  if(present(q2m)) then
302  error=nf90_inq_varid(ncid, "q2m", id_var)
303  call netcdf_err(error, 'reading q2m id' )
304  dum2d = reshape(q2m, (/idim,jdim/))
305  error = nf90_put_var( ncid, id_var, dum2d)
306  call netcdf_err(error, 'writing q2m record' )
307  call remove_checksum(ncid, id_var)
308  endif
309 
310  if(present(vetfcs)) then
311  error=nf90_inq_varid(ncid, "vtype", id_var)
312  call netcdf_err(error, 'reading vtype id' )
313  dum2d = reshape(vetfcs, (/idim,jdim/))
314  error = nf90_put_var( ncid, id_var, dum2d)
315  call netcdf_err(error, 'writing vtype record' )
316  call remove_checksum(ncid, id_var)
317  endif
318 
319  if(present(sotfcs)) then
320  error=nf90_inq_varid(ncid, "stype", id_var)
321  call netcdf_err(error, 'reading stype id' )
322  dum2d = reshape(sotfcs, (/idim,jdim/))
323  error = nf90_put_var( ncid, id_var, dum2d)
324  call netcdf_err(error, 'writing stype record' )
325  call remove_checksum(ncid, id_var)
326  endif
327 
328  if(present(ustar)) then
329  error=nf90_inq_varid(ncid, "uustar", id_var)
330  call netcdf_err(error, 'reading uustar id' )
331  dum2d = reshape(ustar, (/idim,jdim/))
332  error = nf90_put_var( ncid, id_var, dum2d)
333  call netcdf_err(error, 'writing uustar record' )
334  call remove_checksum(ncid, id_var)
335  endif
336 
337  if(present(fmm)) then
338  error=nf90_inq_varid(ncid, "ffmm", id_var)
339  call netcdf_err(error, 'reading ffmm id' )
340  dum2d = reshape(fmm, (/idim,jdim/))
341  error = nf90_put_var( ncid, id_var, dum2d)
342  call netcdf_err(error, 'writing ffmm record' )
343  call remove_checksum(ncid, id_var)
344  endif
345 
346  if(present(fhh)) then
347  error=nf90_inq_varid(ncid, "ffhh", id_var)
348  call netcdf_err(error, 'reading ffhh id' )
349  dum2d = reshape(fhh, (/idim,jdim/))
350  error = nf90_put_var( ncid, id_var, dum2d)
351  call netcdf_err(error, 'writing ffhh record' )
352  call remove_checksum(ncid, id_var)
353  endif
354 
355  if(present(sicfcs)) then
356  error=nf90_inq_varid(ncid, "fice", id_var)
357  call netcdf_err(error, 'reading fice id' )
358  dum2d = reshape(sicfcs, (/idim,jdim/))
359  error = nf90_put_var( ncid, id_var, dum2d)
360  call netcdf_err(error, 'writing fice record' )
361  call remove_checksum(ncid, id_var)
362  endif
363 
364  if(present(sihfcs)) then
365  error=nf90_inq_varid(ncid, "hice", id_var)
366  call netcdf_err(error, 'reading hice id' )
367  dum2d = reshape(sihfcs, (/idim,jdim/))
368  error = nf90_put_var( ncid, id_var, dum2d)
369  call netcdf_err(error, 'writing hice record' )
370  call remove_checksum(ncid, id_var)
371  endif
372 
373  if(present(sitfcs)) then
374  error=nf90_inq_varid(ncid, "tisfc", id_var)
375  call netcdf_err(error, 'reading tisfc id' )
376  dum2d = reshape(sitfcs, (/idim,jdim/))
377  error = nf90_put_var( ncid, id_var, dum2d)
378  call netcdf_err(error, 'writing tisfc record' )
379  call remove_checksum(ncid, id_var)
380  endif
381 
382  if(present(tprcp)) then
383  error=nf90_inq_varid(ncid, "tprcp", id_var)
384  call netcdf_err(error, 'reading tprcp id' )
385  dum2d = reshape(tprcp, (/idim,jdim/))
386  error = nf90_put_var( ncid, id_var, dum2d)
387  call netcdf_err(error, 'writing tprcp record' )
388  call remove_checksum(ncid, id_var)
389  endif
390 
391  if(present(srflag)) then
392  error=nf90_inq_varid(ncid, "srflag", id_var)
393  call netcdf_err(error, 'reading srflag id' )
394  dum2d = reshape(srflag, (/idim,jdim/))
395  error = nf90_put_var( ncid, id_var, dum2d)
396  call netcdf_err(error, 'writing srflag record' )
397  call remove_checksum(ncid, id_var)
398  endif
399 
400  if(present(swdfcs)) then
401  error=nf90_inq_varid(ncid, "snwdph", id_var)
402  call netcdf_err(error, 'reading snwdph id' )
403  dum2d = reshape(swdfcs, (/idim,jdim/))
404  error = nf90_put_var( ncid, id_var, dum2d)
405  call netcdf_err(error, 'writing snwdph record' )
406  call remove_checksum(ncid, id_var)
407  endif
408 
409  if(present(vmnfcs)) then
410  error=nf90_inq_varid(ncid, "shdmin", id_var)
411  call netcdf_err(error, 'reading shdmin id' )
412  dum2d = reshape(vmnfcs, (/idim,jdim/))
413  error = nf90_put_var( ncid, id_var, dum2d)
414  call netcdf_err(error, 'writing shdmin record' )
415  call remove_checksum(ncid, id_var)
416  endif
417 
418  if(present(vmxfcs)) then
419  error=nf90_inq_varid(ncid, "shdmax", id_var)
420  call netcdf_err(error, 'reading shdmax id' )
421  dum2d = reshape(vmxfcs, (/idim,jdim/))
422  error = nf90_put_var( ncid, id_var, dum2d)
423  call netcdf_err(error, 'writing shdmax record' )
424  call remove_checksum(ncid, id_var)
425  endif
426 
427  if(present(slpfcs)) then
428  error=nf90_inq_varid(ncid, "slope", id_var)
429  call netcdf_err(error, 'reading slope id' )
430  dum2d = reshape(slpfcs, (/idim,jdim/))
431  error = nf90_put_var( ncid, id_var, dum2d)
432  call netcdf_err(error, 'writing slope record' )
433  call remove_checksum(ncid, id_var)
434  endif
435 
436  if(present(absfcs)) then
437  error=nf90_inq_varid(ncid, "snoalb", id_var)
438  call netcdf_err(error, 'reading snoalb id' )
439  dum2d = reshape(absfcs, (/idim,jdim/))
440  error = nf90_put_var( ncid, id_var, dum2d)
441  call netcdf_err(error, 'writing snoalb record' )
442  call remove_checksum(ncid, id_var)
443  endif
444 
445  if(present(slcfcs)) then
446  error=nf90_inq_varid(ncid, "slc", id_var)
447  call netcdf_err(error, 'reading slc id' )
448  dum3d = reshape(slcfcs, (/idim,jdim,lsoil/))
449  error = nf90_put_var( ncid, id_var, dum3d)
450  call netcdf_err(error, 'writing slc record' )
451  call remove_checksum(ncid, id_var)
452  endif
453 
454  if(present(smcfcs)) then
455  error=nf90_inq_varid(ncid, "smc", id_var)
456  call netcdf_err(error, 'reading smc id' )
457  dum3d = reshape(smcfcs, (/idim,jdim,lsoil/))
458  error = nf90_put_var( ncid, id_var, dum3d)
459  call netcdf_err(error, 'writing smc record' )
460  call remove_checksum(ncid, id_var)
461  endif
462 
463  if(present(stcfcs)) then
464  error=nf90_inq_varid(ncid, "stc", id_var)
465  call netcdf_err(error, 'reading stc id' )
466  dum3d = reshape(stcfcs, (/idim,jdim,lsoil/))
467  error = nf90_put_var( ncid, id_var, dum3d)
468  call netcdf_err(error, 'writing stc record' )
469  call remove_checksum(ncid, id_var)
470  endif
471 
472  if(do_nsst) then
473 
474  error=nf90_inq_varid(ncid, "tref", id_var)
475  call netcdf_err(error, 'reading tref id' )
476  dum2d = reshape(nsst%tref, (/idim,jdim/))
477  error = nf90_put_var( ncid, id_var, dum2d)
478  call netcdf_err(error, 'WRITING TREF RECORD' )
479  call remove_checksum(ncid, id_var)
480 
481  error=nf90_inq_varid(ncid, "z_c", id_var)
482  call netcdf_err(error, 'reading z_c id' )
483  dum2d = reshape(nsst%z_c, (/idim,jdim/))
484  error = nf90_put_var( ncid, id_var, dum2d)
485  call netcdf_err(error, 'WRITING Z_C RECORD' )
486  call remove_checksum(ncid, id_var)
487 
488  error=nf90_inq_varid(ncid, "c_0", id_var)
489  call netcdf_err(error, 'reading c_0 id' )
490  dum2d = reshape(nsst%c_0, (/idim,jdim/))
491  error = nf90_put_var( ncid, id_var, dum2d)
492  call netcdf_err(error, 'WRITING C_0 RECORD' )
493  call remove_checksum(ncid, id_var)
494 
495  error=nf90_inq_varid(ncid, "c_d", id_var)
496  call netcdf_err(error, 'reading c_d id' )
497  dum2d = reshape(nsst%c_d, (/idim,jdim/))
498  error = nf90_put_var( ncid, id_var, dum2d)
499  call netcdf_err(error, 'WRITING C_D RECORD' )
500  call remove_checksum(ncid, id_var)
501 
502  error=nf90_inq_varid(ncid, "w_0", id_var)
503  call netcdf_err(error, 'reading w_0 id' )
504  dum2d = reshape(nsst%w_0, (/idim,jdim/))
505  error = nf90_put_var( ncid, id_var, dum2d)
506  call netcdf_err(error, 'WRITING W_0 RECORD' )
507  call remove_checksum(ncid, id_var)
508 
509  error=nf90_inq_varid(ncid, "w_d", id_var)
510  call netcdf_err(error, 'reading w_d id' )
511  dum2d = reshape(nsst%w_d, (/idim,jdim/))
512  error = nf90_put_var( ncid, id_var, dum2d)
513  call netcdf_err(error, 'WRITING W_D RECORD' )
514  call remove_checksum(ncid, id_var)
515 
516  error=nf90_inq_varid(ncid, "xt", id_var)
517  call netcdf_err(error, 'reading xt id' )
518  dum2d = reshape(nsst%xt, (/idim,jdim/))
519  error = nf90_put_var( ncid, id_var, dum2d)
520  call netcdf_err(error, 'WRITING XT RECORD' )
521  call remove_checksum(ncid, id_var)
522 
523  error=nf90_inq_varid(ncid, "xs", id_var)
524  call netcdf_err(error, 'reading xs id' )
525  dum2d = reshape(nsst%xs, (/idim,jdim/))
526  error = nf90_put_var( ncid, id_var, dum2d)
527  call netcdf_err(error, 'WRITING XS RECORD' )
528  call remove_checksum(ncid, id_var)
529 
530  error=nf90_inq_varid(ncid, "xu", id_var)
531  call netcdf_err(error, 'reading xu id' )
532  dum2d = reshape(nsst%xu, (/idim,jdim/))
533  error = nf90_put_var( ncid, id_var, dum2d)
534  call netcdf_err(error, 'WRITING XU RECORD' )
535  call remove_checksum(ncid, id_var)
536 
537  error=nf90_inq_varid(ncid, "xv", id_var)
538  call netcdf_err(error, 'reading xv id' )
539  dum2d = reshape(nsst%xv, (/idim,jdim/))
540  error = nf90_put_var( ncid, id_var, dum2d)
541  call netcdf_err(error, 'WRITING XV RECORD' )
542  call remove_checksum(ncid, id_var)
543 
544  error=nf90_inq_varid(ncid, "xz", id_var)
545  call netcdf_err(error, 'reading xz id' )
546  dum2d = reshape(nsst%xz, (/idim,jdim/))
547  error = nf90_put_var( ncid, id_var, dum2d)
548  call netcdf_err(error, 'WRITING XZ RECORD' )
549  call remove_checksum(ncid, id_var)
550 
551  error=nf90_inq_varid(ncid, "zm", id_var)
552  call netcdf_err(error, 'reading zm id' )
553  dum2d = reshape(nsst%zm, (/idim,jdim/))
554  error = nf90_put_var( ncid, id_var, dum2d)
555  call netcdf_err(error, 'WRITING ZM RECORD' )
556  call remove_checksum(ncid, id_var)
557 
558  error=nf90_inq_varid(ncid, "xtts", id_var)
559  call netcdf_err(error, 'reading xtts id' )
560  dum2d = reshape(nsst%xtts, (/idim,jdim/))
561  error = nf90_put_var( ncid, id_var, dum2d)
562  call netcdf_err(error, 'WRITING XTTS RECORD' )
563  call remove_checksum(ncid, id_var)
564 
565  error=nf90_inq_varid(ncid, "xzts", id_var)
566  call netcdf_err(error, 'reading xzts id' )
567  dum2d = reshape(nsst%xzts, (/idim,jdim/))
568  error = nf90_put_var( ncid, id_var, dum2d)
569  call netcdf_err(error, 'WRITING XZTS RECORD' )
570  call remove_checksum(ncid, id_var)
571 
572  error=nf90_inq_varid(ncid, "d_conv", id_var)
573  call netcdf_err(error, 'reading d_conv id' )
574  dum2d = reshape(nsst%d_conv, (/idim,jdim/))
575  error = nf90_put_var( ncid, id_var, dum2d)
576  call netcdf_err(error, 'WRITING D_CONV RECORD' )
577  call remove_checksum(ncid, id_var)
578 
579  error=nf90_inq_varid(ncid, "ifd", id_var)
580  call netcdf_err(error, 'reading idf id' )
581  dum2d = reshape(nsst%ifd, (/idim,jdim/))
582  error = nf90_put_var( ncid, id_var, dum2d)
583  call netcdf_err(error, 'WRITING IFD RECORD' )
584  call remove_checksum(ncid, id_var)
585 
586  error=nf90_inq_varid(ncid, "dt_cool", id_var)
587  call netcdf_err(error, 'reading dt_cool id' )
588  dum2d = reshape(nsst%dt_cool, (/idim,jdim/))
589  error = nf90_put_var( ncid, id_var, dum2d)
590  call netcdf_err(error, 'WRITING DT_COOL RECORD' )
591  call remove_checksum(ncid, id_var)
592 
593  error=nf90_inq_varid(ncid, "qrain", id_var)
594  call netcdf_err(error, 'reading qrain id' )
595  dum2d = reshape(nsst%qrain, (/idim,jdim/))
596  error = nf90_put_var( ncid, id_var, dum2d)
597  call netcdf_err(error, 'WRITING QRAIN RECORD' )
598  call remove_checksum(ncid, id_var)
599 
600 ! Some files don't include 'tfinc', which is just diagnostic.
601 ! If missing, then add it to the restart file.
602  error=nf90_inq_varid(ncid, "tfinc", id_var)
603  if (error /= 0) then
604  error=nf90_inq_dimid(ncid, "xaxis_1", dim_x)
605  call netcdf_err(error, 'finding xaxis_1' )
606  error=nf90_inq_dimid(ncid, "yaxis_1", dim_y)
607  call netcdf_err(error, 'finding yaxis_1' )
608  error=nf90_inq_dimid(ncid, "Time", dim_time)
609  call netcdf_err(error, 'finding Time' )
610  dims_3d(1) = dim_x
611  dims_3d(2) = dim_y
612  dims_3d(3) = dim_time
613  error=nf90_redef(ncid)
614  error = nf90_def_var(ncid, 'tfinc', nf90_double, dims_3d, id_var)
615  call netcdf_err(error, 'DEFINING tfinc' )
616  error = nf90_put_att(ncid, id_var, "long_name", "tfinc")
617  call netcdf_err(error, 'DEFINING tfinc LONG NAME' )
618  error = nf90_put_att(ncid, id_var, "units", "none")
619  call netcdf_err(error, 'DEFINING tfinc UNITS' )
620  error=nf90_enddef(ncid)
621  endif
622  dum2d = reshape(nsst%tfinc, (/idim,jdim/))
623  error = nf90_put_var( ncid, id_var, dum2d)
624  call netcdf_err(error, 'WRITING TFINC RECORD' )
625 
626  endif
627 
628  error = nf90_close(ncid)
629 
630  end subroutine write_data
631 
638  subroutine remove_checksum(ncid, id_var)
640  implicit none
641 
642  integer, intent(in) :: ncid, id_var
643 
644  integer :: error
645 
646  error=nf90_inquire_attribute(ncid, id_var, 'checksum')
647 
648  if (error == 0) then ! attribute was found
649 
650  error = nf90_redef(ncid)
651  call netcdf_err(error, 'entering define mode' )
652 
653  error=nf90_del_att(ncid, id_var, 'checksum')
654  call netcdf_err(error, 'deleting checksum' )
655 
656  error= nf90_enddef(ncid)
657  call netcdf_err(error, 'ending define mode' )
658 
659  endif
660 
661  end subroutine remove_checksum
662 
676  SUBROUTINE read_lat_lon_orog(RLA,RLO,OROG,OROG_UF,&
677  TILE_NUM,IDIM,JDIM,IJDIM)
679  USE mpi
680 
681  IMPLICIT NONE
682 
683  INTEGER, INTENT(IN) :: idim, jdim, ijdim
684 
685  CHARACTER(LEN=5), INTENT(OUT) :: tile_num
686 
687  REAL, INTENT(OUT) :: rla(ijdim),rlo(ijdim)
688  REAL, INTENT(OUT) :: orog(ijdim),orog_uf(ijdim)
689 
690  CHARACTER(LEN=50) :: fnorog, fngrid
691  CHARACTER(LEN=3) :: rankch
692 
693  INTEGER :: error, ncid, ncid_orog
694  INTEGER :: i, ii, j, jj, myrank
695  INTEGER :: id_dim, id_var, nx, ny
696 
697  REAL, ALLOCATABLE :: dummy(:,:), geolat(:,:), geolon(:,:)
698  REAL(KIND=4), ALLOCATABLE :: dummy4(:,:)
699 
700  CALL mpi_comm_rank(mpi_comm_world, myrank, error)
701 
702  WRITE(rankch, '(I3.3)') (myrank+1)
703 
704  fngrid = "./fngrid." // rankch
705 
706  print*
707  print*, "READ FV3 GRID INFO FROM: "//trim(fngrid)
708 
709  error=nf90_open(trim(fngrid),nf90_nowrite,ncid)
710  CALL netcdf_err(error, 'OPENING FILE: '//trim(fngrid) )
711 
712  error=nf90_inq_dimid(ncid, 'nx', id_dim)
713  CALL netcdf_err(error, 'ERROR READING NX ID' )
714 
715  error=nf90_inquire_dimension(ncid,id_dim,len=nx)
716  CALL netcdf_err(error, 'ERROR READING NX' )
717 
718  error=nf90_inq_dimid(ncid, 'ny', id_dim)
719  CALL netcdf_err(error, 'ERROR READING NY ID' )
720 
721  error=nf90_inquire_dimension(ncid,id_dim,len=ny)
722  CALL netcdf_err(error, 'ERROR READING NY' )
723 
724  IF ((nx/2) /= idim .OR. (ny/2) /= jdim) THEN
725  print*,'FATAL ERROR: DIMENSIONS IN FILE: ',(nx/2),(ny/2)
726  print*,'DO NOT MATCH GRID DIMENSIONS: ',idim,jdim
727  CALL mpi_abort(mpi_comm_world, 130, error)
728  ENDIF
729 
730  ALLOCATE(geolon(nx+1,ny+1))
731  ALLOCATE(geolat(nx+1,ny+1))
732 
733  error=nf90_inq_varid(ncid, 'x', id_var)
734  CALL netcdf_err(error, 'ERROR READING X ID' )
735  error=nf90_get_var(ncid, id_var, geolon)
736  CALL netcdf_err(error, 'ERROR READING X RECORD' )
737 
738  error=nf90_inq_varid(ncid, 'y', id_var)
739  CALL netcdf_err(error, 'ERROR READING Y ID' )
740  error=nf90_get_var(ncid, id_var, geolat)
741  CALL netcdf_err(error, 'ERROR READING Y RECORD' )
742 
743  ALLOCATE(dummy(idim,jdim))
744 
745  DO j = 1, jdim
746  DO i = 1, idim
747  ii = 2*i
748  jj = 2*j
749  dummy(i,j) = geolon(ii,jj)
750  ENDDO
751  ENDDO
752 
753  rlo = reshape(dummy, (/ijdim/))
754 
755  DEALLOCATE(geolon)
756 
757  DO j = 1, jdim
758  DO i = 1, idim
759  ii = 2*i
760  jj = 2*j
761  dummy(i,j) = geolat(ii,jj)
762  ENDDO
763  ENDDO
764 
765  rla = reshape(dummy, (/ijdim/))
766 
767  DEALLOCATE(geolat, dummy)
768 
769  error=nf90_inq_varid(ncid, 'tile', id_var)
770  CALL netcdf_err(error, 'ERROR READING TILE ID' )
771  error=nf90_get_var(ncid, id_var, tile_num)
772  CALL netcdf_err(error, 'ERROR READING TILE RECORD' )
773 
774  error = nf90_close(ncid)
775 
776  fnorog = "./fnorog." // rankch
777 
778  print*
779  print*, "READ FV3 OROG INFO FROM: "//trim(fnorog)
780 
781  error=nf90_open(trim(fnorog),nf90_nowrite,ncid_orog)
782  CALL netcdf_err(error, 'OPENING FILE: '//trim(fnorog) )
783 
784  ALLOCATE(dummy4(idim,jdim))
785 
786  error=nf90_inq_varid(ncid_orog, 'orog_raw', id_var)
787  CALL netcdf_err(error, 'ERROR READING orog_raw ID' )
788  error=nf90_get_var(ncid_orog, id_var, dummy4)
789  CALL netcdf_err(error, 'ERROR READING orog_raw RECORD' )
790  orog_uf = reshape(dummy4, (/ijdim/))
791 
792  error=nf90_inq_varid(ncid_orog, 'orog_filt', id_var)
793  CALL netcdf_err(error, 'ERROR READING orog_filt ID' )
794  error=nf90_get_var(ncid_orog, id_var, dummy4)
795  CALL netcdf_err(error, 'ERROR READING orog_filt RECORD' )
796  orog = reshape(dummy4, (/ijdim/))
797 
798  DEALLOCATE(dummy4)
799 
800  error = nf90_close(ncid_orog)
801 
802  END SUBROUTINE read_lat_lon_orog
803 
810  SUBROUTINE netcdf_err( ERR, STRING )
812  USE mpi
813 
814  IMPLICIT NONE
815 
816  INTEGER, INTENT(IN) :: ERR
817  CHARACTER(LEN=*), INTENT(IN) :: STRING
818  CHARACTER(LEN=80) :: ERRMSG
819  INTEGER :: IRET
820 
821  IF( err == nf90_noerr )RETURN
822  errmsg = nf90_strerror(err)
823  print*,''
824  print*,'FATAL ERROR: ', trim(string), ': ', trim(errmsg)
825  print*,'STOP.'
826  CALL mpi_abort(mpi_comm_world, 999, iret)
827 
828  RETURN
829  END SUBROUTINE netcdf_err
830 
843  SUBROUTINE read_gsi_data(GSI_FILE, FILE_TYPE, LSOIL)
845  IMPLICIT NONE
846 
847  CHARACTER(LEN=*), INTENT(IN) :: gsi_file
848  CHARACTER(LEN=3), INTENT(IN) :: file_type
849  INTEGER, INTENT(IN), OPTIONAL :: lsoil
850 
851  INTEGER :: error, id_dim, ncid
852  INTEGER :: id_var, j
853 
854  INTEGER(KIND=1), ALLOCATABLE :: idummy(:,:)
855 
856  REAL(KIND=8), ALLOCATABLE :: dummy(:,:)
857 
858  CHARACTER(LEN=1) :: k_ch
859  CHARACTER(LEN=10) :: incvar
860  CHARACTER(LEN=80) :: err_msg
861  INTEGER :: k, i
862 
863  print*
864  print*, "READ INPUT GSI DATA FROM: "//trim(gsi_file)
865 
866  error=nf90_open(trim(gsi_file),nf90_nowrite,ncid)
867  CALL netcdf_err(error, 'OPENING FILE: '//trim(gsi_file) )
868 
869  error=nf90_inq_dimid(ncid, 'latitude', id_dim)
870  CALL netcdf_err(error, 'READING latitude ID' )
871  error=nf90_inquire_dimension(ncid,id_dim,len=jdim_gaus)
872  CALL netcdf_err(error, 'READING latitude length' )
873  jdim_gaus = jdim_gaus - 2 ! WILL IGNORE POLE POINTS
874 
875  error=nf90_inq_dimid(ncid, 'longitude', id_dim)
876  CALL netcdf_err(error, 'READING longitude ID' )
877  error=nf90_inquire_dimension(ncid,id_dim,len=idim_gaus)
878  CALL netcdf_err(error, 'READING longitude length' )
879 
880  IF (file_type=='NST') then
881  ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
882  ALLOCATE(dtref_gaus(idim_gaus,jdim_gaus))
883 
884  error=nf90_inq_varid(ncid, "dtf", id_var)
885  CALL netcdf_err(error, 'READING dtf ID' )
886  error=nf90_get_var(ncid, id_var, dummy)
887  CALL netcdf_err(error, 'READING dtf' )
888 
889  ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
890  ALLOCATE(slmsk_gaus(idim_gaus,jdim_gaus))
891 
892  error=nf90_inq_varid(ncid, "msk", id_var)
893  CALL netcdf_err(error, 'READING msk ID' )
894  error=nf90_get_var(ncid, id_var, idummy)
895  CALL netcdf_err(error, 'READING msk' )
896 
897  ! REMOVE POLE POINTS.
898 
899  DO j = 1, jdim_gaus
900  slmsk_gaus(:,j) = idummy(:,j+1)
901  dtref_gaus(:,j) = dummy(:,j+1)
902  ENDDO
903 
904  ELSEIF (file_type=='LND') then
905 
906  ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
907  ALLOCATE(stc_inc_gaus(lsoil,idim_gaus,jdim_gaus))
908 
909  ! read in soil temperature increments in each layer
910  DO k = 1, lsoil
911  WRITE(k_ch, '(I1)') k
912 
913  incvar = "soilt"//k_ch//"_inc"
914  error=nf90_inq_varid(ncid, incvar, id_var)
915  err_msg = "reading "//incvar//" ID"
916  CALL netcdf_err(error, trim(err_msg))
917  error=nf90_get_var(ncid, id_var, dummy)
918  err_msg = "reading "//incvar//" data"
919  CALL netcdf_err(error, err_msg)
920 
921  DO j = 1, jdim_gaus
922  stc_inc_gaus(k,:,j) = dummy(:,j+1)
923  ENDDO
924  ENDDO
925 
926  ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
927  ALLOCATE(soilsnow_gaus(idim_gaus,jdim_gaus))
928 
929  error=nf90_inq_varid(ncid, "soilsnow_mask", id_var)
930  CALL netcdf_err(error, 'READING soilsnow_mask ID' )
931  error=nf90_get_var(ncid, id_var, idummy)
932  CALL netcdf_err(error, 'READING soilsnow_mask' )
933 
934  ! REMOVE POLE POINTS.
935 
936  DO j = 1, jdim_gaus
937  soilsnow_gaus(:,j) = idummy(:,j+1)
938  ENDDO
939 
940 
941  ELSE
942  print *, 'WARNING: FILE_TYPE', file_type, 'not recognised.', &
943  ', no increments read in'
944  ENDIF
945 
946  IF(ALLOCATED(dummy)) DEALLOCATE(dummy)
947  IF(ALLOCATED(idummy)) DEALLOCATE(idummy)
948 
949  error = nf90_close(ncid)
950 
951  END SUBROUTINE read_gsi_data
952 
1001  SUBROUTINE read_data(LSOIL,LENSFC,DO_NSST,INC_FILE,IS_NOAHMP, &
1002  TSFFCS,SMCFCS,SWEFCS,STCFCS, &
1003  TG3FCS,ZORFCS, &
1004  CVFCS,CVBFCS,CVTFCS,ALBFCS, &
1005  VEGFCS,SLIFCS,CNPFCS,F10M, &
1006  VETFCS,SOTFCS,ALFFCS, &
1007  USTAR,FMM,FHH, &
1008  SIHFCS,SICFCS,SITFCS, &
1009  TPRCP,SRFLAG,SNDFCS, &
1010  VMNFCS,VMXFCS,SLCFCS, &
1011  SLPFCS,ABSFCS,T2M,Q2M,SLMASK, &
1012  ZSOIL,NSST)
1013  USE mpi
1014 
1015  IMPLICIT NONE
1016 
1017  INTEGER, INTENT(IN) :: lsoil, lensfc
1018  LOGICAL, INTENT(IN) :: do_nsst, inc_file
1019 
1020  LOGICAL, OPTIONAL, INTENT(OUT) :: is_noahmp
1021 
1022  REAL, OPTIONAL, INTENT(OUT) :: cvfcs(lensfc), cvbfcs(lensfc)
1023  REAL, OPTIONAL, INTENT(OUT) :: cvtfcs(lensfc), albfcs(lensfc,4)
1024  REAL, OPTIONAL, INTENT(OUT) :: slifcs(lensfc), cnpfcs(lensfc)
1025  REAL, OPTIONAL, INTENT(OUT) :: vegfcs(lensfc), f10m(lensfc)
1026  REAL, OPTIONAL, INTENT(OUT) :: vetfcs(lensfc), sotfcs(lensfc)
1027  REAL, OPTIONAL, INTENT(OUT) :: tsffcs(lensfc), swefcs(lensfc)
1028  REAL, OPTIONAL, INTENT(OUT) :: tg3fcs(lensfc), zorfcs(lensfc)
1029  REAL, OPTIONAL, INTENT(OUT) :: alffcs(lensfc,2), ustar(lensfc)
1030  REAL, OPTIONAL, INTENT(OUT) :: fmm(lensfc), fhh(lensfc)
1031  REAL, OPTIONAL, INTENT(OUT) :: sihfcs(lensfc), sicfcs(lensfc)
1032  REAL, OPTIONAL, INTENT(OUT) :: sitfcs(lensfc), tprcp(lensfc)
1033  REAL, OPTIONAL, INTENT(OUT) :: srflag(lensfc), sndfcs(lensfc)
1034  REAL, OPTIONAL, INTENT(OUT) :: vmnfcs(lensfc), vmxfcs(lensfc)
1035  REAL, OPTIONAL, INTENT(OUT) :: slpfcs(lensfc), absfcs(lensfc)
1036  REAL, OPTIONAL, INTENT(OUT) :: t2m(lensfc), q2m(lensfc), slmask(lensfc)
1037  REAL, OPTIONAL, INTENT(OUT) :: slcfcs(lensfc,lsoil)
1038  REAL, OPTIONAL, INTENT(OUT) :: smcfcs(lensfc,lsoil)
1039  REAL, OPTIONAL, INTENT(OUT) :: stcfcs(lensfc,lsoil)
1040  REAL(KIND=4), OPTIONAL, INTENT(OUT) :: zsoil(lsoil)
1041 
1042  TYPE(nsst_data), OPTIONAL :: nsst ! intent(out) will crash
1043  ! because subtypes are allocated in main.
1044 
1045  CHARACTER(LEN=50) :: fnbgsi
1046  CHARACTER(LEN=3) :: rankch
1047 
1048  INTEGER :: error, error2, ncid, myrank
1049  INTEGER :: idim, jdim, id_dim
1050  INTEGER :: id_var, ierr
1051 
1052  REAL(KIND=8), ALLOCATABLE :: dummy(:,:), dummy3d(:,:,:)
1053 
1054  CALL mpi_comm_rank(mpi_comm_world, myrank, error)
1055 
1056  WRITE(rankch, '(I3.3)') (myrank+1)
1057 
1058  IF (inc_file) THEN
1059  fnbgsi = "./xainc." // rankch
1060  ELSE
1061  fnbgsi = "./fnbgsi." // rankch
1062  ENDIF
1063 
1064  print*
1065  print*, "READ INPUT SFC DATA FROM: "//trim(fnbgsi)
1066 
1067  error=nf90_open(trim(fnbgsi),nf90_nowrite,ncid)
1068  CALL netcdf_err(error, 'OPENING FILE: '//trim(fnbgsi) )
1069 
1070  error=nf90_inq_dimid(ncid, 'xaxis_1', id_dim)
1071  CALL netcdf_err(error, 'READING xaxis_1' )
1072  error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1073  CALL netcdf_err(error, 'READING xaxis_1' )
1074 
1075  error=nf90_inq_dimid(ncid, 'yaxis_1', id_dim)
1076  CALL netcdf_err(error, 'READING yaxis_1' )
1077  error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1078  CALL netcdf_err(error, 'READING yaxis_1' )
1079 
1080  IF ((idim*jdim) /= lensfc) THEN
1081  print*,'FATAL ERROR: DIMENSIONS WRONG.'
1082  CALL mpi_abort(mpi_comm_world, 88, ierr)
1083  ENDIF
1084 
1085 ! Check for records that indicate the restart file is
1086 ! for the Noah-MP land surface model.
1087 
1088  IF(PRESENT(is_noahmp))THEN
1089  error=nf90_inq_varid(ncid, "canliqxy", id_var)
1090  error2=nf90_inq_varid(ncid, "tsnoxy", id_var)
1091  is_noahmp=.false.
1092  IF(error == 0 .AND. error2 == 0) THEN
1093  is_noahmp=.true.
1094  print*,"- WILL PROCESS FOR NOAH-MP LSM."
1095  ENDIF
1096  ENDIF
1097 
1098  ALLOCATE(dummy(idim,jdim))
1099 
1100  IF (PRESENT(tsffcs)) THEN
1101  error=nf90_inq_varid(ncid, "tsea", id_var)
1102  CALL netcdf_err(error, 'READING tsea ID' )
1103  error=nf90_get_var(ncid, id_var, dummy)
1104  CALL netcdf_err(error, 'READING tsea' )
1105  tsffcs = reshape(dummy, (/lensfc/))
1106  ENDIF
1107 
1108  IF (PRESENT(swefcs)) THEN
1109  error=nf90_inq_varid(ncid, "sheleg", id_var)
1110  CALL netcdf_err(error, 'READING sheleg ID' )
1111  error=nf90_get_var(ncid, id_var, dummy)
1112  CALL netcdf_err(error, 'READING sheleg' )
1113  swefcs = reshape(dummy, (/lensfc/))
1114  ENDIF
1115 
1116  IF (PRESENT(tg3fcs)) THEN
1117  error=nf90_inq_varid(ncid, "tg3", id_var)
1118  CALL netcdf_err(error, 'READING tg3 ID' )
1119  error=nf90_get_var(ncid, id_var, dummy)
1120  CALL netcdf_err(error, 'READING tg3' )
1121  tg3fcs = reshape(dummy, (/lensfc/))
1122  ENDIF
1123 
1124  IF (PRESENT(zorfcs)) THEN
1125  error=nf90_inq_varid(ncid, "zorl", id_var)
1126  CALL netcdf_err(error, 'READING zorl ID' )
1127  error=nf90_get_var(ncid, id_var, dummy)
1128  CALL netcdf_err(error, 'READING zorl' )
1129  zorfcs = reshape(dummy, (/lensfc/))
1130  ENDIF
1131 
1132  IF (PRESENT(albfcs)) THEN
1133 
1134  error=nf90_inq_varid(ncid, "alvsf", id_var)
1135  CALL netcdf_err(error, 'READING alvsf ID' )
1136  error=nf90_get_var(ncid, id_var, dummy)
1137  CALL netcdf_err(error, 'READING alvsf' )
1138  albfcs(:,1) = reshape(dummy, (/lensfc/))
1139 
1140  error=nf90_inq_varid(ncid, "alvwf", id_var)
1141  CALL netcdf_err(error, 'READING alvwf ID' )
1142  error=nf90_get_var(ncid, id_var, dummy)
1143  CALL netcdf_err(error, 'READING alvwf' )
1144  albfcs(:,2) = reshape(dummy, (/lensfc/))
1145 
1146  error=nf90_inq_varid(ncid, "alnsf", id_var)
1147  CALL netcdf_err(error, 'READING alnsf ID' )
1148  error=nf90_get_var(ncid, id_var, dummy)
1149  CALL netcdf_err(error, 'READING alnsf' )
1150  albfcs(:,3) = reshape(dummy, (/lensfc/))
1151 
1152  error=nf90_inq_varid(ncid, "alnwf", id_var)
1153  CALL netcdf_err(error, 'READING alnwf ID' )
1154  error=nf90_get_var(ncid, id_var, dummy)
1155  CALL netcdf_err(error, 'READING alnwf' )
1156  albfcs(:,4) = reshape(dummy, (/lensfc/))
1157 
1158  ENDIF
1159 
1160  IF (PRESENT(slifcs)) THEN
1161  error=nf90_inq_varid(ncid, "slmsk", id_var)
1162  CALL netcdf_err(error, 'READING slmsk ID' )
1163  error=nf90_get_var(ncid, id_var, dummy)
1164  CALL netcdf_err(error, 'READING slmsk' )
1165  slifcs = reshape(dummy, (/lensfc/))
1166  slmask = slifcs
1167  WHERE (slmask > 1.5) slmask=0.0 ! remove sea ice
1168  ENDIF
1169 
1170  IF (PRESENT(cnpfcs)) THEN
1171  error=nf90_inq_varid(ncid, "canopy", id_var)
1172  CALL netcdf_err(error, 'READING canopy ID' )
1173  error=nf90_get_var(ncid, id_var, dummy)
1174  CALL netcdf_err(error, 'READING canopy' )
1175  cnpfcs = reshape(dummy, (/lensfc/))
1176  ENDIF
1177 
1178  IF (PRESENT(vegfcs)) THEN
1179  error=nf90_inq_varid(ncid, "vfrac", id_var)
1180  CALL netcdf_err(error, 'READING vfrac ID' )
1181  error=nf90_get_var(ncid, id_var, dummy)
1182  CALL netcdf_err(error, 'READING vfrac' )
1183  vegfcs = reshape(dummy, (/lensfc/))
1184  ENDIF
1185 
1186  IF (PRESENT(f10m)) THEN
1187  error=nf90_inq_varid(ncid, "f10m", id_var)
1188  CALL netcdf_err(error, 'READING f10m ID' )
1189  error=nf90_get_var(ncid, id_var, dummy)
1190  CALL netcdf_err(error, 'READING f10m' )
1191  f10m = reshape(dummy, (/lensfc/))
1192  ENDIF
1193 
1194  IF (PRESENT(vetfcs)) THEN
1195  error=nf90_inq_varid(ncid, "vtype", id_var)
1196  CALL netcdf_err(error, 'READING vtype ID' )
1197  error=nf90_get_var(ncid, id_var, dummy)
1198  CALL netcdf_err(error, 'READING vtype' )
1199  vetfcs = reshape(dummy, (/lensfc/))
1200  ENDIF
1201 
1202  IF (PRESENT(sotfcs)) THEN
1203  error=nf90_inq_varid(ncid, "stype", id_var)
1204  CALL netcdf_err(error, 'READING stype ID' )
1205  error=nf90_get_var(ncid, id_var, dummy)
1206  CALL netcdf_err(error, 'READING stype' )
1207  sotfcs = reshape(dummy, (/lensfc/))
1208  ENDIF
1209 
1210  IF (PRESENT(alffcs)) THEN
1211  error=nf90_inq_varid(ncid, "facsf", id_var)
1212  CALL netcdf_err(error, 'READING facsf ID' )
1213  error=nf90_get_var(ncid, id_var, dummy)
1214  CALL netcdf_err(error, 'READING facsf' )
1215  alffcs(:,1) = reshape(dummy, (/lensfc/))
1216 
1217  error=nf90_inq_varid(ncid, "facwf", id_var)
1218  CALL netcdf_err(error, 'READING facwf ID' )
1219  error=nf90_get_var(ncid, id_var, dummy)
1220  CALL netcdf_err(error, 'READING facwf' )
1221  alffcs(:,2) = reshape(dummy, (/lensfc/))
1222  ENDIF
1223 
1224  IF (PRESENT(ustar)) THEN
1225  error=nf90_inq_varid(ncid, "uustar", id_var)
1226  CALL netcdf_err(error, 'READING uustar ID' )
1227  error=nf90_get_var(ncid, id_var, dummy)
1228  CALL netcdf_err(error, 'READING uustar' )
1229  ustar = reshape(dummy, (/lensfc/))
1230  ENDIF
1231 
1232  IF (PRESENT(fmm)) THEN
1233  error=nf90_inq_varid(ncid, "ffmm", id_var)
1234  CALL netcdf_err(error, 'READING ffmm ID' )
1235  error=nf90_get_var(ncid, id_var, dummy)
1236  CALL netcdf_err(error, 'READING ffmm' )
1237  fmm = reshape(dummy, (/lensfc/))
1238  ENDIF
1239 
1240  IF (PRESENT(fhh)) THEN
1241  error=nf90_inq_varid(ncid, "ffhh", id_var)
1242  CALL netcdf_err(error, 'READING ffhh ID' )
1243  error=nf90_get_var(ncid, id_var, dummy)
1244  CALL netcdf_err(error, 'READING ffhh' )
1245  fhh = reshape(dummy, (/lensfc/))
1246  ENDIF
1247 
1248  IF (PRESENT(sihfcs)) THEN
1249  error=nf90_inq_varid(ncid, "hice", id_var)
1250  CALL netcdf_err(error, 'READING hice ID' )
1251  error=nf90_get_var(ncid, id_var, dummy)
1252  CALL netcdf_err(error, 'READING hice' )
1253  sihfcs = reshape(dummy, (/lensfc/))
1254  ENDIF
1255 
1256  IF (PRESENT(sicfcs)) THEN
1257  error=nf90_inq_varid(ncid, "fice", id_var)
1258  CALL netcdf_err(error, 'READING fice ID' )
1259  error=nf90_get_var(ncid, id_var, dummy)
1260  CALL netcdf_err(error, 'READING fice' )
1261  sicfcs = reshape(dummy, (/lensfc/))
1262  ENDIF
1263 
1264  IF (PRESENT(sitfcs)) THEN
1265  error=nf90_inq_varid(ncid, "tisfc", id_var)
1266  CALL netcdf_err(error, 'READING tisfc ID' )
1267  error=nf90_get_var(ncid, id_var, dummy)
1268  CALL netcdf_err(error, 'READING tisfc' )
1269  sitfcs = reshape(dummy, (/lensfc/))
1270  ENDIF
1271 
1272  IF (PRESENT(tprcp)) THEN
1273  error=nf90_inq_varid(ncid, "tprcp", id_var)
1274  CALL netcdf_err(error, 'READING tprcp ID' )
1275  error=nf90_get_var(ncid, id_var, dummy)
1276  CALL netcdf_err(error, 'READING tprcp' )
1277  tprcp = reshape(dummy, (/lensfc/))
1278  ENDIF
1279 
1280  IF (PRESENT(srflag)) THEN
1281  error=nf90_inq_varid(ncid, "srflag", id_var)
1282  CALL netcdf_err(error, 'READING srflag ID' )
1283  error=nf90_get_var(ncid, id_var, dummy)
1284  CALL netcdf_err(error, 'READING srflag' )
1285  srflag = reshape(dummy, (/lensfc/))
1286  ENDIF
1287 
1288  IF (PRESENT(sndfcs)) THEN
1289  error=nf90_inq_varid(ncid, "snwdph", id_var)
1290  CALL netcdf_err(error, 'READING snwdph ID' )
1291  error=nf90_get_var(ncid, id_var, dummy)
1292  CALL netcdf_err(error, 'READING snwdph' )
1293  sndfcs = reshape(dummy, (/lensfc/))
1294  ENDIF
1295 
1296  IF (PRESENT(vmnfcs)) THEN
1297  error=nf90_inq_varid(ncid, "shdmin", id_var)
1298  CALL netcdf_err(error, 'READING shdmin ID' )
1299  error=nf90_get_var(ncid, id_var, dummy)
1300  CALL netcdf_err(error, 'READING shdmin' )
1301  vmnfcs = reshape(dummy, (/lensfc/))
1302  ENDIF
1303 
1304  IF (PRESENT(vmxfcs)) THEN
1305  error=nf90_inq_varid(ncid, "shdmax", id_var)
1306  CALL netcdf_err(error, 'READING shdmax ID' )
1307  error=nf90_get_var(ncid, id_var, dummy)
1308  CALL netcdf_err(error, 'READING shdmax' )
1309  vmxfcs = reshape(dummy, (/lensfc/))
1310  ENDIF
1311 
1312  IF (PRESENT(slpfcs)) THEN
1313  error=nf90_inq_varid(ncid, "slope", id_var)
1314  CALL netcdf_err(error, 'READING slope ID' )
1315  error=nf90_get_var(ncid, id_var, dummy)
1316  CALL netcdf_err(error, 'READING slope' )
1317  slpfcs = reshape(dummy, (/lensfc/))
1318  ENDIF
1319 
1320  IF (PRESENT(absfcs)) THEN
1321  error=nf90_inq_varid(ncid, "snoalb", id_var)
1322  CALL netcdf_err(error, 'READING snoalb ID' )
1323  error=nf90_get_var(ncid, id_var, dummy)
1324  CALL netcdf_err(error, 'READING snoalb' )
1325  absfcs = reshape(dummy, (/lensfc/))
1326  ENDIF
1327 
1328  IF (PRESENT(t2m)) THEN
1329  error=nf90_inq_varid(ncid, "t2m", id_var)
1330  CALL netcdf_err(error, 'READING t2m ID' )
1331  error=nf90_get_var(ncid, id_var, dummy)
1332  CALL netcdf_err(error, 'READING t2m' )
1333  t2m = reshape(dummy, (/lensfc/))
1334  ENDIF
1335 
1336  IF (PRESENT(q2m)) THEN
1337  error=nf90_inq_varid(ncid, "q2m", id_var)
1338  CALL netcdf_err(error, 'READING q2m ID' )
1339  error=nf90_get_var(ncid, id_var, dummy)
1340  CALL netcdf_err(error, 'READING q2m' )
1341  q2m = reshape(dummy, (/lensfc/))
1342  ENDIF
1343 
1344  nsst_read : IF(do_nsst) THEN
1345 
1346  print*
1347  print*,"WILL READ NSST RECORDS."
1348 
1349  error=nf90_inq_varid(ncid, "c_0", id_var)
1350  CALL netcdf_err(error, 'READING c_0 ID' )
1351  error=nf90_get_var(ncid, id_var, dummy)
1352  CALL netcdf_err(error, 'READING c_0' )
1353  nsst%C_0 = reshape(dummy, (/lensfc/))
1354 
1355  error=nf90_inq_varid(ncid, "c_d", id_var)
1356  CALL netcdf_err(error, 'READING c_d ID' )
1357  error=nf90_get_var(ncid, id_var, dummy)
1358  CALL netcdf_err(error, 'READING c_d' )
1359  nsst%C_D = reshape(dummy, (/lensfc/))
1360 
1361  error=nf90_inq_varid(ncid, "d_conv", id_var)
1362  CALL netcdf_err(error, 'READING d_conv ID' )
1363  error=nf90_get_var(ncid, id_var, dummy)
1364  CALL netcdf_err(error, 'READING d_conv' )
1365  nsst%D_CONV = reshape(dummy, (/lensfc/))
1366 
1367  error=nf90_inq_varid(ncid, "dt_cool", id_var)
1368  CALL netcdf_err(error, 'READING dt_cool ID' )
1369  error=nf90_get_var(ncid, id_var, dummy)
1370  CALL netcdf_err(error, 'READING dt_cool' )
1371  nsst%DT_COOL = reshape(dummy, (/lensfc/))
1372 
1373  error=nf90_inq_varid(ncid, "ifd", id_var)
1374  CALL netcdf_err(error, 'READING ifd ID' )
1375  error=nf90_get_var(ncid, id_var, dummy)
1376  CALL netcdf_err(error, 'READING ifd' )
1377  nsst%IFD = reshape(dummy, (/lensfc/))
1378 
1379  error=nf90_inq_varid(ncid, "qrain", id_var)
1380  CALL netcdf_err(error, 'READING qrain ID' )
1381  error=nf90_get_var(ncid, id_var, dummy)
1382  CALL netcdf_err(error, 'READING qrain' )
1383  nsst%QRAIN = reshape(dummy, (/lensfc/))
1384 
1385  error=nf90_inq_varid(ncid, "tref", id_var)
1386  CALL netcdf_err(error, 'READING tref ID' )
1387  error=nf90_get_var(ncid, id_var, dummy)
1388  CALL netcdf_err(error, 'READING tref' )
1389  nsst%TREF = reshape(dummy, (/lensfc/))
1390 
1391  error=nf90_inq_varid(ncid, "w_0", id_var)
1392  CALL netcdf_err(error, 'READING w_0 ID' )
1393  error=nf90_get_var(ncid, id_var, dummy)
1394  CALL netcdf_err(error, 'READING w_0' )
1395  nsst%W_0 = reshape(dummy, (/lensfc/))
1396 
1397  error=nf90_inq_varid(ncid, "w_d", id_var)
1398  CALL netcdf_err(error, 'READING w_d ID' )
1399  error=nf90_get_var(ncid, id_var, dummy)
1400  CALL netcdf_err(error, 'READING w_d' )
1401  nsst%W_D = reshape(dummy, (/lensfc/))
1402 
1403  error=nf90_inq_varid(ncid, "xs", id_var)
1404  CALL netcdf_err(error, 'READING xs ID' )
1405  error=nf90_get_var(ncid, id_var, dummy)
1406  CALL netcdf_err(error, 'READING xs' )
1407  nsst%XS = reshape(dummy, (/lensfc/))
1408 
1409  error=nf90_inq_varid(ncid, "xt", id_var)
1410  CALL netcdf_err(error, 'READING xt ID' )
1411  error=nf90_get_var(ncid, id_var, dummy)
1412  CALL netcdf_err(error, 'READING xt' )
1413  nsst%XT = reshape(dummy, (/lensfc/))
1414 
1415  error=nf90_inq_varid(ncid, "xtts", id_var)
1416  CALL netcdf_err(error, 'READING xtts ID' )
1417  error=nf90_get_var(ncid, id_var, dummy)
1418  CALL netcdf_err(error, 'READING xtts' )
1419  nsst%XTTS = reshape(dummy, (/lensfc/))
1420 
1421  error=nf90_inq_varid(ncid, "xu", id_var)
1422  CALL netcdf_err(error, 'READING xu ID' )
1423  error=nf90_get_var(ncid, id_var, dummy)
1424  CALL netcdf_err(error, 'READING xu' )
1425  nsst%XU = reshape(dummy, (/lensfc/))
1426 
1427  error=nf90_inq_varid(ncid, "xv", id_var)
1428  CALL netcdf_err(error, 'READING xv ID' )
1429  error=nf90_get_var(ncid, id_var, dummy)
1430  CALL netcdf_err(error, 'READING xv' )
1431  nsst%XV = reshape(dummy, (/lensfc/))
1432 
1433  error=nf90_inq_varid(ncid, "xz", id_var)
1434  CALL netcdf_err(error, 'READING xz ID' )
1435  error=nf90_get_var(ncid, id_var, dummy)
1436  CALL netcdf_err(error, 'READING xz' )
1437  nsst%XZ = reshape(dummy, (/lensfc/))
1438 
1439  error=nf90_inq_varid(ncid, "xzts", id_var)
1440  CALL netcdf_err(error, 'READING xzts ID' )
1441  error=nf90_get_var(ncid, id_var, dummy)
1442  CALL netcdf_err(error, 'READING xzts' )
1443  nsst%XZTS = reshape(dummy, (/lensfc/))
1444 
1445  error=nf90_inq_varid(ncid, "z_c", id_var)
1446  CALL netcdf_err(error, 'READING z_c ID' )
1447  error=nf90_get_var(ncid, id_var, dummy)
1448  CALL netcdf_err(error, 'READING z_c' )
1449  nsst%Z_C = reshape(dummy, (/lensfc/))
1450 
1451  error=nf90_inq_varid(ncid, "zm", id_var)
1452  CALL netcdf_err(error, 'READING zm ID' )
1453  error=nf90_get_var(ncid, id_var, dummy)
1454  CALL netcdf_err(error, 'READING zm' )
1455  nsst%ZM = reshape(dummy, (/lensfc/))
1456 
1457  END IF nsst_read
1458 
1459  DEALLOCATE(dummy)
1460 
1461  ALLOCATE(dummy3d(idim,jdim,lsoil))
1462 
1463  IF (PRESENT(smcfcs)) THEN
1464  error=nf90_inq_varid(ncid, "smc", id_var)
1465  CALL netcdf_err(error, 'READING smc ID' )
1466  error=nf90_get_var(ncid, id_var, dummy3d)
1467  CALL netcdf_err(error, 'READING smc' )
1468  smcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1469  ENDIF
1470 
1471  IF (PRESENT(slcfcs)) THEN
1472  error=nf90_inq_varid(ncid, "slc", id_var)
1473  CALL netcdf_err(error, 'READING slc ID' )
1474  error=nf90_get_var(ncid, id_var, dummy3d)
1475  CALL netcdf_err(error, 'READING slc' )
1476  slcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1477  ENDIF
1478 
1479  IF (PRESENT(stcfcs)) THEN
1480  error=nf90_inq_varid(ncid, "stc", id_var)
1481  CALL netcdf_err(error, 'READING stc ID' )
1482  error=nf90_get_var(ncid, id_var, dummy3d)
1483  CALL netcdf_err(error, 'READING stc' )
1484  stcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1485  ENDIF
1486 
1487  DEALLOCATE(dummy3d)
1488 
1489 ! cloud fields not in warm restart files. set to zero?
1490 
1491  IF (PRESENT(cvfcs)) cvfcs = 0.0
1492  IF (PRESENT(cvtfcs)) cvtfcs = 0.0
1493  IF (PRESENT(cvbfcs)) cvbfcs = 0.0
1494 
1495 ! soil layer thicknesses not in warm restart files. hardwire
1496 ! for now.
1497 
1498  IF (PRESENT(zsoil)) THEN
1499  zsoil(1) = -0.1
1500  zsoil(2) = -0.4
1501  zsoil(3) = -1.0
1502  zsoil(4) = -2.0
1503  ENDIF
1504 
1505  error = nf90_close(ncid)
1506 
1507  END SUBROUTINE read_data
1508 
1525 subroutine read_tf_clim_grb(file_sst,sst,rlats_sst,rlons_sst,mlat_sst,mlon_sst,mon)
1527  use mpi
1528 
1529  implicit none
1530 
1531 ! declare passed variables and arrays
1532  character(*) , intent(in ) :: file_sst
1533  integer , intent(in ) :: mon,mlat_sst,mlon_sst
1534  real, dimension(mlat_sst) , intent( out) :: rlats_sst
1535  real, dimension(mlon_sst) , intent( out) :: rlons_sst
1536  real, dimension(mlon_sst,mlat_sst) , intent( out) :: sst
1537 
1538 ! declare local parameters
1539  integer,parameter:: lu_sst = 21 ! fortran unit number of grib sst file
1540  real, parameter :: deg2rad = 3.141593/180.0
1541 
1542 ! declare local variables and arrays
1543  logical(1), allocatable, dimension(:) :: lb
1544 
1545  integer :: nlat_sst
1546  integer :: nlon_sst
1547  integer :: iret,ni,nj
1548  integer :: mscan,kb1,ierr
1549  integer :: jincdir,i,iincdir,kb2,kb3,kf,kg,k,j,jf
1550  integer, dimension(22):: jgds,kgds
1551  integer, dimension(25):: jpds,kpds
1552 
1553  real :: xsst0
1554  real :: ysst0
1555  real :: dres
1556  real, allocatable, dimension(:) :: f
1557 
1558 !************+******************************************************************************
1559 !
1560 ! open sst analysis file (grib)
1561  write(*,*) ' sstclm : ',file_sst
1562  call baopenr(lu_sst,trim(file_sst),iret)
1563  if (iret /= 0 ) then
1564  write(6,*)'FATAL ERROR in read_tf_clm_grb: error opening sst file.'
1565  CALL mpi_abort(mpi_comm_world, 111, ierr)
1566  endif
1567 
1568 ! define sst variables for read
1569  j=-1
1570  jpds=-1
1571  jgds=-1
1572  jpds(5)=11 ! sst variable
1573  jpds(6)=1 ! surface
1574  jpds(9) = mon
1575  call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1576 
1577  nlat_sst = kgds(3) ! number points on longitude circle (360)
1578  nlon_sst = kgds(2) ! number points on latitude circle (720)
1579 
1580 ! write(*,*) 'nlat_sst, nlon_sst, mon : ',nlat_sst, nlon_sst, mon
1581 ! write(*,*) 'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1582 
1583 ! allocate arrays
1584  allocate(lb(nlat_sst*nlon_sst))
1585  allocate(f(nlat_sst*nlon_sst))
1586  jf=nlat_sst*nlon_sst
1587 
1588 ! read in the analysis
1589  call getgb(lu_sst,0,jf,j,jpds,jgds,kf,k,kpds,kgds,lb,f,iret)
1590  if (iret /= 0) then
1591  write(6,*)'FATAL ERROR in read_tf_clm_grb: error reading sst analysis data record.'
1592  deallocate(lb,f)
1593  CALL mpi_abort(mpi_comm_world, 111, ierr)
1594  endif
1595 
1596  if ( (nlat_sst /= mlat_sst) .or. (nlon_sst /= mlon_sst) ) then
1597  write(6,*)'FATAL ERROR in read_rtg_org: inconsistent dimensions. mlat_sst,mlon_sst=',&
1598  mlat_sst,mlon_sst,' -versus- nlat_sst,nlon_sst=',nlat_sst,nlon_sst
1599  deallocate(lb,f)
1600  CALL mpi_abort(mpi_comm_world, 111, ierr)
1601  endif
1602 
1603 !
1604 ! get xlats and xlons
1605 !
1606  dres = 180.0/real(nlat_sst)
1607  ysst0 = 0.5*dres-90.0
1608  xsst0 = 0.5*dres
1609 
1610 ! get lat_sst & lon_sst
1611  do j = 1, nlat_sst
1612  rlats_sst(j) = ysst0 + real(j-1)*dres
1613  enddo
1614 
1615  do i = 1, nlon_sst
1616  rlons_sst(i) = (xsst0 + real(i-1)*dres)
1617  enddo
1618 
1619 ! load dimensions and grid specs. check for unusual values
1620  ni=kgds(2) ! 720 for 0.5 x 0.5
1621  nj=kgds(3) ! 360 for 0.5 x 0.5 resolution
1622 
1623  mscan=kgds(11)
1624  kb1=ibits(mscan,7,1) ! i scan direction
1625  kb2=ibits(mscan,6,1) ! j scan direction
1626  kb3=ibits(mscan,5,1) ! (i,j) or (j,i)
1627 
1628 ! get i and j scanning directions from kb1 and kb2.
1629 ! 0 yields +1, 1 yields -1. +1 is west to east, -1 is east to west.
1630  iincdir = 1-kb1*2
1631 
1632 ! 0 yields -1, 1 yields +1. +1 is south to north, -1 is north to south.
1633  jincdir = kb2*2 - 1
1634 
1635 ! write(*,*) 'read_tf_clim_grb,iincdir,jincdir : ',iincdir,jincdir
1636  do k=1,kf
1637 
1638 ! kb3 from scan mode indicates if i points are consecutive
1639 ! or if j points are consecutive
1640  if(kb3==0)then ! (i,j)
1641  i=(ni+1)*kb1+(mod(k-1,ni)+1)*iincdir
1642  j=(nj+1)*(1-kb2)+jincdir*((k-1)/ni+1)
1643  else ! (j,i)
1644  j=(nj+1)*(1-kb2)+(mod(k-1,nj)+1)*jincdir
1645  i=(ni+1)*kb1+iincdir*((k-1)/nj+1)
1646  endif
1647  sst(i,j) = f(k)
1648  end do
1649 
1650  deallocate(lb,f)
1651 
1652  call baclose(lu_sst,iret)
1653  if (iret /= 0 ) then
1654  write(6,*)'FATAL ERROR in read_tf_clm_grb: error closing sst file.'
1655  CALL mpi_abort(mpi_comm_world, 121, ierr)
1656  endif
1657 
1658 end subroutine read_tf_clim_grb
1659 
1667 subroutine get_tf_clm_dim(file_sst,mlat_sst,mlon_sst)
1668  use mpi
1669 
1670  implicit none
1671 
1672 ! declare passed variables and arrays
1673  character(*) , intent(in ) :: file_sst
1674  integer , intent(out) :: mlat_sst,mlon_sst
1675 
1676 ! declare local parameters
1677  integer,parameter:: lu_sst = 21 ! fortran unit number of grib sst file
1678 
1679  integer :: iret
1680  integer :: mscan,kb1
1681  integer :: kf,kg,k,j,ierr
1682  integer, dimension(22):: jgds,kgds
1683  integer, dimension(25):: jpds,kpds
1684 
1685 !************+******************************************************************************
1686 !
1687 ! open sst analysis file (grib)
1688  call baopenr(lu_sst,trim(file_sst),iret)
1689  if (iret /= 0 ) then
1690  write(6,*)'FATAL ERROR in get_tf_clm_dim: error opening sst file.'
1691  CALL mpi_abort(mpi_comm_world, 111, ierr)
1692  endif
1693 
1694 ! define sst variables for read
1695  j=-1
1696  jpds=-1
1697  jgds=-1
1698  jpds(5)=11 ! sst variable
1699  jpds(6)=1 ! surface
1700  jpds(9) = 1
1701  call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1702 
1703  mlat_sst = kgds(3) ! number points on longitude circle (360)
1704  mlon_sst = kgds(2) ! number points on latitude circle (720)
1705 
1706  write(*,*) 'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1707 
1708  call baclose(lu_sst,iret)
1709  if (iret /= 0 ) then
1710  write(6,*)'FATAL ERROR in get_tf_clm_dim: error closing sst file.'
1711  CALL mpi_abort(mpi_comm_world, 121, ierr)
1712  endif
1713 end subroutine get_tf_clm_dim
1714 
1726 subroutine read_salclm_gfs_nc(filename,sal,xlats,xlons,nlat,nlon,itime)
1727  use netcdf
1728  implicit none
1729 
1730  ! This is the name of the data file we will read.
1731  character (len=*), intent(in) :: filename
1732  integer, intent(in) :: nlat,nlon
1733  integer, intent(in) :: itime
1734  real, dimension(nlat), intent(out) :: xlats
1735  real, dimension(nlon), intent(out) :: xlons
1736  real, dimension(nlon,nlat), intent(out) :: sal
1737 ! Local variables
1738  integer :: ncid,ntime
1739 
1740  integer, parameter :: ndims = 3
1741  character (len = *), parameter :: lat_name = "latitude"
1742  character (len = *), parameter :: lon_name = "longitude"
1743  character (len = *), parameter :: t_name = "time"
1744  character (len = *), parameter :: sal_name="sal"
1745  integer :: no_fill,fill_value
1746  integer :: time_varid,lon_varid, lat_varid, z_varid, sal_varid
1747 
1748  ! The start and count arrays will tell the netCDF library where to read our data.
1749  integer, dimension(ndims) :: start, count
1750 
1751  character (len = *), parameter :: units = "units"
1752  character (len = *), parameter :: sal_units = "psu"
1753  ! PSU (Practical SalinitUnit). 1 PSU = 1g/kg
1754  character (len = *), parameter :: time_units = "months"
1755  character (len = *), parameter :: lat_units = "degrees_north"
1756  character (len = *), parameter :: lon_units = "degrees_east"
1757 
1758  integer :: missv
1759 ! Loop indices
1760  integer :: i,j
1761 
1762 ! Open the file.
1763  call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1764 
1765 ! Get the varids of time, latitude, longitude & depth coordinate variables.
1766  call nc_check( nf90_inq_varid(ncid, t_name, time_varid) )
1767  call nc_check( nf90_inq_varid(ncid, lat_name, lat_varid) )
1768  call nc_check( nf90_inq_varid(ncid, lon_name, lon_varid) )
1769 
1770 ! Read the time, latitude and longitude data.
1771 ! call nc_check( nf90_get_var(ncid, time_varid, ntime) )
1772  call nc_check( nf90_get_var(ncid, lat_varid, xlats) )
1773  call nc_check( nf90_get_var(ncid, lon_varid, xlons) )
1774 
1775 ! Get the varids of the sal netCDF variables.
1776  call nc_check( nf90_inq_varid(ncid, sal_name,sal_varid) )
1777 
1778 ! Read 1 record of nlat*nlon values, starting at the beginning
1779 ! of the record (the (1, 1, 1, rec) element in the netCDF file).
1780  start = (/ 1, 1, itime /)
1781  count = (/ nlon, nlat, 1 /)
1782 
1783 ! write(*,*) 'read_salclm_gfs_nc itime : ',itime
1784 ! Read the sal data from the file, one record at a time.
1785  call nc_check( nf90_get_var(ncid, sal_varid, sal, start, count) )
1786 
1787 ! Close the file. This frees up any internal netCDF resources
1788 ! associated with the file.
1789  call nc_check( nf90_close(ncid) )
1790 
1791 ! If we got this far, everything worked as expected. Yipee!
1792 ! print *,"*** SUCCESS reading file ", filename, "!"
1793 
1794 end subroutine read_salclm_gfs_nc
1795 
1802 subroutine get_dim_nc(filename,nlat,nlon)
1803  use netcdf
1804  implicit none
1805 
1806  character (len=*), intent(in) :: filename
1807  integer, intent(out) :: nlat,nlon
1808 ! Local variables
1809  character (len = *), parameter :: lat_name = "latitude"
1810  character (len = *), parameter :: lon_name = "longitude"
1811  integer :: ncid
1812  integer :: latdimid,londimid
1813 
1814 ! Open the file.
1815  call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1816 
1817 ! Get dimensions
1818  call nc_check( nf90_inq_dimid(ncid,lat_name,latdimid) )
1819  call nc_check( nf90_inq_dimid(ncid,lon_name,londimid) )
1820  call nc_check( nf90_inquire_dimension(ncid,latdimid,len=nlat) )
1821  call nc_check( nf90_inquire_dimension(ncid,londimid,len=nlon) )
1822 
1823 ! write(*,'(a,1x,a6,2I8)') 'get_dim_nc, file, nlat, nlon : ',filename,nlat,nlon
1824 
1825 ! Close the file. This frees up any internal netCDF resources
1826 ! associated with the file.
1827  call nc_check( nf90_close(ncid) )
1828 
1829 ! If we got this far, everything worked as expected. Yipee!
1830 ! print *,"*** SUCCESS get dimensions from nc file ", filename, "!"
1831 
1832 end subroutine get_dim_nc
1833 
1839 subroutine nc_check(status)
1841  use mpi
1842  use netcdf
1843 
1844  integer, intent ( in) :: status
1845  integer :: ierr
1846 
1847  if(status /= nf90_noerr) then
1848  print *, 'FATAL ERROR:'
1849  print *, trim(nf90_strerror(status))
1850  CALL mpi_abort(mpi_comm_world, 122, ierr)
1851  end if
1852 end subroutine nc_check
1853 
1854  END MODULE read_write_data