WAVEWATCH III  beta 0.0.1
scrip_remap_read.f
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 !
3 ! This routine reads remapping information from files written
4 ! by remap_setup. If remapping in both directions are required,
5 ! two input files must be specified.
6 !
7 !-----------------------------------------------------------------------
8 !
9 ! CVS:$Id: remap_read.f,v 1.6 2000/04/19 21:56:26 pwjones Exp $
10 !
11 ! Copyright (c) 1997, 1998 the Regents of the University of
12 ! California.
13 !
14 ! This software and ancillary information (herein called software)
15 ! called SCRIP is made available under the terms described here.
16 ! The software has been approved for release with associated
17 ! LA-CC Number 98-45.
18 !
19 ! Unless otherwise indicated, this software has been authored
20 ! by an employee or employees of the University of California,
21 ! operator of the Los Alamos National Laboratory under Contract
22 ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
23 ! Government has rights to use, reproduce, and distribute this
24 ! software. The public may copy and use this software without
25 ! charge, provided that this Notice and any statement of authorship
26 ! are reproduced on all copies. Neither the Government nor the
27 ! University makes any warranty, express or implied, or assumes
28 ! any liability or responsibility for the use of this software.
29 !
30 ! If software is modified to produce derivative works, such modified
31 ! software should be clearly marked, so as not to confuse it with
32 ! the version available from Los Alamos National Laboratory.
33 !
34 ! This code has been modified from the version available from
35 ! Los Alamos National Laboratory, for the purpose of running it
36 ! within WW3.
37 !
38 !***********************************************************************
39 
41 
42 !-----------------------------------------------------------------------
43 !
44 ! contains routines for reading a remap file
45 !
46 !-----------------------------------------------------------------------
47 
48  use scrip_kindsmod ! defines common data types
49  use scrip_errormod ! SCRIP error handler
50  use scrip_netcdfmod ! netCDF error handler
51  use netcdf ! netCDF library
52 
53  use scrip_constants ! defines common scalar constants
54  use scrip_grids ! module containing grid information
55  use scrip_remap_vars ! module containing remap information
56 
57  implicit none
58 
59 !-----------------------------------------------------------------------
60 !
61 ! module variables
62 !
63 !-----------------------------------------------------------------------
64 
65 
66  character(SCRIP_charLength), private ::
67  & map_method ! character string for map_type
68  &, normalize_opt ! character string for normalization option
69  &, convention ! character string for output convention
70 
71 !-----------------------------------------------------------------------
72 !
73 ! various netCDF ids for files variables
74 !
75 !-----------------------------------------------------------------------
76 
77  integer (SCRIP_i4), private ::
78  & ncstat ! error flag for netCDF calls
79  &, nc_file_id ! id for netCDF file
80  &, nc_srcgrdsize_id ! id for source grid size
81  &, nc_dstgrdsize_id ! id for destination grid size
82  &, nc_srcgrdcorn_id ! id for number of source grid corners
83  &, nc_dstgrdcorn_id ! id for number of dest grid corners
84  &, nc_srcgrdrank_id ! id for source grid rank
85  &, nc_dstgrdrank_id ! id for dest grid rank
86  &, nc_numlinks_id ! id for number of links in mapping
87  &, nc_numwgts_id ! id for number of weights for mapping
88  &, nc_srcgrddims_id ! id for source grid dimensions
89  &, nc_dstgrddims_id ! id for dest grid dimensions,
90  &, nc_srcgrdcntrlat_id ! id for source grid center latitude
91  &, nc_dstgrdcntrlat_id ! id for dest grid center latitude
92  &, nc_srcgrdcntrlon_id ! id for source grid center longitude
93  &, nc_dstgrdcntrlon_id ! id for dest grid center longitude
94  &, nc_srcgrdimask_id ! id for source grid mask
95  &, nc_dstgrdimask_id ! id for dest grid mask
96  &, nc_srcgrdcrnrlat_id ! id for latitude of source grid corners
97  &, nc_srcgrdcrnrlon_id ! id for longitude of source grid corners
98  &, nc_dstgrdcrnrlat_id ! id for latitude of dest grid corners
99  &, nc_dstgrdcrnrlon_id ! id for longitude of dest grid corners
100  &, nc_srcgrdarea_id ! id for area of source grid cells
101  &, nc_dstgrdarea_id ! id for area of dest grid cells
102  &, nc_srcgrdfrac_id ! id for area fraction on source grid
103  &, nc_dstgrdfrac_id ! id for area fraction on dest grid
104  &, nc_srcgrdadd_id ! id for map source address
105  &, nc_dstgrdadd_id ! id for map dest address
106  &, nc_rmpmatrix_id ! id for remapping matrix
107 
108 !***********************************************************************
109 
110  contains
111 
112 !***********************************************************************
113 
114  subroutine read_remap_ww3(map_name, interp_file, errorCode)
115 
116 !-----------------------------------------------------------------------
117 !
118 ! this driver routine reads some global attributes and then
119 ! calls a specific read routine based on file conventions
120 !
121 !-----------------------------------------------------------------------
122 
123 !-----------------------------------------------------------------------
124 !
125 ! input variables
126 !
127 !-----------------------------------------------------------------------
128 
129  character(SCRIP_CharLength), intent(in) ::
130  & interp_file ! filename for remap data
131 
132 !-----------------------------------------------------------------------
133 !
134 ! output variables
135 !
136 !-----------------------------------------------------------------------
137 
138  character(SCRIP_CharLength), intent(out) ::
139  & map_name ! name for mapping
140 
141 !-----------------------------------------------------------------------
142 !
143 ! local variables
144 !
145 !-----------------------------------------------------------------------
146 
147  integer (SCRIP_i4) ::
148  & errorCode ! error code for SCRIP routine
149 
150  character (14), parameter ::
151  & rtnName = 'read_remap_ww3'
152 
153 !-----------------------------------------------------------------------
154 !
155 ! open file and read some global information
156 !
157 !-----------------------------------------------------------------------
158 
159  ncstat = nf90_open(interp_file, nf90_noclobber, nc_file_id)
160  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
161  & opening interpolation file')) return
162 
163  !***
164  !*** map name
165  !***
166  map_name = ' '
167  ncstat = nf90_get_att(nc_file_id, nf90_global, 'title',
168  & map_name)
169  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
170  & reading map name')) return
171 
172 
173  !***
174  !*** normalization option
175  !***
176  normalize_opt = ' '
177  ncstat = nf90_get_att(nc_file_id, nf90_global, 'normalization',
178  & normalize_opt)
179  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
180  & reading interpolation option')) return
181 
182  select case(normalize_opt)
183  case ('none')
185  case ('fracarea')
187  case ('destarea')
189  case default
190  print *,'normalize_opt = ',normalize_opt
191  stop 'Invalid normalization option'
192  end select
193 
194  !***
195  !*** map method
196  !***
197  map_method = ' '
198  ncstat = nf90_get_att(nc_file_id, nf90_global, 'map_method',
199  & map_method)
200  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
201  & reading map method')) return
202 
203  select case(map_method)
204  case('Conservative remapping')
206  case('Bilinear remapping')
208  case('Distance weighted avg of nearest neighbors')
210  case('Bicubic remapping')
212  case default
213  print *,'map_type = ',map_method
214  stop 'Invalid Map Type'
215  end select
216 
217  !***
218  !*** file convention
219  !***
220  convention = ' '
221  ncstat = nf90_get_att(nc_file_id, nf90_global, 'conventions',
222  & convention)
223  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
224  & reading file convention')) return
225 
226 !-----------------------------------------------------------------------
227 !
228 ! call appropriate read routine based on output convention
229 !
230 !-----------------------------------------------------------------------
231 
232  select case(convention)
233  case ('SCRIP')
235  case ('NCAR-CSM')
236  print *,'convention = NCAR-CSM not supported for WW3'
237  stop 'unsupported file convention'
238 ! call read_remap_csm
239  case default
240  print *,'convention = ',convention
241  stop 'unknown output file convention'
242  end select
243 
244 !-----------------------------------------------------------------------
245 
246  end subroutine read_remap_ww3
247 
248 !***********************************************************************
249 
250  subroutine read_remap_scrip_ww3
251 
252 !-----------------------------------------------------------------------
253 !
254 ! the routine reads a netCDF file to extract remapping info
255 ! in SCRIP format
256 !
257 ! Only read variables needed by WW3
258 !
259 !-----------------------------------------------------------------------
260 !-----------------------------------------------------------------------
261 !
262 ! local variables
263 !
264 !-----------------------------------------------------------------------
265 
266  character (SCRIP_charLength) ::
267  & grid1_name ! grid name for source grid
268  &, grid2_name ! grid name for dest grid
269 
270  integer (SCRIP_i4) ::
271  & n ! dummy index
272  &, errorCode ! error code for SCRIP routine
273 
274  integer (SCRIP_i4), dimension(:), allocatable ::
275  & grid1_mask_int, ! integer masks to determine
276  & grid2_mask_int ! cells that participate in map
277 
278  character (20), parameter ::
279  & rtnName = 'read_remap_scrip_ww3'
280 
281 !-----------------------------------------------------------------------
282 !
283 ! read some additional global attributes
284 !
285 !-----------------------------------------------------------------------
286 
287  !***
288  !*** source and destination grid names
289  !***
290 
291  grid1_name = ' '
292  grid2_name = ' '
293  errorcode = scrip_success
294  ncstat = nf90_get_att(nc_file_id, nf90_global, 'source_grid',
295  & grid1_name)
296  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
297  & reading source grid name')) return
298 
299  ncstat = nf90_get_att(nc_file_id, nf90_global, 'dest_grid',
300  & grid2_name)
301  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
302  & reading destination grid name')) return
303 
304 ! Let's not include routine write statements w/out check on improc, l_master, etc.
305 ! print *,' '
306 ! print *,'Remapping between:',trim(grid1_name)
307 ! print *,'and ',trim(grid2_name)
308 ! print *,' '
309 
310 !-----------------------------------------------------------------------
311 !
312 ! read dimension information
313 !
314 !-----------------------------------------------------------------------
315 
316  ncstat = nf90_inq_dimid(nc_file_id, 'dst_grid_size',
317  & nc_dstgrdsize_id)
318  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
319  & reading destination grid id')) return
320 ! ncstat = nf90_inq_dimlen(nc_file_id, nc_dstgrdsize_id, grid2_size)
321  ncstat = nf90_inquire_dimension(nc_file_id, nc_dstgrdsize_id,
322  & len = grid2_size)
323  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
324  & reading destination grid size')) return
325 
326  ncstat = nf90_inq_dimid(nc_file_id, 'num_links',
327  & nc_numlinks_id)
328  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
329  & reading number of links id')) return
330 ! ncstat = nf90_inq_dimlen(nc_file_id, nc_numlinks_id,
331 ! & num_links_map1)
332  ncstat = nf90_inquire_dimension(nc_file_id, nc_numlinks_id,
333  & len = num_links_map1)
334  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
335  & reading number of links')) return
337 
338  ncstat = nf90_inq_dimid(nc_file_id, 'num_wgts',
339  & nc_numwgts_id)
340  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
341  & reading number of weights id')) return
342 ! ncstat = nf90_inq_dimlen(nc_file_id, nc_numwgts_id, num_wts)
343  ncstat = nf90_inquire_dimension(nc_file_id, nc_numwgts_id,
344  & len = num_wts)
345  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
346  & reading number of weights')) return
347 
348 !-----------------------------------------------------------------------
349 !
350 ! allocate arrays
351 !
352 !-----------------------------------------------------------------------
353 
354 ! grid2_frac is allocated in grid_init
355  if(allocated(grid2_frac))deallocate(grid2_frac)
356  allocate( grid2_frac(grid2_size) )
357 
358 ! grid1_add_map1, grid2_add_map1, wts_map1 are allocated in init_remap_vars
359  if(allocated(grid1_add_map1))deallocate(grid1_add_map1)
360  if(allocated(grid2_add_map1))deallocate(grid2_add_map1)
361  if(allocated(wts_map1))deallocate(wts_map1)
362  allocate( grid1_add_map1(num_links_map1),
365 
366 !-----------------------------------------------------------------------
367 !
368 ! get variable ids
369 !
370 !-----------------------------------------------------------------------
371 
372  ncstat = nf90_inq_varid(nc_file_id, 'dst_grid_frac',
373  & nc_dstgrdfrac_id)
374  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
375  & reading destination grid fraction id')) return
376 
377  ncstat = nf90_inq_varid(nc_file_id, 'src_address',
378  & nc_srcgrdadd_id)
379  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
380  & reading source grid address id')) return
381 
382  ncstat = nf90_inq_varid(nc_file_id, 'dst_address',
383  & nc_dstgrdadd_id)
384  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
385  & reading destination grid address id')) return
386 
387  ncstat = nf90_inq_varid(nc_file_id, 'remap_matrix',
388  & nc_rmpmatrix_id)
389  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
390  & reading remap matrix id')) return
391 
392 !-----------------------------------------------------------------------
393 !
394 ! read all variables
395 !
396 !-----------------------------------------------------------------------
397 
398  ncstat = nf90_get_var(nc_file_id, nc_dstgrdfrac_id,
399  & grid2_frac)
400  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
401  & reading destination grid fraction')) return
402 
403  ncstat = nf90_get_var(nc_file_id, nc_srcgrdadd_id,
404  & grid1_add_map1)
405  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
406  & reading source grid address')) return
407 
408  ncstat = nf90_get_var(nc_file_id, nc_dstgrdadd_id,
409  & grid2_add_map1)
410  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
411  & reading destination grid address')) return
412 
413  ncstat = nf90_get_var(nc_file_id, nc_rmpmatrix_id,
414  & wts_map1)
415  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
416  & reading remap weights')) return
417 
418 
419  ncstat = nf90_close(nc_file_id)
420  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname, 'error
421  & closing netCDF file')) return
422 
423 !-----------------------------------------------------------------------
424 
425  end subroutine read_remap_scrip_ww3
426 
427 !***********************************************************************
428 
429  end module scrip_remap_read
430 
431 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
scrip_remap_vars::norm_opt_none
integer(scrip_i4), parameter norm_opt_none
Definition: scrip_remap_vars.f:54
scrip_remap_vars::grid1_add_map1
integer(scrip_i4), dimension(:), allocatable, save grid1_add_map1
Definition: scrip_remap_vars.f:77
scrip_remap_read
Definition: scrip_remap_read.f:40
scrip_remap_vars::map_type_bilinear
integer(scrip_i4), parameter map_type_bilinear
Definition: scrip_remap_vars.f:59
scrip_netcdfmod
Definition: scrip_netcdfmod.f90:3
scrip_remap_vars::num_links_map1
integer(scrip_i4), save num_links_map1
Definition: scrip_remap_vars.f:66
scrip_remap_vars::norm_opt
integer(scrip_i4), save norm_opt
Definition: scrip_remap_vars.f:66
scrip_remap_vars::map_type
integer(scrip_i4), save map_type
Definition: scrip_remap_vars.f:66
scrip_grids
Definition: scrip_grids.f:49
scrip_remap_read::read_remap_scrip_ww3
subroutine read_remap_scrip_ww3
Definition: scrip_remap_read.f:251
scrip_remap_vars::norm_opt_frcarea
integer(scrip_i4), parameter norm_opt_frcarea
Definition: scrip_remap_vars.f:54
scrip_remap_vars::max_links_map1
integer(scrip_i4), save max_links_map1
Definition: scrip_remap_vars.f:66
scrip_errormod::scrip_success
integer(scrip_i4), parameter, public scrip_success
Definition: scrip_errormod.f90:42
scrip_grids::grid2_size
integer(scrip_i4), save grid2_size
Definition: scrip_grids.f:68
scrip_remap_vars::wts_map1
real(scrip_r8), dimension(:,:), allocatable, save wts_map1
Definition: scrip_remap_vars.f:83
scrip_constants
Definition: scrip_constants.f:38
scrip_remap_vars::map_type_distwgt
integer(scrip_i4), parameter map_type_distwgt
Definition: scrip_remap_vars.f:59
scrip_remap_read::read_remap_ww3
subroutine read_remap_ww3(map_name, interp_file, errorCode)
Definition: scrip_remap_read.f:115
scrip_kindsmod
Definition: scrip_kindsmod.f90:3
scrip_remap_vars::map_type_conserv
integer(scrip_i4), parameter map_type_conserv
Definition: scrip_remap_vars.f:59
scrip_errormod
Definition: scrip_errormod.f90:3
scrip_remap_vars::grid2_add_map1
integer(scrip_i4), dimension(:), allocatable, save grid2_add_map1
Definition: scrip_remap_vars.f:77
scrip_remap_vars::map_type_bicubic
integer(scrip_i4), parameter map_type_bicubic
Definition: scrip_remap_vars.f:59
scrip_remap_vars
Definition: scrip_remap_vars.f:40
scrip_netcdfmod::scrip_netcdferrorcheck
logical(scrip_logical) function, public scrip_netcdferrorcheck(netcdfStat, errorCode, rtnName, errorMsg)
Definition: scrip_netcdfmod.f90:45
scrip_remap_vars::norm_opt_dstarea
integer(scrip_i4), parameter norm_opt_dstarea
Definition: scrip_remap_vars.f:54
scrip_remap_vars::num_wts
integer(scrip_i4), save num_wts
Definition: scrip_remap_vars.f:66
scrip_grids::grid2_frac
real(scrip_r8), dimension(:), allocatable, target, save grid2_frac
Definition: scrip_grids.f:103