grid_tools 1.14.0
Loading...
Searching...
No Matches
shave_nc.F90
Go to the documentation of this file.
1
4
17 program shave_nc
18 use netcdf
19 implicit none
20 integer,parameter :: kdbl=selected_real_kind(p=13,r=200)
21 character(len=255) :: filename_full ! the netcdf input file
22 character(len=255) :: filename_shaved ! the NetCDF file for shaved data
23 integer :: idim_compute,jdim_compute,halo
24 integer :: i_count_compute,j_count_compute &
25 ,i_count_super,j_count_super
26 integer :: i_start ! starting i of 2-D data with halo rows
27 integer :: j_start ! starting j of 2-D data with halo rows
28 integer :: i_count ! i extent of 2-D data with halo rows
29 integer :: j_count ! j extent of 2-D data with halo rows
30 integer :: n_count
31 integer :: n_shave ! number of rows to shave off full data to leave halo rows
32 integer :: n_start ! location where Start reading incoming character data
33
34 integer :: n,na,ncid_in,ncid_out,nd,ndims,ngatts &
35 ,nvars,unlimdimid
36 integer :: natts ! # of attributes
37 integer :: nctype ! type of the nth variable
38 integer :: dim_id,len_dim,len_x,len_y,var_id,xdim_id,xdim_id_out &
39 ,ydim_id,ydim_id_out
40 integer :: istat
41 integer,dimension(1:2) :: dimids=(/0,0/)
42 real,dimension(:) ,allocatable :: var_1d_with_halo
43 real,dimension(:,:),allocatable :: var_2d_with_halo
44 real(kind=kdbl),dimension(:,:),allocatable :: var_2d_dbl_with_halo
45! real*8,dimension(:,:),allocatable :: var_2d_dbl_with_halo
46 character(len=50) :: file,name_dim,name_xdim,name_ydim &
47 ,xdim,ydim
48 character(len=50) :: name_att ! the attribute's name
49 character(len=50) :: name_var ! the name of the nth variable
50 character(len=255) :: att=' '
51 character(len=255),dimension(:),allocatable :: var_1d_char
52
53!-----------------------------------------------------------------------
54!*** Read in the required compute dimensions, halo size, and filenames.
55!-----------------------------------------------------------------------
56 read(5,*)idim_compute,jdim_compute,halo,filename_full,filename_shaved
57 write(6,*)' id ',idim_compute,' jd ',jdim_compute,' halo ',halo
58 write(6,*)' fn_f ',trim(filename_full)
59 write(6,*)' fn_s ',trim(filename_shaved)
60 i_count_compute=idim_compute+2*halo
61 j_count_compute=jdim_compute+2*halo
62 i_count_super =2*i_count_compute
63 j_count_super =2*j_count_compute
64
65!-----------------------------------------------------------------------
66!*** Open the netcdf file with the incoming data to be shaved.
67!-----------------------------------------------------------------------
68 call check(nf90_open(filename_full,nf90_nowrite,ncid_in)) ! Open the netcdf file; get the file ID.
69 call check(nf90_inquire(ncid_in,ndims,nvars,ngatts,unlimdimid)) ! Find the number of variables in the file.
70
71!-----------------------------------------------------------------------
72!*** Create the NetCDF file that the shaved data will be written into.
73!*** Match the NetCDF format of the input file.
74!-----------------------------------------------------------------------
75 call check(nf90_create(filename_shaved & !
76 ,or(nf90_classic_model,nf90_netcdf4) & ! Create NetCDF file for shaved data.
77 ,ncid_out)) !
78
79!-----------------------------------------------------------------------
80!*** Replicate the dimensions from the input to the output file
81!*** but change values as needed to account for the shaving. We
82!*** know the grid file and orog file have given names for their
83!*** x and y dimensions so use those to adjust to shaved values.
84!*** NOTE: Gridpoints in the grid file are on the domain's
85!*** supergrid while points in the orog file are on the
86!*** compute grid.
87!-----------------------------------------------------------------------
88 do nd=1,ndims
89 call check(nf90_inquire_dimension(ncid_in,nd,name_dim,len_dim)) ! Get this dimension's name and value.
90 select case (name_dim)
91 case ('nx') !---
92 len_dim=i_count_super !
93 xdim=name_dim !
94 file='grid_file' !
95 case ('nxp') ! Used by the
96 len_dim=i_count_super+1 ! grid file.
97 case ('ny') !
98 len_dim=j_count_super !
99 ydim=name_dim !
100 case ('nyp') !
101 len_dim=j_count_super+1 !---
102 case ('lon') !---
103 len_dim=i_count_compute !
104 xdim=name_dim ! Used by the
105 file='orog_file' ! orog file.
106 case ('lat') !
107 len_dim=j_count_compute !
108 ydim=name_dim !---
109 end select
110 dim_id=nd
111 call check(nf90_def_dim(ncid_out,name_dim,len_dim,dim_id)) !-- Insert dimension into the output file.
112 enddo
113
114!-----------------------------------------------------------------------
115!*** The output file's variables must be defined while that file
116!*** is still in define mode. Loop through the variables in the
117!*** input file and define each of them in the output file.
118!-----------------------------------------------------------------------
119 do n=1,nvars
120 var_id=n
121 call check(nf90_inquire_variable(ncid_in,var_id,name_var,nctype & !-- name and type of this variable
122 ,ndims,dimids,natts)) !-- # of dimensions, ID, and attributes in this variable
123 if(ndims==1)then
124 call check(nf90_def_var(ncid_out,name_var,nctype,dimids(1),var_id)) !-- Define this 1-D variable in the output file.
125 elseif(ndims==2)then
126 call check(nf90_def_var(ncid_out,name_var,nctype,dimids,var_id)) !-- Define this 2-D variable in the output file.
127 endif
128
129!-----------------------------------------------------------------------
130!*** Copy this variable's attributes to the output file's variable.
131!-----------------------------------------------------------------------
132 if(natts>0)then
133 do na=1,natts
134 call check(nf90_inq_attname(ncid_in,var_id,na,name_att)) !-- Get the attribute's name and ID from input file.
135 call check(nf90_copy_att(ncid_in,var_id,name_att,ncid_out,var_id)) !-- Copy to output file.
136 enddo
137 endif
138 enddo
139
140!-----------------------------------------------------------------------
141!*** Copy the global attributes to the output file.
142!-----------------------------------------------------------------------
143 do n=1,ngatts
144 call check(nf90_inq_attname(ncid_in,nf90_global,n,name_att))
145 call check(nf90_copy_att(ncid_in,nf90_global,name_att,ncid_out,nf90_global))
146 enddo
147 call check(nf90_enddef(ncid_out)) !-- Put the output file into data mode.
148
149!-----------------------------------------------------------------------
150!*** Get the x and y extents of the incoming grid with extra rows
151!*** so we can find determine how many rows to shave off.
152!-----------------------------------------------------------------------
153 call check(nf90_inq_dimid(ncid_in,xdim,xdim_id)) !-- Find the ID of the x dimension.
154 call check(nf90_inq_dimid(ncid_in,ydim,ydim_id)) !-- Find the ID of the y dimension.
155 call check(nf90_inquire_dimension(ncid_in,xdim_id,name_xdim,len_x)) !-- Length of x dimension of vars in incoming file.
156 call check(nf90_inquire_dimension(ncid_in,ydim_id,name_ydim,len_y)) !-- Length of y dimension of vars in incoming file.
157 if(trim(file)=='orog_file')then
158 i_start=(len_x-idim_compute)/2-halo+1 !-- Starting i of 2-D data with halo rows on compute grid.
159 j_start=(len_y-jdim_compute)/2-halo+1 !-- Starting j of 2-D data with halo rows on compute grid.
160
161 elseif(trim(file)=='grid_file')then
162 i_start=(len_x-2*idim_compute)/2-2*halo+1 !-- Starting i of 2-D data with halo rows on supergrid.
163 j_start=(len_y-2*jdim_compute)/2-2*halo+1 !-- Starting j of 2-D data with halo rows on supergrid.
164 endif
165
166!-----------------------------------------------------------------------
167!*** We assume the # of extra rows on the incoming data is the same
168!*** in both x and y so the # of rows to shave off to leave the
169!*** halo rows is also the same in x and y. So consider only the
170!*** values from the x dimension.
171!-----------------------------------------------------------------------
172 n_shave=i_start-1 !-- # of rows to shave off full data to leave halo rows.
173
174!-----------------------------------------------------------------------
175!*** Now loop through all the variables in the input netcdf file,
176!*** read in the data excluding all extra rows except for halo rows,
177!*** and then write that out to the output file.
178!-----------------------------------------------------------------------
179 var_loop: do n=1,nvars
180 var_id=n
181 call check(nf90_inquire_variable(ncid_in,var_id,name_var,nctype & !-- The name and type of the nth variable
182 ,ndims,dimids,natts)) !-- The dimensions, ID, and attributes in the nth variable
183 call check(nf90_inquire_dimension(ncid_in,dimids(1),name_xdim,len_x)) !-- Get the length of the input 1st dimension.
184 if(ndims==2)then
185 call check(nf90_inquire_dimension(ncid_in,dimids(2),name_ydim,len_y)) !-- Get the length of the input y dimension.
186 endif
187
188!-------------------
189!*** 1-D variables
190!-------------------
191 if(ndims==1)then
192
193!---------------
194!*** Character
195!---------------
196 if(nctype==nf90_char)then
197 n_start=1 !-- Start reading incoming character data at this location.
198 n_count=len_x !-- Character data is not gridded so not shaved.
199 allocate(var_1d_char(1:n_count),stat=istat)
200 call check(nf90_get_var(ncid_in,var_id,var_1d_char(:) & !-- Fill the 1-D character variable.
201 ,start=(/n_start/) &
202 ,count=(/n_count/)))
203 call check(nf90_put_var(ncid_out,var_id,var_1d_char)) !-- Write out the 1-D character variable.
204 deallocate(var_1d_char)
205
206!---------------
207!*** Numerical
208!---------------
209 else
210 n_start=n_shave+1 !-- Start reading incoming data at this location.
211 n_count=len_dim-2*n_shave !-- # of datapoints to fill in the shaved 1-D variable.
212 allocate(var_1d_with_halo(1:n_count),stat=istat)
213 call check(nf90_get_var(ncid_in,var_id,var_1d_with_halo(:) & !-- Fill the shaved 1-D variable.
214 ,start=(/n_start/) &
215 ,count=(/n_count/)))
216 call check(nf90_put_var(ncid_out,var_id,var_1d_with_halo)) !-- Write out the shaved 1-D variable.
217 deallocate(var_1d_with_halo)
218 endif
219
220!-------------------
221!*** 2-D variables
222!-------------------
223 elseif(ndims==2)then
224 if(trim(file)=='orog_file')then
225 i_start=(len_x-idim_compute)/2-halo+1 !-- Starting i of 2-D data with halo rows on compute grid.
226 j_start=(len_y-jdim_compute)/2-halo+1 !-- Starting j of 2-D data with halo rows on compute grid.
227 i_count=i_count_compute !-- i extent of 2-D data with halo rows on compute grid.
228 j_count=j_count_compute !-- j extent of 2-D data with halo rows on compute grid.
229 elseif(trim(file)=='grid_file')then
230 i_start=(len_x-2*idim_compute)/2-2*halo+1 !-- Starting i of 2-D data with halo rows on supergrid.
231 j_start=(len_y-2*jdim_compute)/2-2*halo+1 !-- Starting j of 2-D data with halo rows on supergrid.
232 i_count=i_count_super !-- i extent of 2-D data with halo rows on supergrid.
233 j_count=j_count_super !-- j extent of 2-D data with halo rows on supergrid.
234 if(trim(name_xdim)=='nxp')then
235 i_count=i_count+1 !-- nxp is # of cell corners in x, not centers.
236 endif
237 if(trim(name_ydim)=='nyp')then
238 j_count=j_count+1 !-- nyp is # of cell corners in y, not centers.
239 endif
240 endif
241 if(nctype==nf90_float)then !-- Single precision real variables
242 allocate(var_2d_with_halo(i_count,j_count),stat=istat)
243 call check(nf90_get_var(ncid_in,var_id,var_2d_with_halo(:,:) & !-- Fill array with compute data plus halo rows.
244 ,start=(/i_start,j_start/) &
245 ,count=(/i_count,j_count/)))
246 call check(nf90_put_var(ncid_out,var_id,var_2d_with_halo)) !-- Write out the shaved 2-D single precision variable.
247 deallocate(var_2d_with_halo)
248 elseif(nctype==nf90_double)then !-- Double precision real variables
249 allocate(var_2d_dbl_with_halo(i_count,j_count),stat=istat)
250 call check(nf90_get_var(ncid_in,var_id,var_2d_dbl_with_halo(:,:) & !-- Fill array with compute data plus halo rows.
251 ,start=(/i_start,j_start/) &
252 ,count=(/i_count,j_count/)))
253 call check(nf90_put_var(ncid_out,var_id,var_2d_dbl_with_halo)) !-- Write out the shaved 2-D double precision variable.
254 deallocate(var_2d_dbl_with_halo)
255 endif
256 endif
257 enddo var_loop
258 call check(nf90_close(ncid_out))
259 call check(nf90_close(ncid_in))
260
261 contains
262
266 subroutine check(status)
267 integer,intent(in) :: status
268 if(status /= nf90_noerr) then
269 print *, trim(nf90_strerror(status))
270 stop "Stopped"
271 end if
272 end subroutine check
273 end program shave_nc
subroutine check(status)
Check results of netCDF call.
program shave_nc
The grid driver step in FV3 preprocessing generates a grid_tile file and an oro_tile file for the reg...
Definition shave_nc.F90:17