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