global_cycle  1.7.0
 All Data Structures Files Functions Variables
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 
56 
57  PUBLIC :: read_data
58  PUBLIC :: read_gsi_data
59  PUBLIC :: read_lat_lon_orog
60  PUBLIC :: write_data
63 
64  CONTAINS
65 
117  subroutine write_data(slifcs,tsffcs,swefcs,tg3fcs,zorfcs, &
118  albfcs,alffcs,vegfcs,cnpfcs,f10m, &
119  t2m,q2m,vetfcs,sotfcs,ustar,fmm,fhh, &
120  sicfcs,sihfcs,sitfcs,tprcp,srflag, &
121  swdfcs,vmnfcs,vmxfcs,slpfcs, &
122  absfcs,slcfcs,smcfcs,stcfcs,&
123  idim,jdim,lensfc,lsoil,do_nsst,nsst)
124 
125 
126  use mpi
127 
128  implicit none
129 
130  integer, intent(in) :: idim, jdim, lensfc, lsoil
131 
132  logical, intent(in) :: do_nsst
133 
134  real, intent(in) :: slifcs(lensfc), tsffcs(lensfc)
135  real, intent(in) :: swefcs(lensfc), tg3fcs(lensfc)
136  real, intent(in) :: vegfcs(lensfc), cnpfcs(lensfc)
137  real, intent(in) :: zorfcs(lensfc), albfcs(lensfc,4)
138  real, intent(in) :: f10m(lensfc), alffcs(lensfc,2)
139  real, intent(in) :: t2m(lensfc), q2m(lensfc)
140  real, intent(in) :: vetfcs(lensfc), sotfcs(lensfc)
141  real, intent(in) :: ustar(lensfc), fmm(lensfc)
142  real, intent(in) :: fhh(lensfc), sicfcs(lensfc)
143  real, intent(in) :: sihfcs(lensfc), sitfcs(lensfc)
144  real, intent(in) :: tprcp(lensfc), srflag(lensfc)
145  real, intent(in) :: swdfcs(lensfc), vmnfcs(lensfc)
146  real, intent(in) :: vmxfcs(lensfc), slpfcs(lensfc)
147  real, intent(in) :: absfcs(lensfc), slcfcs(lensfc,lsoil)
148  real, intent(in) :: smcfcs(lensfc,lsoil), stcfcs(lensfc,lsoil)
149 
150  type(nsst_data) :: nsst
151 
152  character(len=3) :: rankch
153  character(len=50) :: fnbgso
154 
155  integer :: fsize=65536, inital=0
156  integer :: header_buffer_val = 16384
157  integer :: dims_3d(3), dims_strt(3), dims_end(3)
158  integer :: dims_4d(4), dims4_strt(4), dims4_end(4)
159  integer :: error, i, ncid
160  integer :: dim_x, dim_y, dim_lsoil, dim_time
161  integer :: id_x, id_y, id_lsoil, id_time
162  integer :: id_slmsk, id_tsea, id_sheleg
163  integer :: id_alnwf, id_alvwf, id_alnsf, id_alvsf
164  integer :: id_tg3, id_zorl, id_facsf, id_facwf
165  integer :: id_vfrac, id_canopy, id_f10m, id_t2m
166  integer :: id_q2m, id_stype, id_vtype, id_uustar
167  integer :: id_ffmm, id_ffhh, id_fice, id_hice
168  integer :: id_tisfc, id_tprcp, id_srflag
169  integer :: id_snwdph, id_shdmin, id_shdmax
170  integer :: id_slope, id_snoalb, id_qrain
171  integer :: id_dt_cool, id_ifd, id_d_conv
172  integer :: id_xzts, id_xtts, id_zm, id_xz
173  integer :: id_xv, id_xu, id_xs, id_xt
174  integer :: id_w_d, id_w_0, id_c_d, id_tfinc
175  integer :: id_c_0, id_z_c, id_tref
176  integer :: id_stc, id_smc, id_slc
177  integer :: myrank
178 
179  real(kind=4) :: times
180  real(kind=4), allocatable :: lsoil_data(:), x_data(:), y_data(:)
181  real(kind=8), allocatable :: dum2d(:,:), dum3d(:,:,:)
182 
183  call mpi_comm_rank(mpi_comm_world, myrank, error)
184 
185  write(rankch, '(i3.3)') (myrank+1)
186 
187  fnbgso = "./fnbgso." // rankch
188 
189  print*
190  print*,"WRITE OUTPUT SFC DATA TO: ",trim(fnbgso)
191 
192 !--- open the file
193  error = nf90_create(fnbgso, ior(nf90_netcdf4,nf90_classic_model), ncid, initialsize=inital, chunksize=fsize)
194  call netcdf_err(error, 'CREATING FILE='//trim(fnbgso) )
195 
196 !--- define dimensions
197  error = nf90_def_dim(ncid, 'xaxis_1', idim, dim_x)
198  call netcdf_err(error, 'DEFINING XAXIS DIMENSION' )
199  error = nf90_def_dim(ncid, 'yaxis_1', jdim, dim_y)
200  call netcdf_err(error, 'DEFINING YAXIS DIMENSION' )
201  error = nf90_def_dim(ncid, 'zaxis_1', lsoil, dim_lsoil)
202  call netcdf_err(error, 'DEFINING ZAXIS DIMENSION' )
203  error = nf90_def_dim(ncid, 'Time', 1, dim_time)
204  call netcdf_err(error, 'DEFINING TIME DIMENSION' )
205 
206  !--- define fields
207  error = nf90_def_var(ncid, 'xaxis_1', nf90_float, dim_x, id_x)
208  call netcdf_err(error, 'DEFINING XAXIS_1 FIELD' )
209  error = nf90_put_att(ncid, id_x, "long_name", "xaxis_1")
210  call netcdf_err(error, 'DEFINING XAXIS_1 LONG NAME' )
211  error = nf90_put_att(ncid, id_x, "units", "none")
212  call netcdf_err(error, 'DEFINING XAXIS_1 UNITS' )
213  error = nf90_put_att(ncid, id_x, "cartesian_axis", "X")
214  call netcdf_err(error, 'WRITING XAXIS_1 FIELD' )
215 
216  error = nf90_def_var(ncid, 'yaxis_1', nf90_float, dim_y, id_y)
217  call netcdf_err(error, 'DEFINING YAXIS_1 FIELD' )
218  error = nf90_put_att(ncid, id_y, "long_name", "yaxis_1")
219  call netcdf_err(error, 'DEFINING YAXIS_1 LONG NAME' )
220  error = nf90_put_att(ncid, id_y, "units", "none")
221  call netcdf_err(error, 'DEFINING YAXIS_1 UNITS' )
222  error = nf90_put_att(ncid, id_y, "cartesian_axis", "Y")
223  call netcdf_err(error, 'WRITING YAXIS_1 FIELD' )
224 
225  error = nf90_def_var(ncid, 'zaxis_1', nf90_float, dim_lsoil, id_lsoil)
226  call netcdf_err(error, 'DEFINING ZAXIS_1 FIELD' )
227  error = nf90_put_att(ncid, id_lsoil, "long_name", "zaxis_1")
228  call netcdf_err(error, 'DEFINING ZAXIS_1 LONG NAME' )
229  error = nf90_put_att(ncid, id_lsoil, "units", "none")
230  call netcdf_err(error, 'DEFINING ZAXIS_1 UNITS' )
231  error = nf90_put_att(ncid, id_lsoil, "cartesian_axis", "Z")
232  call netcdf_err(error, 'WRITING ZAXIS_1 FIELD' )
233 
234  error = nf90_def_var(ncid, 'Time', nf90_float, dim_time, id_time)
235  call netcdf_err(error, 'DEFINING TIME FIELD' )
236  error = nf90_put_att(ncid, id_time, "long_name", "Time")
237  call netcdf_err(error, 'DEFINING TIME LONG NAME' )
238  error = nf90_put_att(ncid, id_time, "units", "time level")
239  call netcdf_err(error, 'DEFINING TIME UNITS' )
240  error = nf90_put_att(ncid, id_time, "cartesian_axis", "T")
241  call netcdf_err(error, 'WRITING TIME FIELD' )
242 
243  dims_3d(1) = dim_x
244  dims_3d(2) = dim_y
245  dims_3d(3) = dim_time
246 
247  error = nf90_def_var(ncid, 'slmsk', nf90_double, dims_3d, id_slmsk)
248  call netcdf_err(error, 'DEFINING SLMSK' )
249  error = nf90_put_att(ncid, id_slmsk, "long_name", "slmsk")
250  call netcdf_err(error, 'DEFINING SLMSK LONG NAME' )
251  error = nf90_put_att(ncid, id_slmsk, "units", "none")
252  call netcdf_err(error, 'DEFINING SLMSK UNITS' )
253 
254  error = nf90_def_var(ncid, 'tsea', nf90_double, dims_3d, id_tsea)
255  call netcdf_err(error, 'DEFINING TSEA' )
256  error = nf90_put_att(ncid, id_tsea, "long_name", "tsea")
257  call netcdf_err(error, 'DEFINING TSEA LONG NAME' )
258  error = nf90_put_att(ncid, id_tsea, "units", "none")
259  call netcdf_err(error, 'DEFINING TSEA UNITS' )
260 
261  error = nf90_def_var(ncid, 'sheleg', nf90_double, dims_3d, id_sheleg)
262  call netcdf_err(error, 'DEFINING SHELEG' )
263  error = nf90_put_att(ncid, id_sheleg, "long_name", "sheleg")
264  call netcdf_err(error, 'DEFINING SHELEG LONG NAME' )
265  error = nf90_put_att(ncid, id_sheleg, "units", "none")
266  call netcdf_err(error, 'DEFINING SHELEG UNITS' )
267 
268  error = nf90_def_var(ncid, 'tg3', nf90_double, dims_3d, id_tg3)
269  call netcdf_err(error, 'DEFINING TG3' )
270  error = nf90_put_att(ncid, id_tg3, "long_name", "tg3")
271  call netcdf_err(error, 'DEFINING TG3 LONG NAME' )
272  error = nf90_put_att(ncid, id_tg3, "units", "none")
273  call netcdf_err(error, 'DEFINING TG3 UNITS' )
274 
275  error = nf90_def_var(ncid, 'zorl', nf90_double, dims_3d, id_zorl)
276  call netcdf_err(error, 'DEFINING ZORL' )
277  error = nf90_put_att(ncid, id_zorl, "long_name", "zorl")
278  call netcdf_err(error, 'DEFINING ZORL LONG NAME' )
279  error = nf90_put_att(ncid, id_zorl, "units", "none")
280  call netcdf_err(error, 'DEFINING ZORL UNITS' )
281 
282  error = nf90_def_var(ncid, 'alvsf', nf90_double, dims_3d, id_alvsf)
283  call netcdf_err(error, 'DEFINING ALVSF' )
284  error = nf90_put_att(ncid, id_alvsf, "long_name", "alvsf")
285  call netcdf_err(error, 'DEFINING ALVSF LONG NAME' )
286  error = nf90_put_att(ncid, id_alvsf, "units", "none")
287  call netcdf_err(error, 'DEFINING ALVSF UNITS' )
288 
289  error = nf90_def_var(ncid, 'alvwf', nf90_double, dims_3d, id_alvwf)
290  call netcdf_err(error, 'DEFINING ALVWF' )
291  error = nf90_put_att(ncid, id_alvwf, "long_name", "alvwf")
292  call netcdf_err(error, 'DEFINING ALVWF LONG NAME' )
293  error = nf90_put_att(ncid, id_alvwf, "units", "none")
294  call netcdf_err(error, 'DEFINING ALVWF UNITS' )
295 
296  error = nf90_def_var(ncid, 'alnsf', nf90_double, dims_3d, id_alnsf)
297  call netcdf_err(error, 'DEFINING ALNSF' )
298  error = nf90_put_att(ncid, id_alnsf, "long_name", "alnsf")
299  call netcdf_err(error, 'DEFINING ALNSF LONG NAME' )
300  error = nf90_put_att(ncid, id_alnsf, "units", "none")
301  call netcdf_err(error, 'DEFINING ALNSF UNITS' )
302 
303  error = nf90_def_var(ncid, 'alnwf', nf90_double, dims_3d, id_alnwf)
304  call netcdf_err(error, 'DEFINING ALNWF' )
305  error = nf90_put_att(ncid, id_alnwf, "long_name", "alnwf")
306  call netcdf_err(error, 'DEFINING ALNWF LONG NAME' )
307  error = nf90_put_att(ncid, id_alnwf, "units", "none")
308  call netcdf_err(error, 'DEFINING ALNWF UNITS' )
309 
310  error = nf90_def_var(ncid, 'facsf', nf90_double, dims_3d, id_facsf)
311  call netcdf_err(error, 'DEFINING FACSF' )
312  error = nf90_put_att(ncid, id_facsf, "long_name", "facsf")
313  call netcdf_err(error, 'DEFINING FACSF LONG NAME' )
314  error = nf90_put_att(ncid, id_facsf, "units", "none")
315  call netcdf_err(error, 'DEFINING FACSF UNITS' )
316 
317  error = nf90_def_var(ncid, 'facwf', nf90_double, dims_3d, id_facwf)
318  call netcdf_err(error, 'DEFINING FACWF' )
319  error = nf90_put_att(ncid, id_facwf, "long_name", "facwf")
320  call netcdf_err(error, 'DEFINING FACWF LONG NAME' )
321  error = nf90_put_att(ncid, id_facwf, "units", "none")
322  call netcdf_err(error, 'DEFINING FACWF UNITS' )
323 
324  error = nf90_def_var(ncid, 'vfrac', nf90_double, dims_3d, id_vfrac)
325  call netcdf_err(error, 'DEFINING VFRAC' )
326  error = nf90_put_att(ncid, id_vfrac, "long_name", "vfrac")
327  call netcdf_err(error, 'DEFINING FACWF LONG NAME' )
328  error = nf90_put_att(ncid, id_vfrac, "units", "none")
329  call netcdf_err(error, 'DEFINING VFRAC UNITS' )
330 
331  error = nf90_def_var(ncid, 'canopy', nf90_double, dims_3d, id_canopy)
332  call netcdf_err(error, 'DEFINING CANOPY' )
333  error = nf90_put_att(ncid, id_canopy, "long_name", "canopy")
334  call netcdf_err(error, 'DEFINING CANOPY LONG NAME' )
335  error = nf90_put_att(ncid, id_canopy, "units", "none")
336  call netcdf_err(error, 'DEFINING CANOPY UNITS' )
337 
338  error = nf90_def_var(ncid, 'f10m', nf90_double, dims_3d, id_f10m)
339  call netcdf_err(error, 'DEFINING F10M' )
340  error = nf90_put_att(ncid, id_f10m, "long_name", "f10m")
341  call netcdf_err(error, 'DEFINING F10M LONG NAME' )
342  error = nf90_put_att(ncid, id_f10m, "units", "none")
343  call netcdf_err(error, 'DEFINING F10M UNITS' )
344 
345  error = nf90_def_var(ncid, 't2m', nf90_double, dims_3d, id_t2m)
346  call netcdf_err(error, 'DEFINING T2M' )
347  error = nf90_put_att(ncid, id_t2m, "long_name", "t2m")
348  call netcdf_err(error, 'DEFINING T2M LONG NAME' )
349  error = nf90_put_att(ncid, id_t2m, "units", "none")
350  call netcdf_err(error, 'DEFINING T2M UNITS' )
351 
352  error = nf90_def_var(ncid, 'q2m', nf90_double, dims_3d, id_q2m)
353  call netcdf_err(error, 'DEFINING Q2M' )
354  error = nf90_put_att(ncid, id_q2m, "long_name", "q2m")
355  call netcdf_err(error, 'DEFINING Q2M LONG NAME' )
356  error = nf90_put_att(ncid, id_q2m, "units", "none")
357  call netcdf_err(error, 'DEFINING Q2M UNITS' )
358 
359  error = nf90_def_var(ncid, 'vtype', nf90_double, dims_3d, id_vtype)
360  call netcdf_err(error, 'DEFINING VTYPE' )
361  error = nf90_put_att(ncid, id_vtype, "long_name", "vtype")
362  call netcdf_err(error, 'DEFINING VTYPE LONG NAME' )
363  error = nf90_put_att(ncid, id_vtype, "units", "none")
364  call netcdf_err(error, 'DEFINING VTYPE UNITS' )
365 
366  error = nf90_def_var(ncid, 'stype', nf90_double, dims_3d, id_stype)
367  call netcdf_err(error, 'DEFINING STYPE' )
368  error = nf90_put_att(ncid, id_stype, "long_name", "stype")
369  call netcdf_err(error, 'DEFINING STYPE LONG NAME' )
370  error = nf90_put_att(ncid, id_stype, "units", "none")
371  call netcdf_err(error, 'DEFINING STYPE UNITS' )
372 
373  error = nf90_def_var(ncid, 'uustar', nf90_double, dims_3d, id_uustar)
374  call netcdf_err(error, 'DEFINING UUSTAR' )
375  error = nf90_put_att(ncid, id_uustar, "long_name", "uustar")
376  call netcdf_err(error, 'DEFINING UUSTAR LONG NAME' )
377  error = nf90_put_att(ncid, id_uustar, "units", "none")
378  call netcdf_err(error, 'DEFINING UUSTAR UNITS' )
379 
380  error = nf90_def_var(ncid, 'ffmm', nf90_double, dims_3d, id_ffmm)
381  call netcdf_err(error, 'DEFINING FFMM' )
382  error = nf90_put_att(ncid, id_ffmm, "long_name", "ffmm")
383  call netcdf_err(error, 'DEFINING FFMM LONG NAME' )
384  error = nf90_put_att(ncid, id_ffmm, "units", "none")
385  call netcdf_err(error, 'DEFINING FFMM UNITS' )
386 
387  error = nf90_def_var(ncid, 'ffhh', nf90_double, dims_3d, id_ffhh)
388  call netcdf_err(error, 'DEFINING FFHH' )
389  error = nf90_put_att(ncid, id_ffhh, "long_name", "ffhh")
390  call netcdf_err(error, 'DEFINING FFHH LONG NAME' )
391  error = nf90_put_att(ncid, id_ffhh, "units", "none")
392  call netcdf_err(error, 'DEFINING FFHH UNITS' )
393 
394  error = nf90_def_var(ncid, 'hice', nf90_double, dims_3d, id_hice)
395  call netcdf_err(error, 'DEFINING HICE' )
396  error = nf90_put_att(ncid, id_hice, "long_name", "hice")
397  call netcdf_err(error, 'DEFINING HICE LONG NAME' )
398  error = nf90_put_att(ncid, id_hice, "units", "none")
399  call netcdf_err(error, 'DEFINING HICE UNITS' )
400 
401  error = nf90_def_var(ncid, 'fice', nf90_double, dims_3d, id_fice)
402  call netcdf_err(error, 'DEFINING FICE' )
403  error = nf90_put_att(ncid, id_fice, "long_name", "fice")
404  call netcdf_err(error, 'DEFINING FICE LONG NAME' )
405  error = nf90_put_att(ncid, id_fice, "units", "none")
406  call netcdf_err(error, 'DEFINING FICE UNITS' )
407 
408  error = nf90_def_var(ncid, 'tisfc', nf90_double, dims_3d, id_tisfc)
409  call netcdf_err(error, 'DEFINING TISFC' )
410  error = nf90_put_att(ncid, id_tisfc, "long_name", "tisfc")
411  call netcdf_err(error, 'DEFINING TISFC LONG NAME' )
412  error = nf90_put_att(ncid, id_tisfc, "units", "none")
413  call netcdf_err(error, 'DEFINING TISFC UNITS' )
414 
415  error = nf90_def_var(ncid, 'tprcp', nf90_double, dims_3d, id_tprcp)
416  call netcdf_err(error, 'DEFINING TPRCP' )
417  error = nf90_put_att(ncid, id_tprcp, "long_name", "tprcp")
418  call netcdf_err(error, 'DEFINING TPRCP LONG NAME' )
419  error = nf90_put_att(ncid, id_tprcp, "units", "none")
420  call netcdf_err(error, 'DEFINING TPRCP UNITS' )
421 
422  error = nf90_def_var(ncid, 'srflag', nf90_double, dims_3d, id_srflag)
423  call netcdf_err(error, 'DEFINING SRFLAG' )
424  error = nf90_put_att(ncid, id_srflag, "long_name", "srflag")
425  call netcdf_err(error, 'DEFINING SRFLAG LONG NAME' )
426  error = nf90_put_att(ncid, id_srflag, "units", "none")
427  call netcdf_err(error, 'DEFINING SRFLAG UNITS' )
428 
429  error = nf90_def_var(ncid, 'snwdph', nf90_double, dims_3d, id_snwdph)
430  call netcdf_err(error, 'DEFINING SNWDPH' )
431  error = nf90_put_att(ncid, id_snwdph, "long_name", "snwdph")
432  call netcdf_err(error, 'DEFINING SNWDPH LONG NAME' )
433  error = nf90_put_att(ncid, id_snwdph, "units", "none")
434  call netcdf_err(error, 'DEFINING SNWDPH UNITS' )
435 
436  error = nf90_def_var(ncid, 'shdmin', nf90_double, dims_3d, id_shdmin)
437  call netcdf_err(error, 'DEFINING SHDMIN' )
438  error = nf90_put_att(ncid, id_shdmin, "long_name", "shdmin")
439  call netcdf_err(error, 'DEFINING SHDMIN LONG NAME' )
440  error = nf90_put_att(ncid, id_shdmin, "units", "none")
441  call netcdf_err(error, 'DEFINING SHDMIN UNITS' )
442 
443  error = nf90_def_var(ncid, 'shdmax', nf90_double, dims_3d, id_shdmax)
444  call netcdf_err(error, 'DEFINING SHDMAX' )
445  error = nf90_put_att(ncid, id_shdmax, "long_name", "shdmax")
446  call netcdf_err(error, 'DEFINING SHDMAX LONG NAME' )
447  error = nf90_put_att(ncid, id_shdmax, "units", "none")
448  call netcdf_err(error, 'DEFINING SHDMAX UNITS' )
449 
450  error = nf90_def_var(ncid, 'slope', nf90_double, dims_3d, id_slope)
451  call netcdf_err(error, 'DEFINING SLOPE' )
452  error = nf90_put_att(ncid, id_slope, "long_name", "slope")
453  call netcdf_err(error, 'DEFINING SLOPE LONG NAME' )
454  error = nf90_put_att(ncid, id_slope, "units", "none")
455  call netcdf_err(error, 'DEFINING SLOPE UNITS' )
456 
457  error = nf90_def_var(ncid, 'snoalb', nf90_double, dims_3d, id_snoalb)
458  call netcdf_err(error, 'DEFINING SNOALB' )
459  error = nf90_put_att(ncid, id_snoalb, "long_name", "snoalb")
460  call netcdf_err(error, 'DEFINING SNOALB LONG NAME' )
461  error = nf90_put_att(ncid, id_snoalb, "units", "none")
462  call netcdf_err(error, 'DEFINING SNOALB UNITS' )
463 
464  nsst_header : if (do_nsst) then
465 
466  print*
467  print*,"WRITE NSST RECORDS."
468 
469  error = nf90_def_var(ncid, 'tref', nf90_double, dims_3d, id_tref)
470  call netcdf_err(error, 'DEFINING TREF' )
471  error = nf90_put_att(ncid, id_tref, "long_name", "tref")
472  call netcdf_err(error, 'DEFINING TREF LONG NAME' )
473  error = nf90_put_att(ncid, id_tref, "units", "none")
474  call netcdf_err(error, 'DEFINING TREF UNITS' )
475 
476  error = nf90_def_var(ncid, 'z_c', nf90_double, dims_3d, id_z_c)
477  call netcdf_err(error, 'DEFINING Z_C' )
478  error = nf90_put_att(ncid, id_z_c, "long_name", "z_c")
479  call netcdf_err(error, 'DEFINING Z_C LONG NAME' )
480  error = nf90_put_att(ncid, id_z_c, "units", "none")
481  call netcdf_err(error, 'DEFINING Z_C UNITS' )
482 
483  error = nf90_def_var(ncid, 'c_0', nf90_double, dims_3d, id_c_0)
484  call netcdf_err(error, 'DEFINING C_0' )
485  error = nf90_put_att(ncid, id_c_0, "long_name", "c_0")
486  call netcdf_err(error, 'DEFINING C_0 LONG NAME' )
487  error = nf90_put_att(ncid, id_c_0, "units", "none")
488  call netcdf_err(error, 'DEFINING C_0 UNITS' )
489 
490  error = nf90_def_var(ncid, 'c_d', nf90_double, dims_3d, id_c_d)
491  call netcdf_err(error, 'DEFINING C_D' )
492  error = nf90_put_att(ncid, id_c_d, "long_name", "c_d")
493  call netcdf_err(error, 'DEFINING C_D LONG NAME' )
494  error = nf90_put_att(ncid, id_c_d, "units", "none")
495  call netcdf_err(error, 'DEFINING C_D UNITS' )
496 
497  error = nf90_def_var(ncid, 'w_0', nf90_double, dims_3d, id_w_0)
498  call netcdf_err(error, 'DEFINING W_0' )
499  error = nf90_put_att(ncid, id_w_0, "long_name", "w_0")
500  call netcdf_err(error, 'DEFINING W_0 LONG NAME' )
501  error = nf90_put_att(ncid, id_w_0, "units", "none")
502  call netcdf_err(error, 'DEFINING W_0 UNITS' )
503 
504  error = nf90_def_var(ncid, 'w_d', nf90_double, dims_3d, id_w_d)
505  call netcdf_err(error, 'DEFINING W_D' )
506  error = nf90_put_att(ncid, id_w_d, "long_name", "w_d")
507  call netcdf_err(error, 'DEFINING W_D LONG NAME' )
508  error = nf90_put_att(ncid, id_w_d, "units", "none")
509  call netcdf_err(error, 'DEFINING W_D UNITS' )
510 
511  error = nf90_def_var(ncid, 'xt', nf90_double, dims_3d, id_xt)
512  call netcdf_err(error, 'DEFINING XT' )
513  error = nf90_put_att(ncid, id_xt, "long_name", "xt")
514  call netcdf_err(error, 'DEFINING XT LONG NAME' )
515  error = nf90_put_att(ncid, id_xt, "units", "none")
516  call netcdf_err(error, 'DEFINING XT UNITS' )
517 
518  error = nf90_def_var(ncid, 'xs', nf90_double, dims_3d, id_xs)
519  call netcdf_err(error, 'DEFINING XS' )
520  error = nf90_put_att(ncid, id_xs, "long_name", "xs")
521  call netcdf_err(error, 'DEFINING XS LONG NAME' )
522  error = nf90_put_att(ncid, id_xs, "units", "none")
523  call netcdf_err(error, 'DEFINING XS UNITS' )
524 
525  error = nf90_def_var(ncid, 'xu', nf90_double, dims_3d, id_xu)
526  call netcdf_err(error, 'DEFINING XU' )
527  error = nf90_put_att(ncid, id_xu, "long_name", "xu")
528  call netcdf_err(error, 'DEFINING XU LONG NAME' )
529  error = nf90_put_att(ncid, id_xu, "units", "none")
530  call netcdf_err(error, 'DEFINING XU UNITS' )
531 
532  error = nf90_def_var(ncid, 'xv', nf90_double, dims_3d, id_xv)
533  call netcdf_err(error, 'DEFINING XV' )
534  error = nf90_put_att(ncid, id_xv, "long_name", "xv")
535  call netcdf_err(error, 'DEFINING XV LONG NAME' )
536  error = nf90_put_att(ncid, id_xv, "units", "none")
537  call netcdf_err(error, 'DEFINING XV UNITS' )
538 
539  error = nf90_def_var(ncid, 'xz', nf90_double, dims_3d, id_xz)
540  call netcdf_err(error, 'DEFINING XZ' )
541  error = nf90_put_att(ncid, id_xz, "long_name", "xz")
542  call netcdf_err(error, 'DEFINING XZ LONG NAME' )
543  error = nf90_put_att(ncid, id_xz, "units", "none")
544  call netcdf_err(error, 'DEFINING XZ UNITS' )
545 
546  error = nf90_def_var(ncid, 'zm', nf90_double, dims_3d, id_zm)
547  call netcdf_err(error, 'DEFINING ZM' )
548  error = nf90_put_att(ncid, id_zm, "long_name", "zm")
549  call netcdf_err(error, 'DEFINING ZM LONG NAME' )
550  error = nf90_put_att(ncid, id_zm, "units", "none")
551  call netcdf_err(error, 'DEFINING ZM UNITS' )
552 
553  error = nf90_def_var(ncid, 'xtts', nf90_double, dims_3d, id_xtts)
554  call netcdf_err(error, 'DEFINING XTTS' )
555  error = nf90_put_att(ncid, id_xtts, "long_name", "xtts")
556  call netcdf_err(error, 'DEFINING XTTS LONG NAME' )
557  error = nf90_put_att(ncid, id_xtts, "units", "none")
558  call netcdf_err(error, 'DEFINING XTTS UNITS' )
559 
560  error = nf90_def_var(ncid, 'xzts', nf90_double, dims_3d, id_xzts)
561  call netcdf_err(error, 'DEFINING XZTS' )
562  error = nf90_put_att(ncid, id_xzts, "long_name", "xzts")
563  call netcdf_err(error, 'DEFINING XZTS LONG NAME' )
564  error = nf90_put_att(ncid, id_xzts, "units", "none")
565  call netcdf_err(error, 'DEFINING XZTS UNITS' )
566 
567  error = nf90_def_var(ncid, 'd_conv', nf90_double, dims_3d, id_d_conv)
568  call netcdf_err(error, 'DEFINING D_CONV' )
569  error = nf90_put_att(ncid, id_d_conv, "long_name", "d_conv")
570  call netcdf_err(error, 'DEFINING D_CONV LONG NAME' )
571  error = nf90_put_att(ncid, id_d_conv, "units", "none")
572  call netcdf_err(error, 'DEFINING D_CONV UNITS' )
573 
574  error = nf90_def_var(ncid, 'ifd', nf90_double, dims_3d, id_ifd)
575  call netcdf_err(error, 'DEFINING IFD' )
576  error = nf90_put_att(ncid, id_ifd, "long_name", "ifd")
577  call netcdf_err(error, 'DEFINING IFD LONG NAME' )
578  error = nf90_put_att(ncid, id_ifd, "units", "none")
579  call netcdf_err(error, 'DEFINING IFD UNITS' )
580 
581  error = nf90_def_var(ncid, 'dt_cool', nf90_double, dims_3d, id_dt_cool)
582  call netcdf_err(error, 'DEFINING DT_COOL' )
583  error = nf90_put_att(ncid, id_dt_cool, "long_name", "dt_cool")
584  call netcdf_err(error, 'DEFINING DT_COOL LONG NAME' )
585  error = nf90_put_att(ncid, id_dt_cool, "units", "none")
586  call netcdf_err(error, 'DEFINING DT_COOL UNITS' )
587 
588  error = nf90_def_var(ncid, 'qrain', nf90_double, dims_3d, id_qrain)
589  call netcdf_err(error, 'DEFINING QRAIN' )
590  error = nf90_put_att(ncid, id_qrain, "long_name", "qrain")
591  call netcdf_err(error, 'DEFINING QRAIN LONG NAME' )
592  error = nf90_put_att(ncid, id_qrain, "units", "none")
593  call netcdf_err(error, 'DEFINING QRAIN UNITS' )
594 
595  error = nf90_def_var(ncid, 'tfinc', nf90_double, dims_3d, id_tfinc)
596  call netcdf_err(error, 'DEFINING TFINC' )
597  error = nf90_put_att(ncid, id_tfinc, "long_name", "tfinc")
598  call netcdf_err(error, 'DEFINING TFINC LONG NAME' )
599  error = nf90_put_att(ncid, id_tfinc, "units", "none")
600  call netcdf_err(error, 'DEFINING TFINC UNITS' )
601 
602  endif nsst_header
603 
604  dims_4d(1) = dim_x
605  dims_4d(2) = dim_y
606  dims_4d(3) = dim_lsoil
607  dims_4d(4) = dim_time
608 
609  error = nf90_def_var(ncid, 'stc', nf90_double, dims_4d, id_stc)
610  call netcdf_err(error, 'DEFINING STC' )
611  error = nf90_put_att(ncid, id_stc, "long_name", "stc")
612  call netcdf_err(error, 'DEFINING STC LONG NAME' )
613  error = nf90_put_att(ncid, id_stc, "units", "none")
614  call netcdf_err(error, 'DEFINING STC UNITS' )
615 
616  error = nf90_def_var(ncid, 'smc', nf90_double, dims_4d, id_smc)
617  call netcdf_err(error, 'DEFINING SMC' )
618  error = nf90_put_att(ncid, id_smc, "long_name", "smc")
619  call netcdf_err(error, 'DEFINING SMC LONG NAME' )
620  error = nf90_put_att(ncid, id_smc, "units", "none")
621  call netcdf_err(error, 'DEFINING SMC UNITS' )
622 
623  error = nf90_def_var(ncid, 'slc', nf90_double, dims_4d, id_slc)
624  call netcdf_err(error, 'DEFINING SLC' )
625  error = nf90_put_att(ncid, id_slc, "long_name", "slc")
626  call netcdf_err(error, 'DEFINING SLC LONG NAME' )
627  error = nf90_put_att(ncid, id_slc, "units", "none")
628  call netcdf_err(error, 'DEFINING SLC UNITS' )
629 
630  error = nf90_enddef(ncid, header_buffer_val,4,0,4)
631  call netcdf_err(error, 'DEFINING HEADER' )
632 
633 !---------------------------------------------------------------------------------
634 ! Write data
635 !---------------------------------------------------------------------------------
636 
637  allocate(lsoil_data(lsoil))
638  do i = 1, lsoil
639  lsoil_data(i) = float(i)
640  enddo
641 
642  allocate(x_data(idim))
643  do i = 1, idim
644  x_data(i) = float(i)
645  enddo
646 
647  allocate(y_data(jdim))
648  do i = 1, jdim
649  y_data(i) = float(i)
650  enddo
651 
652  error = nf90_put_var( ncid, id_lsoil, lsoil_data)
653  call netcdf_err(error, 'WRITING ZAXIS RECORD' )
654  error = nf90_put_var( ncid, id_x, x_data)
655  call netcdf_err(error, 'WRITING XAXIS RECORD' )
656  error = nf90_put_var( ncid, id_y, y_data)
657  call netcdf_err(error, 'WRITING YAXIS RECORD' )
658  times = 1.0
659  error = nf90_put_var( ncid, id_time, times)
660  call netcdf_err(error, 'WRITING TIME RECORD' )
661 
662  deallocate(lsoil_data, x_data, y_data)
663 
664  dims_strt(1:3) = 1
665  dims_end(1) = idim
666  dims_end(2) = jdim
667  dims_end(3) = 1
668 
669  allocate(dum2d(idim,jdim))
670 
671  dum2d = reshape(slifcs, (/idim,jdim/))
672  error = nf90_put_var( ncid, id_slmsk, dum2d, dims_strt, dims_end)
673  call netcdf_err(error, 'WRITING LANDMASK RECORD' )
674 
675  dum2d = reshape(tsffcs, (/idim,jdim/))
676  error = nf90_put_var( ncid, id_tsea, dum2d, dims_strt, dims_end)
677  call netcdf_err(error, 'WRITING TSEA RECORD' )
678 
679  dum2d = reshape(swefcs, (/idim,jdim/))
680  error = nf90_put_var( ncid, id_sheleg, dum2d, dims_strt, dims_end)
681  call netcdf_err(error, 'WRITING SHELEG RECORD' )
682 
683  dum2d = reshape(tg3fcs, (/idim,jdim/))
684  error = nf90_put_var( ncid, id_tg3, dum2d, dims_strt, dims_end)
685  call netcdf_err(error, 'WRITING TG3 RECORD' )
686 
687  dum2d = reshape(zorfcs, (/idim,jdim/))
688  error = nf90_put_var( ncid, id_zorl, dum2d, dims_strt, dims_end)
689  call netcdf_err(error, 'WRITING ZORL RECORD' )
690 
691  dum2d = reshape(albfcs(:,1), (/idim,jdim/))
692  error = nf90_put_var( ncid, id_alvsf, dum2d, dims_strt, dims_end)
693  call netcdf_err(error, 'WRITING ALVSF RECORD' )
694 
695  dum2d = reshape(albfcs(:,2), (/idim,jdim/))
696  error = nf90_put_var( ncid, id_alvwf, dum2d, dims_strt, dims_end)
697  call netcdf_err(error, 'WRITING ALVWF RECORD' )
698 
699  dum2d = reshape(albfcs(:,3), (/idim,jdim/))
700  error = nf90_put_var( ncid, id_alnsf, dum2d, dims_strt, dims_end)
701  call netcdf_err(error, 'WRITING ALNSF RECORD' )
702 
703  dum2d = reshape(albfcs(:,4), (/idim,jdim/))
704  error = nf90_put_var( ncid, id_alnwf, dum2d, dims_strt, dims_end)
705  call netcdf_err(error, 'WRITING ALNWF RECORD' )
706 
707  dum2d = reshape(alffcs(:,1), (/idim,jdim/))
708  error = nf90_put_var( ncid, id_facsf, dum2d, dims_strt, dims_end)
709  call netcdf_err(error, 'WRITING FACSF RECORD' )
710 
711  dum2d = reshape(alffcs(:,2), (/idim,jdim/))
712  error = nf90_put_var( ncid, id_facwf, dum2d, dims_strt, dims_end)
713  call netcdf_err(error, 'WRITING FACWF RECORD' )
714 
715  dum2d = reshape(vegfcs, (/idim,jdim/))
716  error = nf90_put_var( ncid, id_vfrac, dum2d, dims_strt, dims_end)
717  call netcdf_err(error, 'WRITING VFRAC RECORD' )
718 
719  dum2d = reshape(cnpfcs, (/idim,jdim/))
720  error = nf90_put_var( ncid, id_canopy, dum2d, dims_strt, dims_end)
721  call netcdf_err(error, 'WRITING CANOPY RECORD' )
722 
723  dum2d = reshape(f10m, (/idim,jdim/))
724  error = nf90_put_var( ncid, id_f10m, dum2d, dims_strt, dims_end)
725  call netcdf_err(error, 'WRITING F10M RECORD' )
726 
727  dum2d = reshape(t2m, (/idim,jdim/))
728  error = nf90_put_var( ncid, id_t2m, dum2d, dims_strt, dims_end)
729  call netcdf_err(error, 'WRITING T2M RECORD' )
730 
731  dum2d = reshape(q2m, (/idim,jdim/))
732  error = nf90_put_var( ncid, id_q2m, dum2d, dims_strt, dims_end)
733  call netcdf_err(error, 'WRITING Q2M RECORD' )
734 
735  dum2d = reshape(vetfcs, (/idim,jdim/))
736  error = nf90_put_var( ncid, id_vtype, dum2d, dims_strt, dims_end)
737  call netcdf_err(error, 'WRITING VTYPE RECORD' )
738 
739  dum2d = reshape(sotfcs, (/idim,jdim/))
740  error = nf90_put_var( ncid, id_stype, dum2d, dims_strt, dims_end)
741  call netcdf_err(error, 'WRITING STYPE RECORD' )
742 
743  dum2d = reshape(ustar, (/idim,jdim/))
744  error = nf90_put_var( ncid, id_uustar, dum2d, dims_strt, dims_end)
745  call netcdf_err(error, 'WRITING UUSTAR RECORD' )
746 
747  dum2d = reshape(fmm, (/idim,jdim/))
748  error = nf90_put_var( ncid, id_ffmm, dum2d, dims_strt, dims_end)
749  call netcdf_err(error, 'WRITING FFMM RECORD' )
750 
751  dum2d = reshape(fhh, (/idim,jdim/))
752  error = nf90_put_var( ncid, id_ffhh, dum2d, dims_strt, dims_end)
753  call netcdf_err(error, 'WRITING FFHH RECORD' )
754 
755  dum2d = reshape(sihfcs, (/idim,jdim/))
756  error = nf90_put_var( ncid, id_hice, dum2d, dims_strt, dims_end)
757  call netcdf_err(error, 'WRITING HICE RECORD' )
758 
759  dum2d = reshape(sicfcs, (/idim,jdim/))
760  error = nf90_put_var( ncid, id_fice, dum2d, dims_strt, dims_end)
761  call netcdf_err(error, 'WRITING FICE RECORD' )
762 
763  dum2d = reshape(sitfcs, (/idim,jdim/))
764  error = nf90_put_var( ncid, id_tisfc, dum2d, dims_strt, dims_end)
765  call netcdf_err(error, 'WRITING TISFC RECORD' )
766 
767  dum2d = reshape(tprcp, (/idim,jdim/))
768  error = nf90_put_var( ncid, id_tprcp, dum2d, dims_strt, dims_end)
769  call netcdf_err(error, 'WRITING TPRCP RECORD' )
770 
771  dum2d = reshape(srflag, (/idim,jdim/))
772  error = nf90_put_var( ncid, id_srflag, dum2d, dims_strt, dims_end)
773  call netcdf_err(error, 'WRITING SRFLAG RECORD' )
774 
775  dum2d = reshape(swdfcs, (/idim,jdim/))
776  error = nf90_put_var( ncid, id_snwdph, dum2d, dims_strt, dims_end)
777  call netcdf_err(error, 'WRITING SNWDPH RECORD' )
778 
779  dum2d = reshape(vmnfcs, (/idim,jdim/))
780  error = nf90_put_var( ncid, id_shdmin, dum2d, dims_strt, dims_end)
781  call netcdf_err(error, 'WRITING SHDMIN RECORD' )
782 
783  dum2d = reshape(vmxfcs, (/idim,jdim/))
784  error = nf90_put_var( ncid, id_shdmax, dum2d, dims_strt, dims_end)
785  call netcdf_err(error, 'WRITING SHDMAX RECORD' )
786 
787  dum2d = reshape(slpfcs, (/idim,jdim/))
788  error = nf90_put_var( ncid, id_slope, dum2d, dims_strt, dims_end)
789  call netcdf_err(error, 'WRITING SLOPE RECORD' )
790 
791  dum2d = reshape(absfcs, (/idim,jdim/))
792  error = nf90_put_var( ncid, id_snoalb, dum2d, dims_strt, dims_end)
793  call netcdf_err(error, 'WRITING SNOALB RECORD' )
794 
795  nsst_write : if (do_nsst) then
796 
797  dum2d = reshape(nsst%tref, (/idim,jdim/))
798  error = nf90_put_var( ncid, id_tref, dum2d, dims_strt, dims_end)
799  call netcdf_err(error, 'WRITING TREF RECORD' )
800 
801  dum2d = reshape(nsst%z_c, (/idim,jdim/))
802  error = nf90_put_var( ncid, id_z_c, dum2d, dims_strt, dims_end)
803  call netcdf_err(error, 'WRITING Z_C RECORD' )
804 
805  dum2d = reshape(nsst%c_0, (/idim,jdim/))
806  error = nf90_put_var( ncid, id_c_0, dum2d, dims_strt, dims_end)
807  call netcdf_err(error, 'WRITING C_0 RECORD' )
808 
809  dum2d = reshape(nsst%c_d, (/idim,jdim/))
810  error = nf90_put_var( ncid, id_c_d, dum2d, dims_strt, dims_end)
811  call netcdf_err(error, 'WRITING C_D RECORD' )
812 
813  dum2d = reshape(nsst%w_0, (/idim,jdim/))
814  error = nf90_put_var( ncid, id_w_0, dum2d, dims_strt, dims_end)
815  call netcdf_err(error, 'WRITING W_0 RECORD' )
816 
817  dum2d = reshape(nsst%w_d, (/idim,jdim/))
818  error = nf90_put_var( ncid, id_w_d, dum2d, dims_strt, dims_end)
819  call netcdf_err(error, 'WRITING W_D RECORD' )
820 
821  dum2d = reshape(nsst%xt, (/idim,jdim/))
822  error = nf90_put_var( ncid, id_xt, dum2d, dims_strt, dims_end)
823  call netcdf_err(error, 'WRITING XT RECORD' )
824 
825  dum2d = reshape(nsst%xs, (/idim,jdim/))
826  error = nf90_put_var( ncid, id_xs, dum2d, dims_strt, dims_end)
827  call netcdf_err(error, 'WRITING XS RECORD' )
828 
829  dum2d = reshape(nsst%xu, (/idim,jdim/))
830  error = nf90_put_var( ncid, id_xu, dum2d, dims_strt, dims_end)
831  call netcdf_err(error, 'WRITING XU RECORD' )
832 
833  dum2d = reshape(nsst%xv, (/idim,jdim/))
834  error = nf90_put_var( ncid, id_xv, dum2d, dims_strt, dims_end)
835  call netcdf_err(error, 'WRITING XV RECORD' )
836 
837  dum2d = reshape(nsst%xz, (/idim,jdim/))
838  error = nf90_put_var( ncid, id_xz, dum2d, dims_strt, dims_end)
839  call netcdf_err(error, 'WRITING XZ RECORD' )
840 
841  dum2d = reshape(nsst%zm, (/idim,jdim/))
842  error = nf90_put_var( ncid, id_zm, dum2d, dims_strt, dims_end)
843  call netcdf_err(error, 'WRITING ZM RECORD' )
844 
845  dum2d = reshape(nsst%zm, (/idim,jdim/))
846  error = nf90_put_var( ncid, id_zm, dum2d, dims_strt, dims_end)
847  call netcdf_err(error, 'WRITING ZM RECORD' )
848 
849  dum2d = reshape(nsst%xtts, (/idim,jdim/))
850  error = nf90_put_var( ncid, id_xtts, dum2d, dims_strt, dims_end)
851  call netcdf_err(error, 'WRITING XTTS RECORD' )
852 
853  dum2d = reshape(nsst%xzts, (/idim,jdim/))
854  error = nf90_put_var( ncid, id_xzts, dum2d, dims_strt, dims_end)
855  call netcdf_err(error, 'WRITING XZTS RECORD' )
856 
857  dum2d = reshape(nsst%d_conv, (/idim,jdim/))
858  error = nf90_put_var( ncid, id_d_conv, dum2d, dims_strt, dims_end)
859  call netcdf_err(error, 'WRITING D_CONV RECORD' )
860 
861  dum2d = reshape(nsst%ifd, (/idim,jdim/))
862  error = nf90_put_var( ncid, id_ifd, dum2d, dims_strt, dims_end)
863  call netcdf_err(error, 'WRITING IFD RECORD' )
864 
865  dum2d = reshape(nsst%dt_cool, (/idim,jdim/))
866  error = nf90_put_var( ncid, id_dt_cool, dum2d, dims_strt, dims_end)
867  call netcdf_err(error, 'WRITING DT_COOL RECORD' )
868 
869  dum2d = reshape(nsst%qrain, (/idim,jdim/))
870  error = nf90_put_var( ncid, id_qrain, dum2d, dims_strt, dims_end)
871  call netcdf_err(error, 'WRITING QRAIN RECORD' )
872 
873  dum2d = reshape(nsst%tfinc, (/idim,jdim/))
874  error = nf90_put_var( ncid, id_tfinc, dum2d, dims_strt, dims_end)
875  call netcdf_err(error, 'WRITING TFINC RECORD' )
876 
877  endif nsst_write
878 
879  deallocate(dum2d)
880 
881  dims4_strt(1:4) = 1
882  dims4_end(1) = idim
883  dims4_end(2) = jdim
884  dims4_end(3) = lsoil
885  dims4_end(4) = 1
886 
887  allocate(dum3d(idim,jdim,lsoil))
888 
889  dum3d = reshape(slcfcs, (/idim,jdim,lsoil/))
890  error = nf90_put_var( ncid, id_slc, dum3d, dims4_strt, dims4_end)
891  call netcdf_err(error, 'WRITING SLC RECORD' )
892 
893  dum3d = reshape(smcfcs, (/idim,jdim,lsoil/))
894  error = nf90_put_var( ncid, id_smc, dum3d, dims4_strt, dims4_end)
895  call netcdf_err(error, 'WRITING SMC RECORD' )
896 
897  dum3d = reshape(stcfcs, (/idim,jdim,lsoil/))
898  error = nf90_put_var( ncid, id_stc, dum3d, dims4_strt, dims4_end)
899  call netcdf_err(error, 'WRITING STC RECORD' )
900 
901  deallocate(dum3d)
902 
903  error = nf90_close(ncid)
904 
905  end subroutine write_data
906 
920  SUBROUTINE read_lat_lon_orog(RLA,RLO,OROG,OROG_UF,&
921  tile_num,idim,jdim,ijdim)
922 
923  USE mpi
924 
925  IMPLICIT NONE
926 
927  INTEGER, INTENT(IN) :: idim, jdim, ijdim
928 
929  CHARACTER(LEN=5), INTENT(OUT) :: tile_num
930 
931  REAL, INTENT(OUT) :: rla(ijdim),rlo(ijdim)
932  REAL, INTENT(OUT) :: orog(ijdim),orog_uf(ijdim)
933 
934  CHARACTER(LEN=50) :: fnorog, fngrid
935  CHARACTER(LEN=3) :: rankch
936 
937  INTEGER :: error, ncid, ncid_orog
938  INTEGER :: i, ii, j, jj, myrank
939  INTEGER :: id_dim, id_var, nx, ny
940 
941  REAL, ALLOCATABLE :: dummy(:,:), geolat(:,:), geolon(:,:)
942  REAL(KIND=4), ALLOCATABLE :: dummy4(:,:)
943 
944  CALL mpi_comm_rank(mpi_comm_world, myrank, error)
945 
946  WRITE(rankch, '(I3.3)') (myrank+1)
947 
948  fngrid = "./fngrid." // rankch
949 
950  print*
951  print*, "READ FV3 GRID INFO FROM: "//trim(fngrid)
952 
953  error=nf90_open(trim(fngrid),nf90_nowrite,ncid)
954  CALL netcdf_err(error, 'OPENING FILE: '//trim(fngrid) )
955 
956  error=nf90_inq_dimid(ncid, 'nx', id_dim)
957  CALL netcdf_err(error, 'ERROR READING NX ID' )
958 
959  error=nf90_inquire_dimension(ncid,id_dim,len=nx)
960  CALL netcdf_err(error, 'ERROR READING NX' )
961 
962  error=nf90_inq_dimid(ncid, 'ny', id_dim)
963  CALL netcdf_err(error, 'ERROR READING NY ID' )
964 
965  error=nf90_inquire_dimension(ncid,id_dim,len=ny)
966  CALL netcdf_err(error, 'ERROR READING NY' )
967 
968  IF ((nx/2) /= idim .OR. (ny/2) /= jdim) THEN
969  print*,'FATAL ERROR: DIMENSIONS IN FILE: ',(nx/2),(ny/2)
970  print*,'DO NOT MATCH GRID DIMENSIONS: ',idim,jdim
971  CALL mpi_abort(mpi_comm_world, 130, error)
972  ENDIF
973 
974  ALLOCATE(geolon(nx+1,ny+1))
975  ALLOCATE(geolat(nx+1,ny+1))
976 
977  error=nf90_inq_varid(ncid, 'x', id_var)
978  CALL netcdf_err(error, 'ERROR READING X ID' )
979  error=nf90_get_var(ncid, id_var, geolon)
980  CALL netcdf_err(error, 'ERROR READING X RECORD' )
981 
982  error=nf90_inq_varid(ncid, 'y', id_var)
983  CALL netcdf_err(error, 'ERROR READING Y ID' )
984  error=nf90_get_var(ncid, id_var, geolat)
985  CALL netcdf_err(error, 'ERROR READING Y RECORD' )
986 
987  ALLOCATE(dummy(idim,jdim))
988 
989  DO j = 1, jdim
990  DO i = 1, idim
991  ii = 2*i
992  jj = 2*j
993  dummy(i,j) = geolon(ii,jj)
994  ENDDO
995  ENDDO
996 
997  rlo = reshape(dummy, (/ijdim/))
998 
999  DEALLOCATE(geolon)
1000 
1001  DO j = 1, jdim
1002  DO i = 1, idim
1003  ii = 2*i
1004  jj = 2*j
1005  dummy(i,j) = geolat(ii,jj)
1006  ENDDO
1007  ENDDO
1008 
1009  rla = reshape(dummy, (/ijdim/))
1010 
1011  DEALLOCATE(geolat, dummy)
1012 
1013  error=nf90_inq_varid(ncid, 'tile', id_var)
1014  CALL netcdf_err(error, 'ERROR READING TILE ID' )
1015  error=nf90_get_var(ncid, id_var, tile_num)
1016  CALL netcdf_err(error, 'ERROR READING TILE RECORD' )
1017 
1018  error = nf90_close(ncid)
1019 
1020  fnorog = "./fnorog." // rankch
1021 
1022  print*
1023  print*, "READ FV3 OROG INFO FROM: "//trim(fnorog)
1024 
1025  error=nf90_open(trim(fnorog),nf90_nowrite,ncid_orog)
1026  CALL netcdf_err(error, 'OPENING FILE: '//trim(fnorog) )
1027 
1028  ALLOCATE(dummy4(idim,jdim))
1029 
1030  error=nf90_inq_varid(ncid_orog, 'orog_raw', id_var)
1031  CALL netcdf_err(error, 'ERROR READING orog_raw ID' )
1032  error=nf90_get_var(ncid_orog, id_var, dummy4)
1033  CALL netcdf_err(error, 'ERROR READING orog_raw RECORD' )
1034  orog_uf = reshape(dummy4, (/ijdim/))
1035 
1036  error=nf90_inq_varid(ncid_orog, 'orog_filt', id_var)
1037  CALL netcdf_err(error, 'ERROR READING orog_filt ID' )
1038  error=nf90_get_var(ncid_orog, id_var, dummy4)
1039  CALL netcdf_err(error, 'ERROR READING orog_filt RECORD' )
1040  orog = reshape(dummy4, (/ijdim/))
1041 
1042  DEALLOCATE(dummy4)
1043 
1044  error = nf90_close(ncid_orog)
1045 
1046  END SUBROUTINE read_lat_lon_orog
1047 
1054  SUBROUTINE netcdf_err( ERR, STRING )
1055 
1056  USE mpi
1057 
1058  IMPLICIT NONE
1059 
1060  INTEGER, INTENT(IN) :: err
1061  CHARACTER(LEN=*), INTENT(IN) :: string
1062  CHARACTER(LEN=80) :: errmsg
1063  INTEGER :: iret
1064 
1065  IF( err == nf90_noerr )RETURN
1066  errmsg = nf90_strerror(err)
1067  print*,''
1068  print*,'FATAL ERROR: ', trim(string), ': ', trim(errmsg)
1069  print*,'STOP.'
1070  CALL mpi_abort(mpi_comm_world, 999, iret)
1071 
1072  RETURN
1073  END SUBROUTINE netcdf_err
1074 
1087  SUBROUTINE read_gsi_data(GSI_FILE, FILE_TYPE, LSOIL)
1088 
1089  IMPLICIT NONE
1090 
1091  CHARACTER(LEN=*), INTENT(IN) :: gsi_file
1092  CHARACTER(LEN=3), INTENT(IN) :: file_type
1093  INTEGER, INTENT(IN), OPTIONAL :: lsoil
1094 
1095  INTEGER :: error, id_dim, ncid
1096  INTEGER :: id_var, j
1097 
1098  INTEGER(KIND=1), ALLOCATABLE :: idummy(:,:)
1099 
1100  REAL(KIND=8), ALLOCATABLE :: dummy(:,:)
1101 
1102  CHARACTER(LEN=1) :: k_ch
1103  CHARACTER(LEN=10) :: incvar
1104  CHARACTER(LEN=80) :: err_msg
1105  INTEGER :: k, i
1106 
1107  print*
1108  print*, "READ INPUT GSI DATA FROM: "//trim(gsi_file)
1109 
1110  error=nf90_open(trim(gsi_file),nf90_nowrite,ncid)
1111  CALL netcdf_err(error, 'OPENING FILE: '//trim(gsi_file) )
1112 
1113  error=nf90_inq_dimid(ncid, 'latitude', id_dim)
1114  CALL netcdf_err(error, 'READING latitude ID' )
1115  error=nf90_inquire_dimension(ncid,id_dim,len=jdim_gaus)
1116  CALL netcdf_err(error, 'READING latitude length' )
1117  jdim_gaus = jdim_gaus - 2 ! WILL IGNORE POLE POINTS
1118 
1119  error=nf90_inq_dimid(ncid, 'longitude', id_dim)
1120  CALL netcdf_err(error, 'READING longitude ID' )
1121  error=nf90_inquire_dimension(ncid,id_dim,len=idim_gaus)
1122  CALL netcdf_err(error, 'READING longitude length' )
1123 
1124  IF (file_type=='NST') then
1125  ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
1126  ALLOCATE(dtref_gaus(idim_gaus,jdim_gaus))
1127 
1128  error=nf90_inq_varid(ncid, "dtf", id_var)
1129  CALL netcdf_err(error, 'READING dtf ID' )
1130  error=nf90_get_var(ncid, id_var, dummy)
1131  CALL netcdf_err(error, 'READING dtf' )
1132 
1133  ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
1134  ALLOCATE(slmsk_gaus(idim_gaus,jdim_gaus))
1135 
1136  error=nf90_inq_varid(ncid, "msk", id_var)
1137  CALL netcdf_err(error, 'READING msk ID' )
1138  error=nf90_get_var(ncid, id_var, idummy)
1139  CALL netcdf_err(error, 'READING msk' )
1140 
1141  ! REMOVE POLE POINTS.
1142 
1143  DO j = 1, jdim_gaus
1144  slmsk_gaus(:,j) = idummy(:,j+1)
1145  dtref_gaus(:,j) = dummy(:,j+1)
1146  ENDDO
1147 
1148  ELSEIF (file_type=='LND') then
1149 
1150  ALLOCATE(dummy(idim_gaus,jdim_gaus+2))
1151  ALLOCATE(stc_inc_gaus(lsoil,idim_gaus,jdim_gaus))
1152 
1153  ! read in soil temperature increments in each layer
1154  DO k = 1, lsoil
1155  WRITE(k_ch, '(I1)') k
1156 
1157  incvar = "soilt"//k_ch//"_inc"
1158  error=nf90_inq_varid(ncid, incvar, id_var)
1159  err_msg = "reading "//incvar//" ID"
1160  CALL netcdf_err(error, trim(err_msg))
1161  error=nf90_get_var(ncid, id_var, dummy)
1162  err_msg = "reading "//incvar//" data"
1163  CALL netcdf_err(error, err_msg)
1164 
1165  DO j = 1, jdim_gaus
1166  stc_inc_gaus(k,:,j) = dummy(:,j+1)
1167  ENDDO
1168  ENDDO
1169 
1170  ALLOCATE(idummy(idim_gaus,jdim_gaus+2))
1171  ALLOCATE(soilsnow_gaus(idim_gaus,jdim_gaus))
1172 
1173  error=nf90_inq_varid(ncid, "soilsnow_mask", id_var)
1174  CALL netcdf_err(error, 'READING soilsnow_mask ID' )
1175  error=nf90_get_var(ncid, id_var, idummy)
1176  CALL netcdf_err(error, 'READING soilsnow_mask' )
1177 
1178  ! REMOVE POLE POINTS.
1179 
1180  DO j = 1, jdim_gaus
1181  soilsnow_gaus(:,j) = idummy(:,j+1)
1182  ENDDO
1183 
1184 
1185  ELSE
1186  print *, 'WARNING: FILE_TYPE', file_type, 'not recognised.', &
1187  ', no increments read in'
1188  ENDIF
1189 
1190  IF(ALLOCATED(dummy)) DEALLOCATE(dummy)
1191  IF(ALLOCATED(idummy)) DEALLOCATE(idummy)
1192 
1193  error = nf90_close(ncid)
1194 
1195  END SUBROUTINE read_gsi_data
1196 
1244  SUBROUTINE read_data(LSOIL,LENSFC,DO_NSST,INC_FILE,TSFFCS,SMCFCS,SWEFCS,STCFCS, &
1245  tg3fcs,zorfcs, &
1246  cvfcs,cvbfcs,cvtfcs,albfcs, &
1247  vegfcs,slifcs,cnpfcs,f10m, &
1248  vetfcs,sotfcs,alffcs, &
1249  ustar,fmm,fhh, &
1250  sihfcs,sicfcs,sitfcs, &
1251  tprcp,srflag,sndfcs, &
1252  vmnfcs,vmxfcs,slcfcs, &
1253  slpfcs,absfcs,t2m,q2m,slmask, &
1254  zsoil,nsst)
1255  USE mpi
1256 
1257  IMPLICIT NONE
1258 
1259  INTEGER, INTENT(IN) :: lsoil, lensfc
1260  LOGICAL, INTENT(IN) :: do_nsst, inc_file
1261 
1262  REAL, OPTIONAL, INTENT(OUT) :: cvfcs(lensfc), cvbfcs(lensfc)
1263  REAL, OPTIONAL, INTENT(OUT) :: cvtfcs(lensfc), albfcs(lensfc,4)
1264  REAL, OPTIONAL, INTENT(OUT) :: slifcs(lensfc), cnpfcs(lensfc)
1265  REAL, OPTIONAL, INTENT(OUT) :: vegfcs(lensfc), f10m(lensfc)
1266  REAL, OPTIONAL, INTENT(OUT) :: vetfcs(lensfc), sotfcs(lensfc)
1267  REAL, OPTIONAL, INTENT(OUT) :: tsffcs(lensfc), swefcs(lensfc)
1268  REAL, OPTIONAL, INTENT(OUT) :: tg3fcs(lensfc), zorfcs(lensfc)
1269  REAL, OPTIONAL, INTENT(OUT) :: alffcs(lensfc,2), ustar(lensfc)
1270  REAL, OPTIONAL, INTENT(OUT) :: fmm(lensfc), fhh(lensfc)
1271  REAL, OPTIONAL, INTENT(OUT) :: sihfcs(lensfc), sicfcs(lensfc)
1272  REAL, OPTIONAL, INTENT(OUT) :: sitfcs(lensfc), tprcp(lensfc)
1273  REAL, OPTIONAL, INTENT(OUT) :: srflag(lensfc), sndfcs(lensfc)
1274  REAL, OPTIONAL, INTENT(OUT) :: vmnfcs(lensfc), vmxfcs(lensfc)
1275  REAL, OPTIONAL, INTENT(OUT) :: slpfcs(lensfc), absfcs(lensfc)
1276  REAL, OPTIONAL, INTENT(OUT) :: t2m(lensfc), q2m(lensfc), slmask(lensfc)
1277  REAL, OPTIONAL, INTENT(OUT) :: slcfcs(lensfc,lsoil)
1278  REAL, OPTIONAL, INTENT(OUT) :: smcfcs(lensfc,lsoil)
1279  REAL, OPTIONAL, INTENT(OUT) :: stcfcs(lensfc,lsoil)
1280  REAL(KIND=4), OPTIONAL, INTENT(OUT) :: zsoil(lsoil)
1281 
1282  TYPE(nsst_data), OPTIONAL :: nsst ! intent(out) will crash
1283  ! because subtypes are allocated in main.
1284 
1285  CHARACTER(LEN=50) :: fnbgsi
1286  CHARACTER(LEN=3) :: rankch
1287 
1288  INTEGER :: error, ncid, myrank
1289  INTEGER :: idim, jdim, id_dim
1290  INTEGER :: id_var, ierr
1291 
1292  REAL(KIND=8), ALLOCATABLE :: dummy(:,:), dummy3d(:,:,:)
1293 
1294  CALL mpi_comm_rank(mpi_comm_world, myrank, error)
1295 
1296  WRITE(rankch, '(I3.3)') (myrank+1)
1297 
1298  IF (inc_file) THEN
1299  fnbgsi = "./xainc." // rankch
1300  ELSE
1301  fnbgsi = "./fnbgsi." // rankch
1302  ENDIF
1303 
1304  print*
1305  print*, "READ INPUT SFC DATA FROM: "//trim(fnbgsi)
1306 
1307  error=nf90_open(trim(fnbgsi),nf90_nowrite,ncid)
1308  CALL netcdf_err(error, 'OPENING FILE: '//trim(fnbgsi) )
1309 
1310  error=nf90_inq_dimid(ncid, 'xaxis_1', id_dim)
1311  CALL netcdf_err(error, 'READING xaxis_1' )
1312  error=nf90_inquire_dimension(ncid,id_dim,len=idim)
1313  CALL netcdf_err(error, 'READING xaxis_1' )
1314 
1315  error=nf90_inq_dimid(ncid, 'yaxis_1', id_dim)
1316  CALL netcdf_err(error, 'READING yaxis_1' )
1317  error=nf90_inquire_dimension(ncid,id_dim,len=jdim)
1318  CALL netcdf_err(error, 'READING yaxis_1' )
1319 
1320  IF ((idim*jdim) /= lensfc) THEN
1321  print*,'FATAL ERROR: DIMENSIONS WRONG.'
1322  CALL mpi_abort(mpi_comm_world, 88, ierr)
1323  ENDIF
1324 
1325  ALLOCATE(dummy(idim,jdim))
1326 
1327  IF (present(tsffcs)) THEN
1328  error=nf90_inq_varid(ncid, "tsea", id_var)
1329  CALL netcdf_err(error, 'READING tsea ID' )
1330  error=nf90_get_var(ncid, id_var, dummy)
1331  CALL netcdf_err(error, 'READING tsea' )
1332  tsffcs = reshape(dummy, (/lensfc/))
1333  ENDIF
1334 
1335  IF (present(swefcs)) THEN
1336  error=nf90_inq_varid(ncid, "sheleg", id_var)
1337  CALL netcdf_err(error, 'READING sheleg ID' )
1338  error=nf90_get_var(ncid, id_var, dummy)
1339  CALL netcdf_err(error, 'READING sheleg' )
1340  swefcs = reshape(dummy, (/lensfc/))
1341  ENDIF
1342 
1343  IF (present(tg3fcs)) THEN
1344  error=nf90_inq_varid(ncid, "tg3", id_var)
1345  CALL netcdf_err(error, 'READING tg3 ID' )
1346  error=nf90_get_var(ncid, id_var, dummy)
1347  CALL netcdf_err(error, 'READING tg3' )
1348  tg3fcs = reshape(dummy, (/lensfc/))
1349  ENDIF
1350 
1351  IF (present(zorfcs)) THEN
1352  error=nf90_inq_varid(ncid, "zorl", id_var)
1353  CALL netcdf_err(error, 'READING zorl ID' )
1354  error=nf90_get_var(ncid, id_var, dummy)
1355  CALL netcdf_err(error, 'READING zorl' )
1356  zorfcs = reshape(dummy, (/lensfc/))
1357  ENDIF
1358 
1359  IF (present(albfcs)) THEN
1360 
1361  error=nf90_inq_varid(ncid, "alvsf", id_var)
1362  CALL netcdf_err(error, 'READING alvsf ID' )
1363  error=nf90_get_var(ncid, id_var, dummy)
1364  CALL netcdf_err(error, 'READING alvsf' )
1365  albfcs(:,1) = reshape(dummy, (/lensfc/))
1366 
1367  error=nf90_inq_varid(ncid, "alvwf", id_var)
1368  CALL netcdf_err(error, 'READING alvwf ID' )
1369  error=nf90_get_var(ncid, id_var, dummy)
1370  CALL netcdf_err(error, 'READING alvwf' )
1371  albfcs(:,2) = reshape(dummy, (/lensfc/))
1372 
1373  error=nf90_inq_varid(ncid, "alnsf", id_var)
1374  CALL netcdf_err(error, 'READING alnsf ID' )
1375  error=nf90_get_var(ncid, id_var, dummy)
1376  CALL netcdf_err(error, 'READING alnsf' )
1377  albfcs(:,3) = reshape(dummy, (/lensfc/))
1378 
1379  error=nf90_inq_varid(ncid, "alnwf", id_var)
1380  CALL netcdf_err(error, 'READING alnwf ID' )
1381  error=nf90_get_var(ncid, id_var, dummy)
1382  CALL netcdf_err(error, 'READING alnwf' )
1383  albfcs(:,4) = reshape(dummy, (/lensfc/))
1384 
1385  ENDIF
1386 
1387  IF (present(slifcs)) THEN
1388  error=nf90_inq_varid(ncid, "slmsk", id_var)
1389  CALL netcdf_err(error, 'READING slmsk ID' )
1390  error=nf90_get_var(ncid, id_var, dummy)
1391  CALL netcdf_err(error, 'READING slmsk' )
1392  slifcs = reshape(dummy, (/lensfc/))
1393  slmask = slifcs
1394  WHERE (slmask > 1.5) slmask=0.0 ! remove sea ice
1395  ENDIF
1396 
1397  IF (present(cnpfcs)) THEN
1398  error=nf90_inq_varid(ncid, "canopy", id_var)
1399  CALL netcdf_err(error, 'READING canopy ID' )
1400  error=nf90_get_var(ncid, id_var, dummy)
1401  CALL netcdf_err(error, 'READING canopy' )
1402  cnpfcs = reshape(dummy, (/lensfc/))
1403  ENDIF
1404 
1405  IF (present(vegfcs)) THEN
1406  error=nf90_inq_varid(ncid, "vfrac", id_var)
1407  CALL netcdf_err(error, 'READING vfrac ID' )
1408  error=nf90_get_var(ncid, id_var, dummy)
1409  CALL netcdf_err(error, 'READING vfrac' )
1410  vegfcs = reshape(dummy, (/lensfc/))
1411  ENDIF
1412 
1413  IF (present(f10m)) THEN
1414  error=nf90_inq_varid(ncid, "f10m", id_var)
1415  CALL netcdf_err(error, 'READING f10m ID' )
1416  error=nf90_get_var(ncid, id_var, dummy)
1417  CALL netcdf_err(error, 'READING f10m' )
1418  f10m = reshape(dummy, (/lensfc/))
1419  ENDIF
1420 
1421  IF (present(vetfcs)) THEN
1422  error=nf90_inq_varid(ncid, "vtype", id_var)
1423  CALL netcdf_err(error, 'READING vtype ID' )
1424  error=nf90_get_var(ncid, id_var, dummy)
1425  CALL netcdf_err(error, 'READING vtype' )
1426  vetfcs = reshape(dummy, (/lensfc/))
1427  ENDIF
1428 
1429  IF (present(sotfcs)) THEN
1430  error=nf90_inq_varid(ncid, "stype", id_var)
1431  CALL netcdf_err(error, 'READING stype ID' )
1432  error=nf90_get_var(ncid, id_var, dummy)
1433  CALL netcdf_err(error, 'READING stype' )
1434  sotfcs = reshape(dummy, (/lensfc/))
1435  ENDIF
1436 
1437  IF (present(alffcs)) THEN
1438  error=nf90_inq_varid(ncid, "facsf", id_var)
1439  CALL netcdf_err(error, 'READING facsf ID' )
1440  error=nf90_get_var(ncid, id_var, dummy)
1441  CALL netcdf_err(error, 'READING facsf' )
1442  alffcs(:,1) = reshape(dummy, (/lensfc/))
1443 
1444  error=nf90_inq_varid(ncid, "facwf", id_var)
1445  CALL netcdf_err(error, 'READING facwf ID' )
1446  error=nf90_get_var(ncid, id_var, dummy)
1447  CALL netcdf_err(error, 'READING facwf' )
1448  alffcs(:,2) = reshape(dummy, (/lensfc/))
1449  ENDIF
1450 
1451  IF (present(ustar)) THEN
1452  error=nf90_inq_varid(ncid, "uustar", id_var)
1453  CALL netcdf_err(error, 'READING uustar ID' )
1454  error=nf90_get_var(ncid, id_var, dummy)
1455  CALL netcdf_err(error, 'READING uustar' )
1456  ustar = reshape(dummy, (/lensfc/))
1457  ENDIF
1458 
1459  IF (present(fmm)) THEN
1460  error=nf90_inq_varid(ncid, "ffmm", id_var)
1461  CALL netcdf_err(error, 'READING ffmm ID' )
1462  error=nf90_get_var(ncid, id_var, dummy)
1463  CALL netcdf_err(error, 'READING ffmm' )
1464  fmm = reshape(dummy, (/lensfc/))
1465  ENDIF
1466 
1467  IF (present(fhh)) THEN
1468  error=nf90_inq_varid(ncid, "ffhh", id_var)
1469  CALL netcdf_err(error, 'READING ffhh ID' )
1470  error=nf90_get_var(ncid, id_var, dummy)
1471  CALL netcdf_err(error, 'READING ffhh' )
1472  fhh = reshape(dummy, (/lensfc/))
1473  ENDIF
1474 
1475  IF (present(sihfcs)) THEN
1476  error=nf90_inq_varid(ncid, "hice", id_var)
1477  CALL netcdf_err(error, 'READING hice ID' )
1478  error=nf90_get_var(ncid, id_var, dummy)
1479  CALL netcdf_err(error, 'READING hice' )
1480  sihfcs = reshape(dummy, (/lensfc/))
1481  ENDIF
1482 
1483  IF (present(sicfcs)) THEN
1484  error=nf90_inq_varid(ncid, "fice", id_var)
1485  CALL netcdf_err(error, 'READING fice ID' )
1486  error=nf90_get_var(ncid, id_var, dummy)
1487  CALL netcdf_err(error, 'READING fice' )
1488  sicfcs = reshape(dummy, (/lensfc/))
1489  ENDIF
1490 
1491  IF (present(sitfcs)) THEN
1492  error=nf90_inq_varid(ncid, "tisfc", id_var)
1493  CALL netcdf_err(error, 'READING tisfc ID' )
1494  error=nf90_get_var(ncid, id_var, dummy)
1495  CALL netcdf_err(error, 'READING tisfc' )
1496  sitfcs = reshape(dummy, (/lensfc/))
1497  ENDIF
1498 
1499  IF (present(tprcp)) THEN
1500  error=nf90_inq_varid(ncid, "tprcp", id_var)
1501  CALL netcdf_err(error, 'READING tprcp ID' )
1502  error=nf90_get_var(ncid, id_var, dummy)
1503  CALL netcdf_err(error, 'READING tprcp' )
1504  tprcp = reshape(dummy, (/lensfc/))
1505  ENDIF
1506 
1507  IF (present(srflag)) THEN
1508  error=nf90_inq_varid(ncid, "srflag", id_var)
1509  CALL netcdf_err(error, 'READING srflag ID' )
1510  error=nf90_get_var(ncid, id_var, dummy)
1511  CALL netcdf_err(error, 'READING srflag' )
1512  srflag = reshape(dummy, (/lensfc/))
1513  ENDIF
1514 
1515  IF (present(sndfcs)) THEN
1516  error=nf90_inq_varid(ncid, "snwdph", id_var)
1517  CALL netcdf_err(error, 'READING snwdph ID' )
1518  error=nf90_get_var(ncid, id_var, dummy)
1519  CALL netcdf_err(error, 'READING snwdph' )
1520  sndfcs = reshape(dummy, (/lensfc/))
1521  ENDIF
1522 
1523  IF (present(vmnfcs)) THEN
1524  error=nf90_inq_varid(ncid, "shdmin", id_var)
1525  CALL netcdf_err(error, 'READING shdmin ID' )
1526  error=nf90_get_var(ncid, id_var, dummy)
1527  CALL netcdf_err(error, 'READING shdmin' )
1528  vmnfcs = reshape(dummy, (/lensfc/))
1529  ENDIF
1530 
1531  IF (present(vmxfcs)) THEN
1532  error=nf90_inq_varid(ncid, "shdmax", id_var)
1533  CALL netcdf_err(error, 'READING shdmax ID' )
1534  error=nf90_get_var(ncid, id_var, dummy)
1535  CALL netcdf_err(error, 'READING shdmax' )
1536  vmxfcs = reshape(dummy, (/lensfc/))
1537  ENDIF
1538 
1539  IF (present(slpfcs)) THEN
1540  error=nf90_inq_varid(ncid, "slope", id_var)
1541  CALL netcdf_err(error, 'READING slope ID' )
1542  error=nf90_get_var(ncid, id_var, dummy)
1543  CALL netcdf_err(error, 'READING slope' )
1544  slpfcs = reshape(dummy, (/lensfc/))
1545  ENDIF
1546 
1547  IF (present(absfcs)) THEN
1548  error=nf90_inq_varid(ncid, "snoalb", id_var)
1549  CALL netcdf_err(error, 'READING snoalb ID' )
1550  error=nf90_get_var(ncid, id_var, dummy)
1551  CALL netcdf_err(error, 'READING snoalb' )
1552  absfcs = reshape(dummy, (/lensfc/))
1553  ENDIF
1554 
1555  IF (present(t2m)) THEN
1556  error=nf90_inq_varid(ncid, "t2m", id_var)
1557  CALL netcdf_err(error, 'READING t2m ID' )
1558  error=nf90_get_var(ncid, id_var, dummy)
1559  CALL netcdf_err(error, 'READING t2m' )
1560  t2m = reshape(dummy, (/lensfc/))
1561  ENDIF
1562 
1563  IF (present(q2m)) THEN
1564  error=nf90_inq_varid(ncid, "q2m", id_var)
1565  CALL netcdf_err(error, 'READING q2m ID' )
1566  error=nf90_get_var(ncid, id_var, dummy)
1567  CALL netcdf_err(error, 'READING q2m' )
1568  q2m = reshape(dummy, (/lensfc/))
1569  ENDIF
1570 
1571  nsst_read : IF(do_nsst) THEN
1572 
1573  print*
1574  print*,"WILL READ NSST RECORDS."
1575 
1576  error=nf90_inq_varid(ncid, "c_0", id_var)
1577  CALL netcdf_err(error, 'READING c_0 ID' )
1578  error=nf90_get_var(ncid, id_var, dummy)
1579  CALL netcdf_err(error, 'READING c_0' )
1580  nsst%C_0 = reshape(dummy, (/lensfc/))
1581 
1582  error=nf90_inq_varid(ncid, "c_d", id_var)
1583  CALL netcdf_err(error, 'READING c_d ID' )
1584  error=nf90_get_var(ncid, id_var, dummy)
1585  CALL netcdf_err(error, 'READING c_d' )
1586  nsst%C_D = reshape(dummy, (/lensfc/))
1587 
1588  error=nf90_inq_varid(ncid, "d_conv", id_var)
1589  CALL netcdf_err(error, 'READING d_conv ID' )
1590  error=nf90_get_var(ncid, id_var, dummy)
1591  CALL netcdf_err(error, 'READING d_conv' )
1592  nsst%D_CONV = reshape(dummy, (/lensfc/))
1593 
1594  error=nf90_inq_varid(ncid, "dt_cool", id_var)
1595  CALL netcdf_err(error, 'READING dt_cool ID' )
1596  error=nf90_get_var(ncid, id_var, dummy)
1597  CALL netcdf_err(error, 'READING dt_cool' )
1598  nsst%DT_COOL = reshape(dummy, (/lensfc/))
1599 
1600  error=nf90_inq_varid(ncid, "ifd", id_var)
1601  CALL netcdf_err(error, 'READING ifd ID' )
1602  error=nf90_get_var(ncid, id_var, dummy)
1603  CALL netcdf_err(error, 'READING ifd' )
1604  nsst%IFD = reshape(dummy, (/lensfc/))
1605 
1606  error=nf90_inq_varid(ncid, "qrain", id_var)
1607  CALL netcdf_err(error, 'READING qrain ID' )
1608  error=nf90_get_var(ncid, id_var, dummy)
1609  CALL netcdf_err(error, 'READING qrain' )
1610  nsst%QRAIN = reshape(dummy, (/lensfc/))
1611 
1612  error=nf90_inq_varid(ncid, "tref", id_var)
1613  CALL netcdf_err(error, 'READING tref ID' )
1614  error=nf90_get_var(ncid, id_var, dummy)
1615  CALL netcdf_err(error, 'READING tref' )
1616  nsst%TREF = reshape(dummy, (/lensfc/))
1617 
1618  error=nf90_inq_varid(ncid, "w_0", id_var)
1619  CALL netcdf_err(error, 'READING w_0 ID' )
1620  error=nf90_get_var(ncid, id_var, dummy)
1621  CALL netcdf_err(error, 'READING w_0' )
1622  nsst%W_0 = reshape(dummy, (/lensfc/))
1623 
1624  error=nf90_inq_varid(ncid, "w_d", id_var)
1625  CALL netcdf_err(error, 'READING w_d ID' )
1626  error=nf90_get_var(ncid, id_var, dummy)
1627  CALL netcdf_err(error, 'READING w_d' )
1628  nsst%W_D = reshape(dummy, (/lensfc/))
1629 
1630  error=nf90_inq_varid(ncid, "xs", id_var)
1631  CALL netcdf_err(error, 'READING xs ID' )
1632  error=nf90_get_var(ncid, id_var, dummy)
1633  CALL netcdf_err(error, 'READING xs' )
1634  nsst%XS = reshape(dummy, (/lensfc/))
1635 
1636  error=nf90_inq_varid(ncid, "xt", id_var)
1637  CALL netcdf_err(error, 'READING xt ID' )
1638  error=nf90_get_var(ncid, id_var, dummy)
1639  CALL netcdf_err(error, 'READING xt' )
1640  nsst%XT = reshape(dummy, (/lensfc/))
1641 
1642  error=nf90_inq_varid(ncid, "xtts", id_var)
1643  CALL netcdf_err(error, 'READING xtts ID' )
1644  error=nf90_get_var(ncid, id_var, dummy)
1645  CALL netcdf_err(error, 'READING xtts' )
1646  nsst%XTTS = reshape(dummy, (/lensfc/))
1647 
1648  error=nf90_inq_varid(ncid, "xu", id_var)
1649  CALL netcdf_err(error, 'READING xu ID' )
1650  error=nf90_get_var(ncid, id_var, dummy)
1651  CALL netcdf_err(error, 'READING xu' )
1652  nsst%XU = reshape(dummy, (/lensfc/))
1653 
1654  error=nf90_inq_varid(ncid, "xv", id_var)
1655  CALL netcdf_err(error, 'READING xv ID' )
1656  error=nf90_get_var(ncid, id_var, dummy)
1657  CALL netcdf_err(error, 'READING xv' )
1658  nsst%XV = reshape(dummy, (/lensfc/))
1659 
1660  error=nf90_inq_varid(ncid, "xz", id_var)
1661  CALL netcdf_err(error, 'READING xz ID' )
1662  error=nf90_get_var(ncid, id_var, dummy)
1663  CALL netcdf_err(error, 'READING xz' )
1664  nsst%XZ = reshape(dummy, (/lensfc/))
1665 
1666  error=nf90_inq_varid(ncid, "xzts", id_var)
1667  CALL netcdf_err(error, 'READING xzts ID' )
1668  error=nf90_get_var(ncid, id_var, dummy)
1669  CALL netcdf_err(error, 'READING xzts' )
1670  nsst%XZTS = reshape(dummy, (/lensfc/))
1671 
1672  error=nf90_inq_varid(ncid, "z_c", id_var)
1673  CALL netcdf_err(error, 'READING z_c ID' )
1674  error=nf90_get_var(ncid, id_var, dummy)
1675  CALL netcdf_err(error, 'READING z_c' )
1676  nsst%Z_C = reshape(dummy, (/lensfc/))
1677 
1678  error=nf90_inq_varid(ncid, "zm", id_var)
1679  CALL netcdf_err(error, 'READING zm ID' )
1680  error=nf90_get_var(ncid, id_var, dummy)
1681  CALL netcdf_err(error, 'READING zm' )
1682  nsst%ZM = reshape(dummy, (/lensfc/))
1683 
1684  END IF nsst_read
1685 
1686  DEALLOCATE(dummy)
1687 
1688  ALLOCATE(dummy3d(idim,jdim,lsoil))
1689 
1690  IF (present(smcfcs)) THEN
1691  error=nf90_inq_varid(ncid, "smc", id_var)
1692  CALL netcdf_err(error, 'READING smc ID' )
1693  error=nf90_get_var(ncid, id_var, dummy3d)
1694  CALL netcdf_err(error, 'READING smc' )
1695  smcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1696  ENDIF
1697 
1698  IF (present(slcfcs)) THEN
1699  error=nf90_inq_varid(ncid, "slc", id_var)
1700  CALL netcdf_err(error, 'READING slc ID' )
1701  error=nf90_get_var(ncid, id_var, dummy3d)
1702  CALL netcdf_err(error, 'READING slc' )
1703  slcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1704  ENDIF
1705 
1706  IF (present(stcfcs)) THEN
1707  error=nf90_inq_varid(ncid, "stc", id_var)
1708  CALL netcdf_err(error, 'READING stc ID' )
1709  error=nf90_get_var(ncid, id_var, dummy3d)
1710  CALL netcdf_err(error, 'READING stc' )
1711  stcfcs = reshape(dummy3d, (/lensfc,lsoil/))
1712  ENDIF
1713 
1714  DEALLOCATE(dummy3d)
1715 
1716 ! cloud fields not in warm restart files. set to zero?
1717 
1718  IF (present(cvfcs)) cvfcs = 0.0
1719  IF (present(cvtfcs)) cvtfcs = 0.0
1720  IF (present(cvbfcs)) cvbfcs = 0.0
1721 
1722 ! soil layer thicknesses not in warm restart files. hardwire
1723 ! for now.
1724 
1725  IF (present(zsoil)) THEN
1726  zsoil(1) = -0.1
1727  zsoil(2) = -0.4
1728  zsoil(3) = -1.0
1729  zsoil(4) = -2.0
1730  ENDIF
1731 
1732  error = nf90_close(ncid)
1733 
1734  END SUBROUTINE read_data
1735 
1752 subroutine read_tf_clim_grb(file_sst,sst,rlats_sst,rlons_sst,mlat_sst,mlon_sst,mon)
1753 
1754  use mpi
1755 
1756  implicit none
1757 
1758 ! declare passed variables and arrays
1759  character(*) , intent(in ) :: file_sst
1760  integer , intent(in ) :: mon,mlat_sst,mlon_sst
1761  real, dimension(mlat_sst) , intent( out) :: rlats_sst
1762  real, dimension(mlon_sst) , intent( out) :: rlons_sst
1763  real, dimension(mlon_sst,mlat_sst) , intent( out) :: sst
1764 
1765 ! declare local parameters
1766  integer,parameter:: lu_sst = 21 ! fortran unit number of grib sst file
1767  real, parameter :: deg2rad = 3.141593/180.0
1768 
1769 ! declare local variables and arrays
1770  logical(1), allocatable, dimension(:) :: lb
1771 
1772  integer :: nlat_sst
1773  integer :: nlon_sst
1774  integer :: iret,ni,nj
1775  integer :: mscan,kb1,ierr
1776  integer :: jincdir,i,iincdir,kb2,kb3,kf,kg,k,j,jf
1777  integer, dimension(22):: jgds,kgds
1778  integer, dimension(25):: jpds,kpds
1779 
1780  real :: xsst0
1781  real :: ysst0
1782  real :: dres
1783  real, allocatable, dimension(:) :: f
1784 
1785 !************+******************************************************************************
1786 !
1787 ! open sst analysis file (grib)
1788  write(*,*) ' sstclm : ',file_sst
1789  call baopenr(lu_sst,trim(file_sst),iret)
1790  if (iret /= 0 ) then
1791  write(6,*)'FATAL ERROR in read_tf_clm_grb: error opening sst file.'
1792  CALL mpi_abort(mpi_comm_world, 111, ierr)
1793  endif
1794 
1795 ! define sst variables for read
1796  j=-1
1797  jpds=-1
1798  jgds=-1
1799  jpds(5)=11 ! sst variable
1800  jpds(6)=1 ! surface
1801  jpds(9) = mon
1802  call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1803 
1804  nlat_sst = kgds(3) ! number points on longitude circle (360)
1805  nlon_sst = kgds(2) ! number points on latitude circle (720)
1806 
1807 ! write(*,*) 'nlat_sst, nlon_sst, mon : ',nlat_sst, nlon_sst, mon
1808 ! write(*,*) 'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1809 
1810 ! allocate arrays
1811  allocate(lb(nlat_sst*nlon_sst))
1812  allocate(f(nlat_sst*nlon_sst))
1813  jf=nlat_sst*nlon_sst
1814 
1815 ! read in the analysis
1816  call getgb(lu_sst,0,jf,j,jpds,jgds,kf,k,kpds,kgds,lb,f,iret)
1817  if (iret /= 0) then
1818  write(6,*)'FATAL ERROR in read_tf_clm_grb: error reading sst analysis data record.'
1819  deallocate(lb,f)
1820  CALL mpi_abort(mpi_comm_world, 111, ierr)
1821  endif
1822 
1823  if ( (nlat_sst /= mlat_sst) .or. (nlon_sst /= mlon_sst) ) then
1824  write(6,*)'FATAL ERROR in read_rtg_org: inconsistent dimensions. mlat_sst,mlon_sst=',&
1825  mlat_sst,mlon_sst,' -versus- nlat_sst,nlon_sst=',nlat_sst,nlon_sst
1826  deallocate(lb,f)
1827  CALL mpi_abort(mpi_comm_world, 111, ierr)
1828  endif
1829 
1830 !
1831 ! get xlats and xlons
1832 !
1833  dres = 180.0/real(nlat_sst)
1834  ysst0 = 0.5*dres-90.0
1835  xsst0 = 0.5*dres
1836 
1837 ! get lat_sst & lon_sst
1838  do j = 1, nlat_sst
1839  rlats_sst(j) = ysst0 + real(j-1)*dres
1840  enddo
1841 
1842  do i = 1, nlon_sst
1843  rlons_sst(i) = (xsst0 + real(i-1)*dres)
1844  enddo
1845 
1846 ! load dimensions and grid specs. check for unusual values
1847  ni=kgds(2) ! 720 for 0.5 x 0.5
1848  nj=kgds(3) ! 360 for 0.5 x 0.5 resolution
1849 
1850  mscan=kgds(11)
1851  kb1=ibits(mscan,7,1) ! i scan direction
1852  kb2=ibits(mscan,6,1) ! j scan direction
1853  kb3=ibits(mscan,5,1) ! (i,j) or (j,i)
1854 
1855 ! get i and j scanning directions from kb1 and kb2.
1856 ! 0 yields +1, 1 yields -1. +1 is west to east, -1 is east to west.
1857  iincdir = 1-kb1*2
1858 
1859 ! 0 yields -1, 1 yields +1. +1 is south to north, -1 is north to south.
1860  jincdir = kb2*2 - 1
1861 
1862 ! write(*,*) 'read_tf_clim_grb,iincdir,jincdir : ',iincdir,jincdir
1863  do k=1,kf
1864 
1865 ! kb3 from scan mode indicates if i points are consecutive
1866 ! or if j points are consecutive
1867  if(kb3==0)then ! (i,j)
1868  i=(ni+1)*kb1+(mod(k-1,ni)+1)*iincdir
1869  j=(nj+1)*(1-kb2)+jincdir*((k-1)/ni+1)
1870  else ! (j,i)
1871  j=(nj+1)*(1-kb2)+(mod(k-1,nj)+1)*jincdir
1872  i=(ni+1)*kb1+iincdir*((k-1)/nj+1)
1873  endif
1874  sst(i,j) = f(k)
1875  end do
1876 
1877  deallocate(lb,f)
1878 
1879  call baclose(lu_sst,iret)
1880  if (iret /= 0 ) then
1881  write(6,*)'FATAL ERROR in read_tf_clm_grb: error closing sst file.'
1882  CALL mpi_abort(mpi_comm_world, 121, ierr)
1883  endif
1884 
1885 end subroutine read_tf_clim_grb
1886 
1894 subroutine get_tf_clm_dim(file_sst,mlat_sst,mlon_sst)
1895  use mpi
1896 
1897  implicit none
1898 
1899 ! declare passed variables and arrays
1900  character(*) , intent(in ) :: file_sst
1901  integer , intent(out) :: mlat_sst,mlon_sst
1902 
1903 ! declare local parameters
1904  integer,parameter:: lu_sst = 21 ! fortran unit number of grib sst file
1905 
1906  integer :: iret
1907  integer :: mscan,kb1
1908  integer :: kf,kg,k,j,ierr
1909  integer, dimension(22):: jgds,kgds
1910  integer, dimension(25):: jpds,kpds
1911 
1912 !************+******************************************************************************
1913 !
1914 ! open sst analysis file (grib)
1915  call baopenr(lu_sst,trim(file_sst),iret)
1916  if (iret /= 0 ) then
1917  write(6,*)'FATAL ERROR in get_tf_clm_dim: error opening sst file.'
1918  CALL mpi_abort(mpi_comm_world, 111, ierr)
1919  endif
1920 
1921 ! define sst variables for read
1922  j=-1
1923  jpds=-1
1924  jgds=-1
1925  jpds(5)=11 ! sst variable
1926  jpds(6)=1 ! surface
1927  jpds(9) = 1
1928  call getgbh(lu_sst,0,j,jpds,jgds,kg,kf,k,kpds,kgds,iret)
1929 
1930  mlat_sst = kgds(3) ! number points on longitude circle (360)
1931  mlon_sst = kgds(2) ! number points on latitude circle (720)
1932 
1933  write(*,*) 'mlat_sst, mlon_sst : ',mlat_sst, mlon_sst
1934 
1935  call baclose(lu_sst,iret)
1936  if (iret /= 0 ) then
1937  write(6,*)'FATAL ERROR in get_tf_clm_dim: error closing sst file.'
1938  CALL mpi_abort(mpi_comm_world, 121, ierr)
1939  endif
1940 end subroutine get_tf_clm_dim
1941 
1953 subroutine read_salclm_gfs_nc(filename,sal,xlats,xlons,nlat,nlon,itime)
1954  use netcdf
1955  implicit none
1956 
1957  ! This is the name of the data file we will read.
1958  character (len=*), intent(in) :: filename
1959  integer, intent(in) :: nlat,nlon
1960  integer, intent(in) :: itime
1961  real, dimension(nlat), intent(out) :: xlats
1962  real, dimension(nlon), intent(out) :: xlons
1963  real, dimension(nlon,nlat), intent(out) :: sal
1964 ! Local variables
1965  integer :: ncid,ntime
1966 
1967  integer, parameter :: ndims = 3
1968  character (len = *), parameter :: lat_name = "latitude"
1969  character (len = *), parameter :: lon_name = "longitude"
1970  character (len = *), parameter :: t_name = "time"
1971  character (len = *), parameter :: sal_name="sal"
1972  integer :: no_fill,fill_value
1973  integer :: time_varid,lon_varid, lat_varid, z_varid, sal_varid
1974 
1975  ! The start and count arrays will tell the netCDF library where to read our data.
1976  integer, dimension(ndims) :: start, count
1977 
1978  character (len = *), parameter :: units = "units"
1979  character (len = *), parameter :: sal_units = "psu"
1980  ! PSU (Practical SalinitUnit). 1 PSU = 1g/kg
1981  character (len = *), parameter :: time_units = "months"
1982  character (len = *), parameter :: lat_units = "degrees_north"
1983  character (len = *), parameter :: lon_units = "degrees_east"
1984 
1985  integer :: missv
1986 ! Loop indices
1987  integer :: i,j
1988 
1989 ! Open the file.
1990  call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
1991 
1992 ! Get the varids of time, latitude, longitude & depth coordinate variables.
1993  call nc_check( nf90_inq_varid(ncid, t_name, time_varid) )
1994  call nc_check( nf90_inq_varid(ncid, lat_name, lat_varid) )
1995  call nc_check( nf90_inq_varid(ncid, lon_name, lon_varid) )
1996 
1997 ! Read the time, latitude and longitude data.
1998 ! call nc_check( nf90_get_var(ncid, time_varid, ntime) )
1999  call nc_check( nf90_get_var(ncid, lat_varid, xlats) )
2000  call nc_check( nf90_get_var(ncid, lon_varid, xlons) )
2001 
2002 ! Get the varids of the sal netCDF variables.
2003  call nc_check( nf90_inq_varid(ncid, sal_name,sal_varid) )
2004 
2005 ! Read 1 record of nlat*nlon values, starting at the beginning
2006 ! of the record (the (1, 1, 1, rec) element in the netCDF file).
2007  start = (/ 1, 1, itime /)
2008  count = (/ nlon, nlat, 1 /)
2009 
2010 ! write(*,*) 'read_salclm_gfs_nc itime : ',itime
2011 ! Read the sal data from the file, one record at a time.
2012  call nc_check( nf90_get_var(ncid, sal_varid, sal, start, count) )
2013 
2014 ! Close the file. This frees up any internal netCDF resources
2015 ! associated with the file.
2016  call nc_check( nf90_close(ncid) )
2017 
2018 ! If we got this far, everything worked as expected. Yipee!
2019 ! print *,"*** SUCCESS reading file ", filename, "!"
2020 
2021 end subroutine read_salclm_gfs_nc
2022 
2029 subroutine get_dim_nc(filename,nlat,nlon)
2030  use netcdf
2031  implicit none
2032 
2033  character (len=*), intent(in) :: filename
2034  integer, intent(out) :: nlat,nlon
2035 ! Local variables
2036  character (len = *), parameter :: lat_name = "latitude"
2037  character (len = *), parameter :: lon_name = "longitude"
2038  integer :: ncid
2039  integer :: latdimid,londimid
2040 
2041 ! Open the file.
2042  call nc_check( nf90_open(filename, nf90_nowrite, ncid) )
2043 
2044 ! Get dimensions
2045  call nc_check( nf90_inq_dimid(ncid,lat_name,latdimid) )
2046  call nc_check( nf90_inq_dimid(ncid,lon_name,londimid) )
2047  call nc_check( nf90_inquire_dimension(ncid,latdimid,len=nlat) )
2048  call nc_check( nf90_inquire_dimension(ncid,londimid,len=nlon) )
2049 
2050 ! write(*,'(a,1x,a6,2I8)') 'get_dim_nc, file, nlat, nlon : ',filename,nlat,nlon
2051 
2052 ! Close the file. This frees up any internal netCDF resources
2053 ! associated with the file.
2054  call nc_check( nf90_close(ncid) )
2055 
2056 ! If we got this far, everything worked as expected. Yipee!
2057 ! print *,"*** SUCCESS get dimensions from nc file ", filename, "!"
2058 
2059 end subroutine get_dim_nc
2060 
2066 subroutine nc_check(status)
2067 
2068  use mpi
2069  use netcdf
2070 
2071  integer, intent ( in) :: status
2072  integer :: ierr
2073 
2074  if(status /= nf90_noerr) then
2075  print *, 'FATAL ERROR:'
2076  print *, trim(nf90_strerror(status))
2077  CALL mpi_abort(mpi_comm_world, 122, ierr)
2078  end if
2079 end subroutine nc_check
2080 
2081  END MODULE read_write_data
subroutine, public read_lat_lon_orog(RLA, RLO, OROG, OROG_UF, TILE_NUM, IDIM, JDIM, IJDIM)
Read latitude and longitude for the cubed-sphere tile from the 'grid' file.
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 get_dim_nc(filename, nlat, nlon)
Get the i/j dimensions of the data from a NetCDF file.
subroutine, public read_salclm_gfs_nc(filename, sal, xlats, xlons, nlat, nlon, itime)
Read the woa05 salinity monthly climatology file.
subroutine, public read_data(LSOIL, LENSFC, DO_NSST, INC_FILE, 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...
subroutine, public get_tf_clm_dim(file_sst, mlat_sst, mlon_sst)
Get the i/j dimensions of RTG SST climatology file.
subroutine, public write_data(slifcs, tsffcs, swefcs, tg3fcs, zorfcs, albfcs, alffcs, vegfcs, cnpfcs, f10m, t2m, q2m, vetfcs, sotfcs, ustar, fmm, fhh, sicfcs, sihfcs, sitfcs, tprcp, srflag, swdfcs, vmnfcs, vmxfcs, slpfcs, absfcs, slcfcs, smcfcs, stcfcs, idim, jdim, lensfc, lsoil, do_nsst, nsst)
Write out all surface records - and nsst records if selected - on a single cubed-sphere tile to a mod...
subroutine, public read_gsi_data(GSI_FILE, FILE_TYPE, LSOIL)
Read increment file from the GSI containing either the foundation temperature increments and mask...
This module contains routines that read and write data.
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.