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 :: mask_field_mdl
33 type(esmf_field),
public :: latitude_field_mdl
35 type(esmf_field),
public :: longitude_field_mdl
37 type(esmf_field),
public :: vegt_field_mdl
63 integer,
intent(in) :: localpet, npets
65 character(len=500) :: the_file
67 integer :: error, id_dim, id_tiles, ncid
68 integer :: id_grid_tiles, ierr
69 integer :: extra, rc, tile
70 integer,
allocatable :: decomptile(:,:)
72 integer(esmf_kind_i4),
allocatable :: mask_mdl_one_tile(:,:)
73 integer(esmf_kind_i4),
pointer :: mask_field_mdl_ptr(:,:)
74 integer(esmf_kind_i4),
pointer :: mask_mdl_ptr(:,:)
76 real(esmf_kind_r4),
allocatable :: latitude_one_tile(:,:)
77 real(esmf_kind_r4),
allocatable :: longitude_one_tile(:,:)
83 print*,
'- OPEN MODEL GRID MOSAIC FILE: ',trim(mosaic_file_mdl)
84 error=nf90_open(trim(mosaic_file_mdl),nf90_nowrite,ncid)
85 call
netcdf_err(error,
"OPENING MODEL GRID MOSAIC FILE")
87 print*,
"- READ NUMBER OF TILES"
88 error=nf90_inq_dimid(ncid,
'ntiles', id_tiles)
90 error=nf90_inquire_dimension(ncid,id_tiles,len=num_tiles)
92 error=nf90_inq_varid(ncid,
'gridtiles', id_grid_tiles)
94 allocate(grid_tiles(num_tiles))
96 print*,
"- READ TILE NAMES"
97 error=nf90_get_var(ncid, id_grid_tiles, grid_tiles)
100 error = nf90_close(ncid)
102 print*,
'- NUMBER OF TILES, MODEL GRID IS ', num_tiles
104 if (mod(npets,num_tiles) /= 0)
then
105 print*,
'- FATAL ERROR: MUST RUN THIS PROGRAM WITH A TASK COUNT THAT'
106 print*,
'- IS A MULTIPLE OF THE NUMBER OF TILES.'
107 call mpi_abort(mpi_comm_world, 44, ierr)
114 orog_dir_mdl = trim(orog_dir_mdl) //
'/'
115 the_file = trim(orog_dir_mdl) // trim(orog_files_mdl(1))
117 print*,
'- OPEN FIRST MODEL GRID OROGRAPHY FILE: ',trim(the_file)
118 error=nf90_open(trim(the_file),nf90_nowrite,ncid)
119 call
netcdf_err(error,
"OPENING MODEL GRID OROGRAPHY FILE")
120 print*,
"- READ GRID DIMENSIONS"
121 error=nf90_inq_dimid(ncid,
'lon', id_dim)
122 call
netcdf_err(error,
"READING MODEL LON ID")
123 error=nf90_inquire_dimension(ncid,id_dim,len=i_mdl)
125 error=nf90_inq_dimid(ncid,
'lat', id_dim)
126 call
netcdf_err(error,
"READING MODEL LAT ID")
127 error=nf90_inquire_dimension(ncid,id_dim,len=j_mdl)
129 error = nf90_close(ncid)
131 print*,
"- I/J DIMENSIONS OF THE MODEL GRID TILES ", i_mdl, j_mdl
133 ij_mdl = i_mdl * j_mdl
139 extra = npets / num_tiles
141 allocate(decomptile(2,num_tiles))
143 do tile = 1, num_tiles
144 decomptile(:,tile)=(/1,extra/)
147 print*,
"- CALL GridCreateMosaic FOR MODEL GRID"
148 grid_mdl = esmf_gridcreatemosaic(filename=trim(mosaic_file_mdl), &
149 regdecompptile=decomptile, &
150 staggerloclist=(/esmf_staggerloc_center, &
151 esmf_staggerloc_corner/), &
152 indexflag=esmf_index_global, &
153 tilefilepath=trim(orog_dir_mdl), &
155 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
158 print*,
"- CALL FieldCreate FOR DATA INTERPOLATED TO MODEL GRID."
159 data_field_mdl = esmf_fieldcreate(grid_mdl, &
160 typekind=esmf_typekind_r4, &
161 staggerloc=esmf_staggerloc_center, &
162 name=
"data on model grid", &
164 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
167 print*,
"- CALL FieldCreate FOR VEGETATION TYPE INTERPOLATED TO MODEL GRID."
168 vegt_field_mdl = esmf_fieldcreate(grid_mdl, &
169 typekind=esmf_typekind_r4, &
170 staggerloc=esmf_staggerloc_center, &
171 name=
"veg type on model grid", &
173 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
176 print*,
"- CALL FieldCreate FOR MODEL GRID LATITUDE."
177 latitude_field_mdl = esmf_fieldcreate(grid_mdl, &
178 typekind=esmf_typekind_r4, &
179 staggerloc=esmf_staggerloc_center, &
180 name=
"latitude on model grid", &
182 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
185 print*,
"- CALL FieldCreate FOR MODEL GRID LONGITUDE."
186 longitude_field_mdl = esmf_fieldcreate(grid_mdl, &
187 typekind=esmf_typekind_r4, &
188 staggerloc=esmf_staggerloc_center, &
189 name=
"longitude on model grid", &
191 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
198 print*,
"- CALL FieldCreate FOR MODEL GRID LANDMASK."
199 mask_field_mdl = esmf_fieldcreate(grid_mdl, &
200 typekind=esmf_typekind_i4, &
201 staggerloc=esmf_staggerloc_center, &
202 name=
"model grid land mask", &
204 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
207 print*,
"- CALL FieldGet FOR MODEL GRID LANDMASK."
208 nullify(mask_field_mdl_ptr)
209 call esmf_fieldget(mask_field_mdl, &
210 farrayptr=mask_field_mdl_ptr, &
212 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
215 if (localpet == 0)
then
216 allocate(mask_mdl_one_tile(i_mdl,j_mdl))
217 allocate(latitude_one_tile(i_mdl,j_mdl))
218 allocate(longitude_one_tile(i_mdl,j_mdl))
220 allocate(mask_mdl_one_tile(0,0))
221 allocate(latitude_one_tile(0,0))
222 allocate(longitude_one_tile(0,0))
225 do tile = 1, num_tiles
226 if (localpet == 0)
then
227 the_file = trim(orog_dir_mdl) // trim(orog_files_mdl(tile))
229 latitude_one_tile, longitude_one_tile, i_mdl, j_mdl)
232 print*,
"- CALL FieldScatter FOR MODEL GRID MASK. TILE IS: ", tile
233 call esmf_fieldscatter(mask_field_mdl, mask_mdl_one_tile, rootpet=0, tile=tile, rc=rc)
234 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
237 print*,
"- CALL FieldScatter FOR MODEL LATITUDE. TILE IS: ", tile
238 call esmf_fieldscatter(latitude_field_mdl, latitude_one_tile, rootpet=0, tile=tile, rc=rc)
239 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
242 print*,
"- CALL FieldScatter FOR MODEL LONGITUDE. TILE IS: ", tile
243 call esmf_fieldscatter(longitude_field_mdl, longitude_one_tile, rootpet=0, tile=tile, rc=rc)
244 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
249 deallocate(mask_mdl_one_tile, latitude_one_tile, longitude_one_tile)
251 print*,
"- CALL GridAddItem FOR MODEL GRID."
252 call esmf_gridadditem(grid_mdl, &
253 itemflag=esmf_griditem_mask, &
254 staggerloc=esmf_staggerloc_center, &
256 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
259 print*,
"- CALL GridGetItem FOR MODEL GRID."
260 nullify(mask_mdl_ptr)
261 call esmf_gridgetitem(grid_mdl, &
262 itemflag=esmf_griditem_mask, &
263 farrayptr=mask_mdl_ptr, &
265 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
268 mask_mdl_ptr = mask_field_mdl_ptr
291 character(len=*),
intent(in) :: orog_file
293 integer,
intent(in) :: idim, jdim
294 integer(esmf_kind_i4),
intent(out) :: mask(idim,jdim)
296 real(esmf_kind_r4),
intent(out) :: lat2d(idim,jdim)
297 real(esmf_kind_r4),
intent(out) :: lon2d(idim,jdim)
299 integer :: error, lat, lon, i, j
300 integer :: ncid, id_dim, id_var
302 real(kind=4),
allocatable :: dummy(:,:)
304 print*,
"- READ MODEL OROGRAPHY FILE"
306 print*,
'- OPEN FILE: ', orog_file
307 error=nf90_open(orog_file,nf90_nowrite,ncid)
308 call
netcdf_err(error,
"OPENING MODEL OROGRAPHY FILE")
310 print*,
"- READ I-DIMENSION"
311 error=nf90_inq_dimid(ncid,
'lon', id_dim)
313 error=nf90_inquire_dimension(ncid,id_dim,len=lon)
316 print*,
"- READ J-DIMENSION"
317 error=nf90_inq_dimid(ncid,
'lat', id_dim)
319 error=nf90_inquire_dimension(ncid,id_dim,len=lat)
322 print*,
"- I/J DIMENSIONS: ", lon, lat
324 if ((lon /= idim) .or. (lat /= jdim))
then
328 allocate(dummy(idim,jdim))
338 error=nf90_inq_varid(ncid,
'lake_frac', id_var)
340 print*,
"- READ LAND MASK (SLMSK)"
341 error=nf90_inq_varid(ncid,
'slmsk', id_var)
343 error=nf90_get_var(ncid, id_var, dummy)
347 print*,
"- READ LAND FRACTION"
348 error=nf90_inq_varid(ncid,
'land_frac', id_var)
349 call
netcdf_err(error,
"READING LAND_FRAC ID")
350 error=nf90_get_var(ncid, id_var, dummy)
355 if (dummy(i,j) > 0.0)
then
362 print*,
"- READ LATITUDE"
363 error=nf90_inq_varid(ncid,
'geolat', id_var)
365 error=nf90_get_var(ncid, id_var, dummy)
369 print*,
"- READ LONGITUDE"
370 error=nf90_inq_varid(ncid,
'geolon', id_var)
372 error=nf90_get_var(ncid, id_var, dummy)
376 error = nf90_close(ncid)
393 print*,
"- CALL GridDestroy FOR MODEL GRID."
394 call esmf_griddestroy(grid_mdl,rc=rc)
396 print*,
"- CALL FieldDestroy FOR MODEL GRID LAND MASK."
397 call esmf_fielddestroy(mask_field_mdl,rc=rc)
399 print*,
"- CALL FieldDestroy FOR MODEL GRID DATA FIELD."
400 call esmf_fielddestroy(data_field_mdl,rc=rc)
402 print*,
"- CALL FieldDestroy FOR MODEL GRID VEGETATION TYPE."
403 call esmf_fielddestroy(vegt_field_mdl,rc=rc)
405 print*,
"- CALL FieldDestroy FOR MODEL GRID LATITUDE."
406 call esmf_fielddestroy(latitude_field_mdl,rc=rc)
408 print*,
"- CALL FieldDestroy FOR MODEL GRID LONGITUDE."
409 call esmf_fielddestroy(longitude_field_mdl,rc=rc)
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.
subroutine get_model_info(orog_file, mask, lat2d, lon2d, idim, jdim)
Get model information.