sfc_climo_gen  1.8.0
output.f90
Go to the documentation of this file.
1 
4 
18  subroutine output(data_one_tile, lat_one_tile, lon_one_tile, i_mdl, j_mdl, &
19  tile, record, time, field_idx)
20 
21  use mpi
22  use esmf
23  use netcdf
24  use utils
25  use source_grid, only : field_names, source, num_fields, &
27  use model_grid, only : missing, grid_tiles
28  use program_setup, only : halo
29 
30  implicit none
31 
32  integer, intent(in) :: i_mdl, j_mdl, tile
33  integer, intent(in) :: record, time, field_idx
34 
35  real(esmf_kind_r4), intent(in) :: data_one_tile(i_mdl,j_mdl)
36  real(esmf_kind_r4), intent(in) :: lat_one_tile(i_mdl,j_mdl)
37  real(esmf_kind_r4), intent(in) :: lon_one_tile(i_mdl,j_mdl)
38 
39  character(len=200) :: out_file
40  character(len=200) :: out_file_with_halo
41 
42  integer :: initialsiz, fsize, error, j
43  integer :: dim_x, dim_y, id_data
44  integer :: dim_time, id_times, ierr
45  integer :: header_buffer_val = 16384
46  integer :: i_out, j_out, id_lat, id_lon
47  integer :: i_start, i_end, j_start, j_end
48  integer, save :: ncid(6), ncid_with_halo
49 
50  select case (field_names(field_idx))
51  case ('substrate_temperature')
52  out_file = "./substrate_temperature." // grid_tiles(tile) // ".nc"
53  out_file_with_halo = "./substrate_temperature." // grid_tiles(tile) // ".halo.nc"
54  case ('vegetation_greenness')
55  out_file = "./vegetation_greenness." // grid_tiles(tile) // ".nc"
56  out_file_with_halo = "./vegetation_greenness." // grid_tiles(tile) // ".halo.nc"
57  case ('maximum_snow_albedo')
58  out_file = "./maximum_snow_albedo." // grid_tiles(tile) // ".nc"
59  out_file_with_halo = "./maximum_snow_albedo." // grid_tiles(tile) // ".halo.nc"
60  case ('leaf_area_index')
61  out_file = "./leaf_area_index." // grid_tiles(tile) // ".nc"
62  out_file_with_halo = "./leaf_area_index." // grid_tiles(tile) // ".halo.nc"
63  case ('visible_black_sky_albedo', 'visible_white_sky_albedo', 'near_IR_black_sky_albedo', 'near_IR_white_sky_albedo')
64  out_file = "./snowfree_albedo." // grid_tiles(tile) // ".nc"
65  out_file_with_halo = "./snowfree_albedo." // grid_tiles(tile) // ".halo.nc"
66  case ('facsf')
67  out_file = "./facsf." // grid_tiles(tile) // ".nc"
68  out_file_with_halo = "./facsf." // grid_tiles(tile) // ".halo.nc"
69  case ('slope_type')
70  out_file = "./slope_type." // grid_tiles(tile) // ".nc"
71  out_file_with_halo = "./slope_type." // grid_tiles(tile) // ".halo.nc"
72  case ('soil_type')
73  out_file = "./soil_type." // grid_tiles(tile) // ".nc"
74  out_file_with_halo = "./soil_type." // grid_tiles(tile) // ".halo.nc"
75  case ('vegetation_type')
76  out_file = "./vegetation_type." // grid_tiles(tile) // ".nc"
77  out_file_with_halo = "./vegetation_type." // grid_tiles(tile) // ".halo.nc"
78  case default
79  print*,'- FATAL ERROR IN ROUTINE OUTPUT. UNIDENTIFIED FIELD : ', field_names(field_idx)
80  call mpi_abort(mpi_comm_world, 67, ierr)
81  end select
82 
83 !----------------------------------------------------------------------
84 ! If user specified a halo (for running stand-alone regional grid),
85 ! remove it.
86 !----------------------------------------------------------------------
87 
88  if (halo > 0) then
89  print*,"- WILL REMOVE HALO REGION OF ", halo, " ROWS/COLS."
90  i_start = 1 + halo
91  i_end = i_mdl - halo
92  j_start = 1 + halo
93  j_end = j_mdl - halo
94  else
95  i_start = 1
96  i_end = i_mdl
97  j_start = 1
98  j_end = j_mdl
99  endif
100 
101  i_out = i_end - i_start + 1
102  j_out = j_end - j_start + 1
103 
104  if (record == 1) then
105 
106  initialsiz = 0
107  fsize = 65536
108  error = nf90_create(out_file, ior(nf90_netcdf4,nf90_classic_model), &
109  ncid(tile), initialsize=initialsiz, chunksize=fsize)
110  call netcdf_err(error, 'ERROR IN NF90_CREATE' )
111  error = nf90_def_dim(ncid(tile), 'nx', i_out, dim_x)
112  call netcdf_err(error, 'DEFINING NX DIMENSION' )
113  error = nf90_def_dim(ncid(tile), 'ny', j_out, dim_y)
114  call netcdf_err(error, 'DEFINING NY DIMENSION' )
115  error = nf90_def_dim(ncid(tile), 'time', num_time_recs, dim_time)
116  call netcdf_err(error, 'DEFINING TIME DIMENSION' )
117  error = nf90_def_var(ncid(tile), 'time', nf90_float, dim_time, id_times)
118  call netcdf_err(error, 'DEFINING TIME VARIABLE' )
119  error = nf90_put_att(ncid(tile), id_times, "units", "days since 2015-1-1")
120  call netcdf_err(error, 'DEFINING TIME ATTRIBUTE' )
121  if (len_trim(source) > 0) then
122  error = nf90_put_att(ncid(tile), nf90_global, 'source', source)
123  call netcdf_err(error, 'DEFINING GLOBAL SOURCE ATTRIBUTE' )
124  endif
125 
126  error = nf90_def_var(ncid(tile), 'geolat', nf90_float, (/dim_x,dim_y/), id_lat)
127  call netcdf_err(error, 'DEFINING GEOLAT FIELD' )
128  error = nf90_put_att(ncid(tile), id_lat, "long_name", "Latitude")
129  call netcdf_err(error, 'DEFINING GEOLAT NAME ATTRIBUTE' )
130  error = nf90_put_att(ncid(tile), id_lat, "units", "degrees_north")
131  call netcdf_err(error, 'DEFINING GEOLAT UNIT ATTRIBUTE' )
132  error = nf90_def_var(ncid(tile), 'geolon', nf90_float, (/dim_x,dim_y/), id_lon)
133  call netcdf_err(error, 'DEFINING GEOLON FIELD' )
134  error = nf90_put_att(ncid(tile), id_lon, "long_name", "Longitude")
135  call netcdf_err(error, 'DEFINING GEOLON NAME ATTRIBUTE' )
136  error = nf90_put_att(ncid(tile), id_lon, "units", "degrees_east")
137  call netcdf_err(error, 'DEFINING GEOLON UNIT ATTRIBUTE' )
138 
139  do j = 1, num_fields
140  error = nf90_def_var(ncid(tile), trim(field_names(j)), nf90_float, (/dim_x,dim_y,dim_time/), id_data)
141  call netcdf_err(error, 'DEFINING FIELD' )
142  error = nf90_put_att(ncid(tile), id_data, "missing_value", missing)
143  call netcdf_err(error, 'DEFINING FIELD ATTRIBUTE' )
144  error = nf90_put_att(ncid(tile), id_data, "coordinates", "geolon geolat")
145  call netcdf_err(error, 'DEFINING COORD ATTRIBUTE' )
146  enddo
147 
148  error = nf90_enddef(ncid(tile), header_buffer_val,4,0,4)
149  call netcdf_err(error, 'IN NF90_ENDDEF' )
150 
151  error = nf90_put_var( ncid(tile), id_times, day_of_rec)
152  call netcdf_err(error, 'WRITING TIME FIELD' )
153 
154  error = nf90_put_var( ncid(tile), id_lat, lat_one_tile(i_start:i_end,j_start:j_end), &
155  start=(/1,1/), count=(/i_out,j_out/))
156  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLAT' )
157 
158  error = nf90_put_var( ncid(tile), id_lon, lon_one_tile(i_start:i_end,j_start:j_end), &
159  start=(/1,1/), count=(/i_out,j_out/))
160  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLON' )
161 
162  endif
163 
164  print*,'- WRITE DATA FOR RECORD: ',record
165  error = nf90_inq_varid( ncid(tile), field_names(field_idx), id_data)
166  call netcdf_err(error, 'IN NF90_INQ_VARID' )
167  error = nf90_put_var( ncid(tile), id_data, data_one_tile(i_start:i_end,j_start:j_end), &
168  start=(/1,1,time/), count=(/i_out,j_out,1/))
169  call netcdf_err(error, 'IN NF90_PUT_VAR' )
170 
171  if (record == num_records) then
172  error = nf90_close(ncid(tile))
173  endif
174 
175 !----------------------------------------------------------------------
176 ! For regional nests, also output files including the halo
177 !----------------------------------------------------------------------
178 
179  if (halo == 0) return
180 
181  print*,"- WRITE OUT FILES THAT INCLUDE HALO REGION."
182 
183  if (record == 1) then
184 
185  initialsiz = 0
186  fsize = 65536
187  error = nf90_create(out_file_with_halo, ior(nf90_netcdf4,nf90_classic_model), &
188  ncid_with_halo, initialsize=initialsiz, chunksize=fsize)
189  call netcdf_err(error, 'IN NF90_CREATE' )
190  error = nf90_def_dim(ncid_with_halo, 'nx', i_mdl, dim_x)
191  call netcdf_err(error, 'DEFINING NX DIMENSION' )
192  error = nf90_def_dim(ncid_with_halo, 'ny', j_mdl, dim_y)
193  call netcdf_err(error, 'DEFINING NY DIMENSION' )
194  error = nf90_def_dim(ncid_with_halo, 'time', num_time_recs, dim_time)
195  call netcdf_err(error, 'DEFINING TIME DIMENSION' )
196  error = nf90_def_var(ncid_with_halo, 'time', nf90_float, dim_time, id_times)
197  call netcdf_err(error, 'DEFINING TIME VARIABLE' )
198  error = nf90_put_att(ncid_with_halo, id_times, "units", "days since 2015-1-1")
199  call netcdf_err(error, 'DEFINING TIME ATTRIBUTE' )
200  if (len_trim(source) > 0) then
201  error = nf90_put_att(ncid_with_halo, nf90_global, 'source', source)
202  call netcdf_err(error, 'DEFINING GLOBAL SOURCE ATTRIBUTE' )
203  endif
204 
205  error = nf90_def_var(ncid_with_halo, 'geolat', nf90_float, (/dim_x,dim_y/), id_lat)
206  call netcdf_err(error, 'DEFINING GEOLAT FIELD' )
207  error = nf90_put_att(ncid_with_halo, id_lat, "long_name", "Latitude")
208  call netcdf_err(error, 'DEFINING GEOLAT NAME ATTRIBUTE' )
209  error = nf90_put_att(ncid_with_halo, id_lat, "units", "degrees_north")
210  call netcdf_err(error, 'DEFINING GEOLAT UNIT ATTRIBUTE' )
211  error = nf90_def_var(ncid_with_halo, 'geolon', nf90_float, (/dim_x,dim_y/), id_lon)
212  call netcdf_err(error, 'DEFINING GEOLON FIELD' )
213  error = nf90_put_att(ncid_with_halo, id_lon, "long_name", "Longitude")
214  call netcdf_err(error, 'DEFINING GEOLON NAME ATTRIBUTE' )
215  error = nf90_put_att(ncid_with_halo, id_lon, "units", "degrees_east")
216  call netcdf_err(error, 'DEFINING GEOLON UNIT ATTRIBUTE' )
217 
218  do j = 1, num_fields
219  error = nf90_def_var(ncid_with_halo, field_names(j), nf90_float, (/dim_x,dim_y,dim_time/), id_data)
220  call netcdf_err(error, 'DEFINING FIELD VARIABLE' )
221  error = nf90_put_att(ncid_with_halo, id_data, "missing_value", missing)
222  call netcdf_err(error, 'DEFINING FIELD ATTRIBUTE' )
223  error = nf90_put_att(ncid_with_halo, id_data, "coordinates", "geolon geolat")
224  call netcdf_err(error, 'DEFINING COORD ATTRIBUTE' )
225  enddo
226 
227  error = nf90_enddef(ncid_with_halo, header_buffer_val,4,0,4)
228  call netcdf_err(error, 'WRITING HEADER ENDDEF' )
229 
230  error = nf90_put_var(ncid_with_halo, id_times, day_of_rec)
231  call netcdf_err(error, 'WRITING TIME VARIABLE' )
232 
233  error = nf90_put_var( ncid_with_halo, id_lat, lat_one_tile, &
234  start=(/1,1/), count=(/i_mdl,j_mdl/))
235  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLAT' )
236 
237  error = nf90_put_var( ncid_with_halo, id_lon, lon_one_tile, &
238  start=(/1,1/), count=(/i_mdl,j_mdl/))
239  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLON' )
240 
241  endif
242 
243  print*,'- WRITE DATA FOR RECORD: ',record
244  error = nf90_inq_varid(ncid_with_halo, field_names(field_idx), id_data)
245  call netcdf_err(error, 'IN NF90_INQ_VARID' )
246 
247  error = nf90_put_var(ncid_with_halo, id_data, data_one_tile, &
248  start=(/1,1,time/), count=(/i_mdl,j_mdl,1/))
249  call netcdf_err(error, 'IN NF90_PUT_VAR' )
250 
251  if (record == num_records) then
252  error = nf90_close(ncid_with_halo)
253  endif
254 
255  return
256 
257  end subroutine output
Set up program execution.
subroutine output(data_one_tile, lat_one_tile, lon_one_tile, i_mdl, j_mdl, tile, record, time, field_idx)
Output model data for a single tile and a single record in netcdf format.
Definition: output.f90:20
character(len=5), dimension(:), allocatable, public grid_tiles
Array of model grid tile names.
Definition: model_grid.F90:18
character(len=50), dimension(:), allocatable, public field_names
Names of fields to be processed.
Definition: source_grid.F90:21
real(kind=4), public missing
Value assigned to undefined points (i.e., ocean points for a land field).
Definition: model_grid.F90:25
This module defines the model grid.
Definition: model_grid.F90:10
Read grid specs, date information and land/sea mask for the source data that will be interpolated to ...
Definition: source_grid.F90:12
integer, public num_time_recs
Number of time records.
Definition: source_grid.F90:27
integer, public halo
Number of row/cols defining the lateral boundary halo.
integer, dimension(:), allocatable, public day_of_rec
Day of each time record with respect to Jan 1.
Definition: source_grid.F90:31
integer, public num_records
Number of fields times time records.
Definition: source_grid.F90:26
integer, public num_fields
Number of fields in the file.
Definition: source_grid.F90:28
subroutine, public netcdf_err(err, string)
Handle netCDF error codes.
Definition: utils.f90:23
Utilities.
Definition: utils.f90:8
character(len=75), public source
Original source of the data.
Definition: source_grid.F90:22