ROMS
Loading...
Searching...
No Matches
get_varcoords_mod Module Reference

Data Types

interface  get_varcoords
 

Functions/Subroutines

subroutine get_varcoords_nf90 (ng, model, ncname, ncid, ncvname, ncvarid, nx, ny, xmin, xmax, x, ymin, ymax, y, rectangular)
 
subroutine get_varcoords_pio (ng, model, ncname, piofile, ncvname, piovar, nx, ny, xmin, xmax, x, ymin, ymax, y, rectangular)
 

Function/Subroutine Documentation

◆ get_varcoords_nf90()

subroutine get_varcoords_mod::get_varcoords_nf90 ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
integer, intent(in) ncid,
character (len=*), intent(in) ncvname,
integer, intent(in) ncvarid,
integer, intent(in) nx,
integer, intent(in) ny,
real(r8), intent(out) xmin,
real(r8), intent(out) xmax,
real(r8), dimension(nx,ny), intent(out) x,
real(r8), intent(out) ymin,
real(r8), intent(out) ymax,
real(r8), dimension(nx,ny), intent(out) y,
logical, intent(out) rectangular )

Definition at line 60 of file get_varcoords.F.

64!***********************************************************************
65!
66 USE mod_netcdf
67!
68! Imported variable declarations.
69!
70 integer, intent(in) :: ng, model, ncid, ncvarid
71 integer, intent(in) :: Nx, Ny
72
73 character (len=*), intent(in) :: ncname
74 character (len=*), intent(in) :: ncvname
75
76 logical, intent(out) :: rectangular
77
78 real(r8), intent(out) :: Xmin, Xmax, Ymin, Ymax
79
80 real(r8), intent(out) :: X(Nx,Ny)
81 real(r8), intent(out) :: Y(Nx,Ny)
82!
83! Local variable declarations
84!
85 logical :: foundit
86!
87 integer :: i, ic, j, jc
88 integer :: alen, blank1, blank2, my_varid
89!
90 real(r8), dimension(Nx) :: Xwrk
91 real(r8), dimension(Ny) :: Ywrk
92!
93 character (len=20) :: Xname, Yname
94
95 character (len=*), parameter :: MyFile = &
96 & __FILE__//", get_varcoords_nf90"
97!
98 sourcefile=myfile
99!
100!-----------------------------------------------------------------------
101! Get coarse variable coordinates.
102!-----------------------------------------------------------------------
103!
104! Inquire the variable attributes.
105!
106 CALL netcdf_inq_var (ng, model, ncname, &
107 & ncid = ncid, &
108 & myvarname = ncvname)
109 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
110!
111! Get names of variable coordinates.
112!
113 foundit=.false.
114 rectangular=.false.
115
116 DO i=1,n_vatt
117 IF (trim(var_aname(i)).eq.'coordinates') THEN
118 alen=len_trim(var_achar(i))
119 blank1=index(var_achar(i)(1:alen),' ')
120 blank2=index(var_achar(i)(blank1+1:alen),' ')
121 xname=var_achar(i)(1:blank1-1)
122 IF (blank2.gt.0) THEN
123 yname=var_achar(i)(blank1+1:blank1+blank2-1)
124 ELSE
125 yname=var_achar(i)(blank1+1:alen)
126 END IF
127 foundit=.true.
128 EXIT
129 END IF
130 END DO
131 IF (.not.foundit) THEN
132 IF (master) WRITE (stdout,10) trim(ncvname), trim(ncname)
133 exit_flag=2
134 ioerror=0
135 END IF
136!
137! Read in X-coordinates.
138!
139 CALL netcdf_inq_var (ng, model, ncname, &
140 & ncid = ncid, &
141 & myvarname = xname, &
142 & varid = my_varid)
143 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
144!
145 IF (n_vdim.eq.1) THEN
146 CALL netcdf_get_fvar (ng, model, ncname, &
147 & xname, xwrk, &
148 & ncid = ncid, &
149 & start = (/1/), &
150 & total = (/nx/))
151 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
152
153 rectangular=.true.
154 jc=0
155 DO j=1,ny
156 DO i=1,nx
157 x(i,j)=xwrk(i)
158 END DO
159 END DO
160 ELSE
161 CALL netcdf_get_fvar (ng, model, ncname, &
162 & xname, x, &
163 & ncid = ncid, &
164 & start = (/1,1/), &
165 & total = (/nx,ny/))
166 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
167
168 jc=1
169 DO j=2,ny
170 IF (x(1,j).eq.x(1,1)) THEN
171 jc=jc+1
172 END IF
173 END DO
174 END IF
175!
176! Read in Y-coordinates.
177!
178 CALL netcdf_inq_var (ng, model, ncname, &
179 & ncid = ncid, &
180 & myvarname = yname, &
181 & varid = my_varid)
182 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
183!
184 IF (n_vdim.eq.1) THEN
185 CALL netcdf_get_fvar (ng, model, ncname, yname, &
186 & ywrk, &
187 & ncid = ncid, &
188 & start = (/1/), &
189 & total = (/ny/))
190 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
191
192 rectangular=.true.
193 ic=0
194 DO j=1,ny
195 DO i=1,nx
196 y(i,j)=ywrk(j)
197 END DO
198 END DO
199 ELSE
200 CALL netcdf_get_fvar (ng, model, ncname, &
201 & yname, y, &
202 & ncid = ncid, &
203 & start = (/1,1/), &
204 & total = (/nx,ny/))
205 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
206
207 ic=1
208 DO i=2,nx
209 IF (y(i,1).eq.y(1,1)) THEN
210 ic=ic+1
211 END IF
212 END DO
213 END IF
214!
215! Determine "rectangular" switch.
216!
217 IF (((ic.ne.0).and.(ic.eq.nx)).and. &
218 & ((jc.ne.0).and.(jc.eq.ny))) THEN
219 rectangular=.true.
220 END IF
221!
222! Determine minimum and maximum positions.
223!
224 xmin= spval
225 xmax=-spval
226 ymin= spval
227 ymax=-spval
228 DO j=1,ny
229 DO i=1,nx
230 xmin=min(xmin,x(i,j))
231 xmax=max(xmax,x(i,j))
232 ymin=min(ymin,y(i,j))
233 ymax=max(ymax,y(i,j))
234 END DO
235 END DO
236!
237 10 FORMAT (/,' GET_VARCOORDS_NF90 - Cannot find "coordinates" ', &
238 & 'attribute for variable:',2x,a,/,22x,'in file:',2x,a,/, &
239 & /,22x,'This attribute is needed to interpolate input data', &
240 & /,22x,'to model grid. Following CF compliance, we need:',/, &
241 & /,22x,'float my_var(time, lat, lon) ;', &
242 & /,22x,' my_var:long_name = "my variable long name" ;', &
243 & /,22x,' my_var:units = "my variable units" ;', &
244 & /,22x,' my_var:coordinates = "lon lat" ;', &
245 & /,22x,' my_var:time = "my_var_time" ;',/)
246!
247 RETURN
integer n_vatt
Definition mod_netcdf.F:174
character(len=1024), dimension(nvara) var_achar
Definition mod_netcdf.F:183
integer n_vdim
Definition mod_netcdf.F:173
character(len=100), dimension(nvara) var_aname
Definition mod_netcdf.F:181
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)

◆ get_varcoords_pio()

subroutine get_varcoords_mod::get_varcoords_pio ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
type (file_desc_t), intent(inout) piofile,
character (len=*), intent(in) ncvname,
type (my_vardesc), intent(inout) piovar,
integer, intent(in) nx,
integer, intent(in) ny,
real(r8), intent(out) xmin,
real(r8), intent(out) xmax,
real(r8), dimension(nx,ny), intent(out) x,
real(r8), intent(out) ymin,
real(r8), intent(out) ymax,
real(r8), dimension(nx,ny), intent(out) y,
logical, intent(out) rectangular )

Definition at line 253 of file get_varcoords.F.

257!***********************************************************************
258!
260!
261! Imported variable declarations.
262!
263 logical, intent(out) :: rectangular
264!
265 integer, intent(in) :: ng, model
266 integer, intent(in) :: Nx, Ny
267!
268 real(r8), intent(out) :: Xmin, Xmax, Ymin, Ymax
269
270 real(r8), intent(out) :: X(Nx,Ny)
271 real(r8), intent(out) :: Y(Nx,Ny)
272!
273 character (len=*), intent(in) :: ncname
274 character (len=*), intent(in) :: ncvname
275!
276 TYPE (File_desc_t), intent(inout) :: pioFile
277 TYPE (My_VarDesc), intent(inout) :: pioVar
278!
279! Local variable declarations
280!
281 logical :: foundit
282!
283 integer :: i, ic, j, jc
284 integer :: alen, blank1, blank2
285!
286 real(r8), dimension(Nx) :: Xwrk
287 real(r8), dimension(Ny) :: Ywrk
288!
289 character (len=20) :: Xname, Yname
290
291 character (len=*), parameter :: MyFile = &
292 & __FILE__//", get_varcoords_pio"
293!
294 TYPE (Var_desc_t) :: my_pioVar
295!
296 sourcefile=myfile
297!
298!-----------------------------------------------------------------------
299! Get coarse variable coordinates.
300!-----------------------------------------------------------------------
301!
302! Inquire the variable attributes.
303!
304 CALL pio_netcdf_inq_var (ng, model, ncname, &
305 & piofile = piofile, &
306 & myvarname = ncvname)
307 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
308!
309! Get names of variable coordinates.
310!
311 foundit=.false.
312 rectangular=.false.
313
314 DO i=1,n_vatt
315 IF (trim(var_aname(i)).eq.'coordinates') THEN
316 alen=len_trim(var_achar(i))
317 blank1=index(var_achar(i)(1:alen),' ')
318 blank2=index(var_achar(i)(blank1+1:alen),' ')
319 xname=var_achar(i)(1:blank1-1)
320 IF (blank2.gt.0) THEN
321 yname=var_achar(i)(blank1+1:blank1+blank2-1)
322 ELSE
323 yname=var_achar(i)(blank1+1:alen)
324 END IF
325 foundit=.true.
326 EXIT
327 END IF
328 END DO
329 IF (.not.foundit) THEN
330 IF (master) WRITE (stdout,10) trim(ncvname), trim(ncname)
331 exit_flag=2
332 ioerror=0
333 END IF
334!
335! Read in X-coordinates.
336!
337 CALL pio_netcdf_inq_var (ng, model, ncname, &
338 & piofile = piofile, &
339 & myvarname = xname, &
340 & piovar = my_piovar)
341 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
342!
343 IF (n_vdim.eq.1) THEN
344 CALL pio_netcdf_get_fvar (ng, model, ncname, &
345 & xname, xwrk, &
346 & piofile = piofile, &
347 & start = (/1/), &
348 & total = (/nx/))
349 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
350
351 rectangular=.true.
352 jc=0
353 DO j=1,ny
354 DO i=1,nx
355 x(i,j)=xwrk(i)
356 END DO
357 END DO
358 ELSE
359 CALL pio_netcdf_get_fvar (ng, model, ncname, &
360 & xname, x, &
361 & piofile = piofile, &
362 & start = (/1,1/), &
363 & total = (/nx,ny/))
364 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
365
366 jc=1
367 DO j=2,ny
368 IF (x(1,j).eq.x(1,1)) THEN
369 jc=jc+1
370 END IF
371 END DO
372 END IF
373!
374! Read in Y-coordinates.
375!
376 CALL pio_netcdf_inq_var (ng, model, ncname, &
377 & piofile = piofile, &
378 & myvarname = yname, &
379 & piovar = my_piovar)
380 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
381!
382 IF (n_vdim.eq.1) THEN
383 CALL pio_netcdf_get_fvar (ng, model, ncname, &
384 & yname, ywrk, &
385 & piofile = piofile, &
386 & start = (/1/), &
387 & total = (/ny/))
388 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
389
390 rectangular=.true.
391 ic=0
392 DO j=1,ny
393 DO i=1,nx
394 y(i,j)=ywrk(j)
395 END DO
396 END DO
397 ELSE
398 CALL pio_netcdf_get_fvar (ng, model, ncname, &
399 & yname, y, &
400 & piofile = piofile, &
401 & start = (/1,1/), &
402 & total = (/nx,ny/))
403 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
404
405 ic=1
406 DO i=2,nx
407 IF (y(i,1).eq.y(1,1)) THEN
408 ic=ic+1
409 END IF
410 END DO
411 END IF
412!
413! Determine "rectangular" switch.
414!
415 IF (((ic.ne.0).and.(ic.eq.nx)).and. &
416 & ((jc.ne.0).and.(jc.eq.ny))) THEN
417 rectangular=.true.
418 END IF
419!
420! Determine minimum and maximum positions.
421!
422 xmin= spval
423 xmax=-spval
424 ymin= spval
425 ymax=-spval
426 DO j=1,ny
427 DO i=1,nx
428 xmin=min(xmin,x(i,j))
429 xmax=max(xmax,x(i,j))
430 ymin=min(ymin,y(i,j))
431 ymax=max(ymax,y(i,j))
432 END DO
433 END DO
434!
435 10 FORMAT (/,' GET_VARCOORDS_PIO - Cannot find "coordinates" ', &
436 & 'attribute for variable:',2x,a,/,21x,'in file:',2x,a,/, &
437 & /,21x,'This attribute is needed to interpolate input data', &
438 & /,21x,'to model grid. Following CF compliance, we need:',/, &
439 & /,21x,'float my_var(time, lat, lon) ;', &
440 & /,21x,' my_var:long_name = "my variable long name" ;', &
441 & /,21x,' my_var:units = "my variable units" ;', &
442 & /,21x,' my_var:coordinates = "lon lat" ;', &
443 & /,21x,' my_var:time = "my_var_time" ;',/)
444!
445 RETURN
subroutine, public pio_netcdf_inq_var(ng, model, ncname, piofile, myvarname, searchvar, piovar, nvardim, nvaratt)