sfc_climo_gen  1.13.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 ('soil_color')
76  out_file = "./soil_color." // grid_tiles(tile) // ".nc"
77  out_file_with_halo = "./soil_color." // grid_tiles(tile) // ".halo.nc"
78  case ('vegetation_type')
79  out_file = "./vegetation_type." // grid_tiles(tile) // ".nc"
80  out_file_with_halo = "./vegetation_type." // grid_tiles(tile) // ".halo.nc"
81  case default
82  print*,'- FATAL ERROR IN ROUTINE OUTPUT. UNIDENTIFIED FIELD : ', field_names(field_idx)
83  call mpi_abort(mpi_comm_world, 67, ierr)
84  end select
85 
86 !----------------------------------------------------------------------
87 ! If user specified a halo (for running stand-alone regional grid),
88 ! remove it.
89 !----------------------------------------------------------------------
90 
91  if (halo > 0) then
92  print*,"- WILL REMOVE HALO REGION OF ", halo, " ROWS/COLS."
93  i_start = 1 + halo
94  i_end = i_mdl - halo
95  j_start = 1 + halo
96  j_end = j_mdl - halo
97  else
98  i_start = 1
99  i_end = i_mdl
100  j_start = 1
101  j_end = j_mdl
102  endif
103 
104  i_out = i_end - i_start + 1
105  j_out = j_end - j_start + 1
106 
107  if (record == 1) then
108 
109  initialsiz = 0
110  fsize = 65536
111  error = nf90_create(out_file, ior(nf90_netcdf4,nf90_classic_model), &
112  ncid(tile), initialsize=initialsiz, chunksize=fsize)
113  call netcdf_err(error, 'ERROR IN NF90_CREATE' )
114  error = nf90_def_dim(ncid(tile), 'nx', i_out, dim_x)
115  call netcdf_err(error, 'DEFINING NX DIMENSION' )
116  error = nf90_def_dim(ncid(tile), 'ny', j_out, dim_y)
117  call netcdf_err(error, 'DEFINING NY DIMENSION' )
118  error = nf90_def_dim(ncid(tile), 'time', num_time_recs, dim_time)
119  call netcdf_err(error, 'DEFINING TIME DIMENSION' )
120  error = nf90_def_var(ncid(tile), 'time', nf90_float, dim_time, id_times)
121  call netcdf_err(error, 'DEFINING TIME VARIABLE' )
122  error = nf90_put_att(ncid(tile), id_times, "units", "days since 2015-1-1")
123  call netcdf_err(error, 'DEFINING TIME ATTRIBUTE' )
124  if (len_trim(source) > 0) then
125  error = nf90_put_att(ncid(tile), nf90_global, 'source', source)
126  call netcdf_err(error, 'DEFINING GLOBAL SOURCE ATTRIBUTE' )
127  endif
128 
129  error = nf90_def_var(ncid(tile), 'geolat', nf90_float, (/dim_x,dim_y/), id_lat)
130  call netcdf_err(error, 'DEFINING GEOLAT FIELD' )
131  error = nf90_put_att(ncid(tile), id_lat, "long_name", "Latitude")
132  call netcdf_err(error, 'DEFINING GEOLAT NAME ATTRIBUTE' )
133  error = nf90_put_att(ncid(tile), id_lat, "units", "degrees_north")
134  call netcdf_err(error, 'DEFINING GEOLAT UNIT ATTRIBUTE' )
135  error = nf90_def_var(ncid(tile), 'geolon', nf90_float, (/dim_x,dim_y/), id_lon)
136  call netcdf_err(error, 'DEFINING GEOLON FIELD' )
137  error = nf90_put_att(ncid(tile), id_lon, "long_name", "Longitude")
138  call netcdf_err(error, 'DEFINING GEOLON NAME ATTRIBUTE' )
139  error = nf90_put_att(ncid(tile), id_lon, "units", "degrees_east")
140  call netcdf_err(error, 'DEFINING GEOLON UNIT ATTRIBUTE' )
141 
142  do j = 1, num_fields
143  error = nf90_def_var(ncid(tile), trim(field_names(j)), nf90_float, (/dim_x,dim_y,dim_time/), id_data)
144  call netcdf_err(error, 'DEFINING FIELD' )
145  error = nf90_put_att(ncid(tile), id_data, "missing_value", missing)
146  call netcdf_err(error, 'DEFINING FIELD ATTRIBUTE' )
147  error = nf90_put_att(ncid(tile), id_data, "coordinates", "geolon geolat")
148  call netcdf_err(error, 'DEFINING COORD ATTRIBUTE' )
149  enddo
150 
151  error = nf90_enddef(ncid(tile), header_buffer_val,4,0,4)
152  call netcdf_err(error, 'IN NF90_ENDDEF' )
153 
154  error = nf90_put_var( ncid(tile), id_times, day_of_rec)
155  call netcdf_err(error, 'WRITING TIME FIELD' )
156 
157  error = nf90_put_var( ncid(tile), id_lat, lat_one_tile(i_start:i_end,j_start:j_end), &
158  start=(/1,1/), count=(/i_out,j_out/))
159  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLAT' )
160 
161  error = nf90_put_var( ncid(tile), id_lon, lon_one_tile(i_start:i_end,j_start:j_end), &
162  start=(/1,1/), count=(/i_out,j_out/))
163  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLON' )
164 
165  endif
166 
167  print*,'- WRITE DATA FOR RECORD: ',record
168  error = nf90_inq_varid( ncid(tile), field_names(field_idx), id_data)
169  call netcdf_err(error, 'IN NF90_INQ_VARID' )
170  error = nf90_put_var( ncid(tile), id_data, data_one_tile(i_start:i_end,j_start:j_end), &
171  start=(/1,1,time/), count=(/i_out,j_out,1/))
172  call netcdf_err(error, 'IN NF90_PUT_VAR' )
173 
174  if (record == num_records) then
175  error = nf90_close(ncid(tile))
176  endif
177 
178 !----------------------------------------------------------------------
179 ! For regional nests, also output files including the halo
180 !----------------------------------------------------------------------
181 
182  if (halo == 0) return
183 
184  print*,"- WRITE OUT FILES THAT INCLUDE HALO REGION."
185 
186  if (record == 1) then
187 
188  initialsiz = 0
189  fsize = 65536
190  error = nf90_create(out_file_with_halo, ior(nf90_netcdf4,nf90_classic_model), &
191  ncid_with_halo, initialsize=initialsiz, chunksize=fsize)
192  call netcdf_err(error, 'IN NF90_CREATE' )
193  error = nf90_def_dim(ncid_with_halo, 'nx', i_mdl, dim_x)
194  call netcdf_err(error, 'DEFINING NX DIMENSION' )
195  error = nf90_def_dim(ncid_with_halo, 'ny', j_mdl, dim_y)
196  call netcdf_err(error, 'DEFINING NY DIMENSION' )
197  error = nf90_def_dim(ncid_with_halo, 'time', num_time_recs, dim_time)
198  call netcdf_err(error, 'DEFINING TIME DIMENSION' )
199  error = nf90_def_var(ncid_with_halo, 'time', nf90_float, dim_time, id_times)
200  call netcdf_err(error, 'DEFINING TIME VARIABLE' )
201  error = nf90_put_att(ncid_with_halo, id_times, "units", "days since 2015-1-1")
202  call netcdf_err(error, 'DEFINING TIME ATTRIBUTE' )
203  if (len_trim(source) > 0) then
204  error = nf90_put_att(ncid_with_halo, nf90_global, 'source', source)
205  call netcdf_err(error, 'DEFINING GLOBAL SOURCE ATTRIBUTE' )
206  endif
207 
208  error = nf90_def_var(ncid_with_halo, 'geolat', nf90_float, (/dim_x,dim_y/), id_lat)
209  call netcdf_err(error, 'DEFINING GEOLAT FIELD' )
210  error = nf90_put_att(ncid_with_halo, id_lat, "long_name", "Latitude")
211  call netcdf_err(error, 'DEFINING GEOLAT NAME ATTRIBUTE' )
212  error = nf90_put_att(ncid_with_halo, id_lat, "units", "degrees_north")
213  call netcdf_err(error, 'DEFINING GEOLAT UNIT ATTRIBUTE' )
214  error = nf90_def_var(ncid_with_halo, 'geolon', nf90_float, (/dim_x,dim_y/), id_lon)
215  call netcdf_err(error, 'DEFINING GEOLON FIELD' )
216  error = nf90_put_att(ncid_with_halo, id_lon, "long_name", "Longitude")
217  call netcdf_err(error, 'DEFINING GEOLON NAME ATTRIBUTE' )
218  error = nf90_put_att(ncid_with_halo, id_lon, "units", "degrees_east")
219  call netcdf_err(error, 'DEFINING GEOLON UNIT ATTRIBUTE' )
220 
221  do j = 1, num_fields
222  error = nf90_def_var(ncid_with_halo, field_names(j), nf90_float, (/dim_x,dim_y,dim_time/), id_data)
223  call netcdf_err(error, 'DEFINING FIELD VARIABLE' )
224  error = nf90_put_att(ncid_with_halo, id_data, "missing_value", missing)
225  call netcdf_err(error, 'DEFINING FIELD ATTRIBUTE' )
226  error = nf90_put_att(ncid_with_halo, id_data, "coordinates", "geolon geolat")
227  call netcdf_err(error, 'DEFINING COORD ATTRIBUTE' )
228  enddo
229 
230  error = nf90_enddef(ncid_with_halo, header_buffer_val,4,0,4)
231  call netcdf_err(error, 'WRITING HEADER ENDDEF' )
232 
233  error = nf90_put_var(ncid_with_halo, id_times, day_of_rec)
234  call netcdf_err(error, 'WRITING TIME VARIABLE' )
235 
236  error = nf90_put_var( ncid_with_halo, id_lat, lat_one_tile, &
237  start=(/1,1/), count=(/i_mdl,j_mdl/))
238  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLAT' )
239 
240  error = nf90_put_var( ncid_with_halo, id_lon, lon_one_tile, &
241  start=(/1,1/), count=(/i_mdl,j_mdl/))
242  call netcdf_err(error, 'IN NF90_PUT_VAR FOR GEOLON' )
243 
244  endif
245 
246  print*,'- WRITE DATA FOR RECORD: ',record
247  error = nf90_inq_varid(ncid_with_halo, field_names(field_idx), id_data)
248  call netcdf_err(error, 'IN NF90_INQ_VARID' )
249 
250  error = nf90_put_var(ncid_with_halo, id_data, data_one_tile, &
251  start=(/1,1,time/), count=(/i_mdl,j_mdl,1/))
252  call netcdf_err(error, 'IN NF90_PUT_VAR' )
253 
254  if (record == num_records) then
255  error = nf90_close(ncid_with_halo)
256  endif
257 
258  return
259 
260  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