18 character(len=5),
allocatable,
public :: grid_tiles(:)
20 integer,
public :: i_mdl
21 integer,
public :: j_mdl
22 integer,
public :: ij_mdl
23 integer,
public :: num_tiles
25 real(kind=4),
public :: missing = -999.
28 type(esmf_grid
),
public :: grid_mdl
29 type(esmf_field),
public :: data_field_mdl
31 type(esmf_field),
public :: land_frac_field_mdl
37 type(esmf_field),
public :: mask_field_mdl
41 type(esmf_field),
public :: latitude_field_mdl
43 type(esmf_field),
public :: longitude_field_mdl
45 type(esmf_field),
public :: vegt_field_mdl
71 integer,
intent(in) :: localpet, npets
73 character(len=500) :: the_file
75 integer :: error, id_dim, id_tiles, ncid
76 integer :: id_grid_tiles, ierr
77 integer :: extra, rc, tile
78 integer,
allocatable :: decomptile(:,:)
80 integer(esmf_kind_i4),
allocatable :: mask_mdl_one_tile(:,:)
81 integer(esmf_kind_i4),
pointer :: mask_field_mdl_ptr(:,:)
82 integer(esmf_kind_i4),
pointer :: mask_mdl_ptr(:,:)
84 real(esmf_kind_r4),
allocatable :: latitude_one_tile(:,:)
85 real(esmf_kind_r4),
allocatable :: longitude_one_tile(:,:)
86 real(esmf_kind_r4),
allocatable :: land_frac_one_tile(:,:)
92 print*,
'- OPEN MODEL GRID MOSAIC FILE: ',trim(mosaic_file_mdl)
93 error=nf90_open(trim(mosaic_file_mdl),nf90_nowrite,ncid)
94 call
netcdf_err(error,
"OPENING MODEL GRID MOSAIC FILE")
96 print*,
"- READ NUMBER OF TILES"
97 error=nf90_inq_dimid(ncid,
'ntiles', id_tiles)
99 error=nf90_inquire_dimension(ncid,id_tiles,len=num_tiles)
101 error=nf90_inq_varid(ncid,
'gridtiles', id_grid_tiles)
102 call
netcdf_err(error,
"READING GRIDTILES ID")
103 allocate(grid_tiles(num_tiles))
105 print*,
"- READ TILE NAMES"
106 error=nf90_get_var(ncid, id_grid_tiles, grid_tiles)
109 error = nf90_close(ncid)
111 print*,
'- NUMBER OF TILES, MODEL GRID IS ', num_tiles
113 if (mod(npets,num_tiles) /= 0)
then
114 print*,
'- FATAL ERROR: MUST RUN THIS PROGRAM WITH A TASK COUNT THAT'
115 print*,
'- IS A MULTIPLE OF THE NUMBER OF TILES.'
116 call mpi_abort(mpi_comm_world, 44, ierr)
123 orog_dir_mdl = trim(orog_dir_mdl) //
'/'
124 the_file = trim(orog_dir_mdl) // trim(orog_files_mdl(1))
126 print*,
'- OPEN FIRST MODEL GRID OROGRAPHY FILE: ',trim(the_file)
127 error=nf90_open(trim(the_file),nf90_nowrite,ncid)
128 call
netcdf_err(error,
"OPENING MODEL GRID OROGRAPHY FILE")
129 print*,
"- READ GRID DIMENSIONS"
130 error=nf90_inq_dimid(ncid,
'lon', id_dim)
131 call
netcdf_err(error,
"READING MODEL LON ID")
132 error=nf90_inquire_dimension(ncid,id_dim,len=i_mdl)
134 error=nf90_inq_dimid(ncid,
'lat', id_dim)
135 call
netcdf_err(error,
"READING MODEL LAT ID")
136 error=nf90_inquire_dimension(ncid,id_dim,len=j_mdl)
138 error = nf90_close(ncid)
140 print*,
"- I/J DIMENSIONS OF THE MODEL GRID TILES ", i_mdl, j_mdl
142 ij_mdl = i_mdl * j_mdl
148 extra = npets / num_tiles
150 allocate(decomptile(2,num_tiles))
152 do tile = 1, num_tiles
153 decomptile(:,tile)=(/1,extra/)
156 print*,
"- CALL GridCreateMosaic FOR MODEL GRID"
157 grid_mdl = esmf_gridcreatemosaic(filename=trim(mosaic_file_mdl), &
158 regdecompptile=decomptile, &
159 staggerloclist=(/esmf_staggerloc_center, &
160 esmf_staggerloc_corner/), &
161 indexflag=esmf_index_global, &
162 tilefilepath=trim(orog_dir_mdl), &
164 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
167 print*,
"- CALL FieldCreate FOR DATA INTERPOLATED TO MODEL GRID."
168 data_field_mdl = esmf_fieldcreate(grid_mdl, &
169 typekind=esmf_typekind_r4, &
170 staggerloc=esmf_staggerloc_center, &
171 name=
"data on model grid", &
173 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
176 if (.not. fract_vegsoil_type)
then
177 print*,
"- CALL FieldCreate FOR VEGETATION TYPE INTERPOLATED TO MODEL GRID."
178 vegt_field_mdl = esmf_fieldcreate(grid_mdl, &
179 typekind=esmf_typekind_r4, &
180 staggerloc=esmf_staggerloc_center, &
181 name=
"veg type on model grid", &
183 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
187 print*,
"- CALL FieldCreate FOR MODEL GRID LATITUDE."
188 latitude_field_mdl = esmf_fieldcreate(grid_mdl, &
189 typekind=esmf_typekind_r4, &
190 staggerloc=esmf_staggerloc_center, &
191 name=
"latitude on model grid", &
193 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
196 print*,
"- CALL FieldCreate FOR MODEL GRID LONGITUDE."
197 longitude_field_mdl = esmf_fieldcreate(grid_mdl, &
198 typekind=esmf_typekind_r4, &
199 staggerloc=esmf_staggerloc_center, &
200 name=
"longitude on model grid", &
202 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
209 print*,
"- CALL FieldCreate FOR MODEL GRID LANDMASK."
210 mask_field_mdl = esmf_fieldcreate(grid_mdl, &
211 typekind=esmf_typekind_i4, &
212 staggerloc=esmf_staggerloc_center, &
213 name=
"model grid land mask", &
215 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
218 print*,
"- CALL FieldGet FOR MODEL GRID LANDMASK."
219 nullify(mask_field_mdl_ptr)
220 call esmf_fieldget(mask_field_mdl, &
221 farrayptr=mask_field_mdl_ptr, &
223 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
226 print*,
"- CALL FieldCreate FOR MODEL GRID LAND FRACTION."
227 land_frac_field_mdl = esmf_fieldcreate(grid_mdl, &
228 typekind=esmf_typekind_r4, &
229 staggerloc=esmf_staggerloc_center, &
230 name=
"model grid land fraction", &
232 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
235 if (localpet == 0)
then
236 allocate(mask_mdl_one_tile(i_mdl,j_mdl))
237 allocate(land_frac_one_tile(i_mdl,j_mdl))
238 allocate(latitude_one_tile(i_mdl,j_mdl))
239 allocate(longitude_one_tile(i_mdl,j_mdl))
241 allocate(mask_mdl_one_tile(0,0))
242 allocate(land_frac_one_tile(0,0))
243 allocate(latitude_one_tile(0,0))
244 allocate(longitude_one_tile(0,0))
247 do tile = 1, num_tiles
248 if (localpet == 0)
then
249 the_file = trim(orog_dir_mdl) // trim(orog_files_mdl(tile))
250 call
get_model_info(trim(the_file), mask_mdl_one_tile, land_frac_one_tile, &
251 latitude_one_tile, longitude_one_tile, i_mdl, j_mdl)
254 print*,
"- CALL FieldScatter FOR MODEL GRID MASK. TILE IS: ", tile
255 call esmf_fieldscatter(mask_field_mdl, mask_mdl_one_tile, rootpet=0, tile=tile, rc=rc)
256 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
259 print*,
"- CALL FieldScatter FOR MODEL GRID LAND FRACTION. TILE IS: ", tile
260 call esmf_fieldscatter(land_frac_field_mdl, land_frac_one_tile, rootpet=0, tile=tile, rc=rc)
261 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
264 print*,
"- CALL FieldScatter FOR MODEL LATITUDE. TILE IS: ", tile
265 call esmf_fieldscatter(latitude_field_mdl, latitude_one_tile, rootpet=0, tile=tile, rc=rc)
266 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
269 print*,
"- CALL FieldScatter FOR MODEL LONGITUDE. TILE IS: ", tile
270 call esmf_fieldscatter(longitude_field_mdl, longitude_one_tile, rootpet=0, tile=tile, rc=rc)
271 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
276 deallocate(mask_mdl_one_tile, latitude_one_tile, longitude_one_tile, land_frac_one_tile)
278 print*,
"- CALL GridAddItem FOR MODEL GRID."
279 call esmf_gridadditem(grid_mdl, &
280 itemflag=esmf_griditem_mask, &
281 staggerloc=esmf_staggerloc_center, &
283 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
286 print*,
"- CALL GridGetItem FOR MODEL GRID."
287 nullify(mask_mdl_ptr)
288 call esmf_gridgetitem(grid_mdl, &
289 itemflag=esmf_griditem_mask, &
290 farrayptr=mask_mdl_ptr, &
292 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
295 mask_mdl_ptr = mask_field_mdl_ptr
319 character(len=*),
intent(in) :: orog_file
321 integer,
intent(in) :: idim, jdim
322 integer(esmf_kind_i4),
intent(out) :: mask(idim,jdim)
324 real(esmf_kind_r4),
intent(out) :: lat2d(idim,jdim)
325 real(esmf_kind_r4),
intent(out) :: lon2d(idim,jdim)
326 real(esmf_kind_r4),
intent(out) :: land_frac(idim,jdim)
328 integer :: error, lat, lon, i, j
329 integer :: ncid, id_dim, id_var
331 real(kind=4),
allocatable :: dummy(:,:)
333 print*,
"- READ MODEL OROGRAPHY FILE"
335 print*,
'- OPEN FILE: ', orog_file
336 error=nf90_open(orog_file,nf90_nowrite,ncid)
337 call
netcdf_err(error,
"OPENING MODEL OROGRAPHY FILE")
339 print*,
"- READ I-DIMENSION"
340 error=nf90_inq_dimid(ncid,
'lon', id_dim)
342 error=nf90_inquire_dimension(ncid,id_dim,len=lon)
345 print*,
"- READ J-DIMENSION"
346 error=nf90_inq_dimid(ncid,
'lat', id_dim)
348 error=nf90_inquire_dimension(ncid,id_dim,len=lat)
351 print*,
"- I/J DIMENSIONS: ", lon, lat
353 if ((lon /= idim) .or. (lat /= jdim))
then
357 allocate(dummy(idim,jdim))
373 error=nf90_inq_varid(ncid,
'lake_frac', id_var)
375 print*,
"- READ LAND MASK (SLMSK)"
376 error=nf90_inq_varid(ncid,
'slmsk', id_var)
378 error=nf90_get_var(ncid, id_var, dummy)
383 print*,
"- READ LAND FRACTION"
384 error=nf90_inq_varid(ncid,
'land_frac', id_var)
385 call
netcdf_err(error,
"READING LAND_FRAC ID")
386 error=nf90_get_var(ncid, id_var, land_frac)
391 if (land_frac(i,j) > 0.0)
then
398 print*,
"- READ LATITUDE"
399 error=nf90_inq_varid(ncid,
'geolat', id_var)
401 error=nf90_get_var(ncid, id_var, dummy)
405 print*,
"- READ LONGITUDE"
406 error=nf90_inq_varid(ncid,
'geolon', id_var)
408 error=nf90_get_var(ncid, id_var, dummy)
412 error = nf90_close(ncid)
429 print*,
"- CALL GridDestroy FOR MODEL GRID."
430 call esmf_griddestroy(grid_mdl,rc=rc)
432 print*,
"- CALL FieldDestroy FOR MODEL GRID LAND MASK."
433 call esmf_fielddestroy(mask_field_mdl,rc=rc)
435 print*,
"- CALL FieldDestroy FOR MODEL GRID LAND FRACTION."
436 call esmf_fielddestroy(land_frac_field_mdl,rc=rc)
438 print*,
"- CALL FieldDestroy FOR MODEL GRID DATA FIELD."
439 call esmf_fielddestroy(data_field_mdl,rc=rc)
441 if (esmf_fieldiscreated(vegt_field_mdl))
then
442 print*,
"- CALL FieldDestroy FOR MODEL GRID VEGETATION TYPE."
443 call esmf_fielddestroy(vegt_field_mdl,rc=rc)
446 print*,
"- CALL FieldDestroy FOR MODEL GRID LATITUDE."
447 call esmf_fielddestroy(latitude_field_mdl,rc=rc)
449 print*,
"- CALL FieldDestroy FOR MODEL GRID LONGITUDE."
450 call esmf_fielddestroy(longitude_field_mdl,rc=rc)
subroutine get_model_info(orog_file, mask, land_frac, lat2d, lon2d, idim, jdim)
Get model information.
subroutine, public model_grid_cleanup
Model grid cleanup.
subroutine, public define_model_grid(localpet, npets)
Define model grid.
This module defines the model grid.
subroutine, public netcdf_err(err, string)
Handle netCDF error codes.
subroutine, public error_handler(string, rc)
Handle errors.
Set up program execution.