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