New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
remap_read.f90 in branches/UKMO/dev_r8183_GC_couple_pkg/NEMOGCM/TOOLS/WEIGHTS/src – NEMO

source: branches/UKMO/dev_r8183_GC_couple_pkg/NEMOGCM/TOOLS/WEIGHTS/src/remap_read.f90 @ 8730

Last change on this file since 8730 was 8730, checked in by dancopsey, 6 years ago

Cleared out SVN keywords.

File size: 36.6 KB
Line 
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$
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!***********************************************************************
35
36      module remap_read
37
38!-----------------------------------------------------------------------
39!
40!     contains routines for reading a remap file
41!
42!-----------------------------------------------------------------------
43
44      use kinds_mod     ! defines common data types
45      use constants     ! defines useful constants
46      use grids         ! includes all grid information
47      use netcdf_mod    ! module with netcdf vars and utilities
48      use remap_vars    ! module for all required remapping variables
49
50      implicit none
51
52!-----------------------------------------------------------------------
53!
54!     module variables
55!
56!-----------------------------------------------------------------------
57!-----------------------------------------------------------------------
58!
59!     various netCDF ids for files variables
60!
61!-----------------------------------------------------------------------
62
63      integer (kind=int_kind), private ::  & ! netCDF ids
64               ncstat, nc_file_id, &
65               nc_srcgrdsize_id, nc_dstgrdsize_id, &
66               nc_srcgrdcorn_id, nc_dstgrdcorn_id, &
67               nc_srcgrdrank_id, nc_dstgrdrank_id, &
68               nc_srcgrddims_id, nc_dstgrddims_id, &
69               nc_numlinks_id, nc_numwgts_id,  &
70               nc_srcgrdimask_id, nc_dstgrdimask_id, &
71               nc_srcgrdcntrlat_id, nc_srcgrdcntrlon_id, &
72               nc_srcgrdcrnrlat_id, nc_srcgrdcrnrlon_id, &
73               nc_srcgrdarea_id, nc_srcgrdfrac_id, &
74               nc_dstgrdcntrlat_id, nc_dstgrdcntrlon_id, &
75               nc_dstgrdcrnrlat_id, nc_dstgrdcrnrlon_id, &
76               nc_dstgrdarea_id, nc_dstgrdfrac_id, &
77               nc_srcgrdadd_id, nc_dstgrdadd_id, nc_rmpmatrix_id
78
79!***********************************************************************
80
81      contains
82
83!***********************************************************************
84
85      subroutine read_remap(map_name, interp_file)
86
87!-----------------------------------------------------------------------
88!
89!     this driver routine reads some global attributes and then
90!     calls a specific read routine based on file conventions
91!
92!-----------------------------------------------------------------------
93
94!-----------------------------------------------------------------------
95!
96!     input variables
97!
98!-----------------------------------------------------------------------
99
100      character(char_len), intent(in) :: &
101        interp_file        ! filename for remap data
102
103!-----------------------------------------------------------------------
104!
105!     output variables
106!
107!-----------------------------------------------------------------------
108
109      character(char_len), intent(out) :: &
110        map_name            ! name for mapping
111
112!-----------------------------------------------------------------------
113!
114!     local variables
115!
116!-----------------------------------------------------------------------
117
118      character(char_len) ::  &
119         map_method        & ! character string for map_type
120      ,  normalize_opt     & ! character string for normalization option
121      ,  convention       ! character string for output convention
122
123!-----------------------------------------------------------------------
124!
125!     open file and read some global information
126!
127!-----------------------------------------------------------------------
128
129      ncstat = nf_open(interp_file, NF_NOWRITE, nc_file_id)
130      call netcdf_error_handler(ncstat)
131
132      !***
133      !*** map name
134      !***
135      map_name = ' '
136      ncstat = nf_get_att_text(nc_file_id, NF_GLOBAL, 'title', &
137                               map_name)
138      call netcdf_error_handler(ncstat)
139
140      print *,'Reading remapping:',trim(map_name)
141      print *,'From file:',trim(interp_file)
142
143      !***
144      !*** normalization option
145      !***
146      normalize_opt = ' '
147      ncstat = nf_get_att_text(nc_file_id, NF_GLOBAL, 'normalization', &
148                               normalize_opt)
149      call netcdf_error_handler(ncstat)
150
151      select case(normalize_opt)
152      case ('none')
153        norm_opt = norm_opt_none
154      case ('fracarea')
155        norm_opt = norm_opt_frcarea
156      case ('destarea')
157        norm_opt = norm_opt_dstarea
158      case default
159        print *,'normalize_opt = ',normalize_opt
160        stop 'Invalid normalization option'
161      end select
162
163      !***
164      !*** map method
165      !***
166      map_method = ' '
167      ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'map_method', &
168                                map_method)
169      call netcdf_error_handler(ncstat)
170
171      select case(map_method)
172      case('Conservative remapping')
173        map_type = map_type_conserv
174      case('Bilinear remapping')
175        map_type = map_type_bilinear
176      case('Distance weighted avg of nearest neighbors')
177        map_type = map_type_distwgt
178      case('Bicubic remapping')
179        map_type = map_type_bicubic
180      case default
181        print *,'map_type = ',map_method
182        stop 'Invalid Map Type'
183      end select
184
185      !***
186      !*** file convention
187      !***
188      convention = ' '
189      ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'conventions', &
190                                convention)
191      call netcdf_error_handler(ncstat)
192
193!-----------------------------------------------------------------------
194!
195!     call appropriate read routine based on output convention
196!
197!-----------------------------------------------------------------------
198
199      select case(convention)
200      case ('SCRIP')
201        call read_remap_scrip
202      case ('NCAR-CSM')
203        call read_remap_csm
204      case default
205        print *,'convention = ',convention
206        stop 'unknown output file convention'
207      end select
208
209!-----------------------------------------------------------------------
210
211      end subroutine read_remap
212
213!***********************************************************************
214
215      subroutine read_remap_scrip
216
217!-----------------------------------------------------------------------
218!
219!     the routine reads a netCDF file to extract remapping info
220!     in SCRIP format
221!
222!-----------------------------------------------------------------------
223!-----------------------------------------------------------------------
224!
225!     local variables
226!
227!-----------------------------------------------------------------------
228
229      character (char_len) :: &
230        grid1_name            & ! grid name for source grid
231      , grid2_name           ! grid name for dest   grid
232
233      integer (kind=int_kind) ::   &
234        n                    ! dummy index
235
236      integer (kind=int_kind), dimension(:), allocatable :: &
237        grid1_mask_int,       & ! integer masks to determine
238        grid2_mask_int       ! cells that participate in map
239
240!-----------------------------------------------------------------------
241!
242!     read some additional global attributes
243!
244!-----------------------------------------------------------------------
245
246      !***
247      !*** source and destination grid names
248      !***
249
250      grid1_name = ' '
251      grid2_name = ' '
252      ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'source_grid', &
253                                grid1_name)
254      call netcdf_error_handler(ncstat)
255
256      ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'dest_grid', &
257                                grid2_name)
258      call netcdf_error_handler(ncstat)
259
260      print *,' '
261      print *,'Remapping between:',trim(grid1_name)
262      print *,'and ',trim(grid2_name)
263      print *,' '
264
265!-----------------------------------------------------------------------
266!
267!     read dimension information
268!
269!-----------------------------------------------------------------------
270
271      ncstat = nf_inq_dimid(nc_file_id, 'src_grid_size',  &
272                            nc_srcgrdsize_id)
273      call netcdf_error_handler(ncstat)
274      ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdsize_id, grid1_size)
275      call netcdf_error_handler(ncstat)
276
277      ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_size',  &
278                            nc_dstgrdsize_id)
279      call netcdf_error_handler(ncstat)
280      ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdsize_id, grid2_size)
281      call netcdf_error_handler(ncstat)
282
283      ncstat = nf_inq_dimid(nc_file_id, 'src_grid_corners',  &
284                            nc_srcgrdcorn_id)
285      call netcdf_error_handler(ncstat)
286      ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdcorn_id,  &
287                             grid1_corners)
288      call netcdf_error_handler(ncstat)
289
290      ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_corners',  &
291                            nc_dstgrdcorn_id)
292      call netcdf_error_handler(ncstat)
293      ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdcorn_id,  &
294                             grid2_corners)
295      call netcdf_error_handler(ncstat)
296
297      ncstat = nf_inq_dimid(nc_file_id, 'src_grid_rank',  &
298                            nc_srcgrdrank_id)
299      call netcdf_error_handler(ncstat)
300      ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdrank_id,  &
301                             grid1_rank)
302      call netcdf_error_handler(ncstat)
303
304      ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_rank',  &
305                            nc_dstgrdrank_id)
306      call netcdf_error_handler(ncstat)
307      ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdrank_id,  &
308                             grid2_rank)
309      call netcdf_error_handler(ncstat)
310
311      ncstat = nf_inq_dimid(nc_file_id, 'num_links',  &
312                            nc_numlinks_id)
313      call netcdf_error_handler(ncstat)
314      ncstat = nf_inq_dimlen(nc_file_id, nc_numlinks_id,  &
315                             num_links_map1)
316      call netcdf_error_handler(ncstat)
317
318      ncstat = nf_inq_dimid(nc_file_id, 'num_wgts',  &
319                            nc_numwgts_id)
320      call netcdf_error_handler(ncstat)
321      ncstat = nf_inq_dimlen(nc_file_id, nc_numwgts_id, num_wts)
322      call netcdf_error_handler(ncstat)
323
324!-----------------------------------------------------------------------
325!
326!     allocate arrays
327!
328!-----------------------------------------------------------------------
329
330      allocate( grid1_dims      (grid1_rank), &
331                grid1_center_lat(grid1_size),  &
332                grid1_center_lon(grid1_size), &
333                grid1_area      (grid1_size), &
334                grid1_frac      (grid1_size), &
335                grid1_mask      (grid1_size), &
336                grid1_mask_int  (grid1_size), &
337                grid1_corner_lat(grid1_corners, grid1_size), &
338                grid1_corner_lon(grid1_corners, grid1_size) )
339
340      allocate( grid2_dims      (grid2_rank), &
341                grid2_center_lat(grid2_size),  &
342                grid2_center_lon(grid2_size), &
343                grid2_area      (grid2_size), &
344                grid2_frac      (grid2_size), &
345                grid2_mask      (grid2_size), &
346                grid2_mask_int  (grid2_size), &
347                grid2_corner_lat(grid2_corners, grid2_size), &
348                grid2_corner_lon(grid2_corners, grid2_size) )
349
350      allocate( grid1_add_map1(num_links_map1), &
351                grid2_add_map1(num_links_map1), &
352                wts_map1(num_wts,num_links_map1) )
353
354!-----------------------------------------------------------------------
355!
356!     get variable ids
357!
358!-----------------------------------------------------------------------
359
360      ncstat = nf_inq_varid(nc_file_id, 'src_grid_dims',  &
361                            nc_srcgrddims_id)
362      call netcdf_error_handler(ncstat)
363
364      ncstat = nf_inq_varid(nc_file_id, 'src_grid_imask',  &
365                            nc_srcgrdimask_id)
366      call netcdf_error_handler(ncstat)
367
368      ncstat = nf_inq_varid(nc_file_id, 'src_grid_center_lat',  &
369                                         nc_srcgrdcntrlat_id)
370      call netcdf_error_handler(ncstat)
371
372      ncstat = nf_inq_varid(nc_file_id, 'src_grid_center_lon',  &
373                                         nc_srcgrdcntrlon_id)
374      call netcdf_error_handler(ncstat)
375
376      ncstat = nf_inq_varid(nc_file_id, 'src_grid_corner_lat',  &
377                                         nc_srcgrdcrnrlat_id)
378      call netcdf_error_handler(ncstat)
379
380      ncstat = nf_inq_varid(nc_file_id, 'src_grid_corner_lon',  &
381                                         nc_srcgrdcrnrlon_id)
382      call netcdf_error_handler(ncstat)
383
384      ncstat = nf_inq_varid(nc_file_id, 'src_grid_area',  &
385                                         nc_srcgrdarea_id)
386      call netcdf_error_handler(ncstat)
387
388      ncstat = nf_inq_varid(nc_file_id, 'src_grid_frac',  &
389                                         nc_srcgrdfrac_id)
390      call netcdf_error_handler(ncstat)
391
392      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_dims',  &
393                            nc_dstgrddims_id)
394      call netcdf_error_handler(ncstat)
395
396      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_imask',  &
397                            nc_dstgrdimask_id)
398      call netcdf_error_handler(ncstat)
399
400      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_center_lat',  &
401                                         nc_dstgrdcntrlat_id)
402      call netcdf_error_handler(ncstat)
403
404      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_center_lon',  &
405                                         nc_dstgrdcntrlon_id)
406      call netcdf_error_handler(ncstat)
407
408      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_corner_lat',  &
409                                         nc_dstgrdcrnrlat_id)
410      call netcdf_error_handler(ncstat)
411
412      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_corner_lon',  &
413                                         nc_dstgrdcrnrlon_id)
414      call netcdf_error_handler(ncstat)
415
416      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_area',  &
417                                         nc_dstgrdarea_id)
418      call netcdf_error_handler(ncstat)
419
420      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_frac',  &
421                                         nc_dstgrdfrac_id)
422      call netcdf_error_handler(ncstat)
423
424      ncstat = nf_inq_varid(nc_file_id, 'src_address',  &
425                                         nc_srcgrdadd_id)
426      call netcdf_error_handler(ncstat)
427
428      ncstat = nf_inq_varid(nc_file_id, 'dst_address',  &
429                                         nc_dstgrdadd_id)
430      call netcdf_error_handler(ncstat)
431
432      ncstat = nf_inq_varid(nc_file_id, 'remap_matrix',  &
433                                         nc_rmpmatrix_id)
434      call netcdf_error_handler(ncstat)
435
436!-----------------------------------------------------------------------
437!
438!     read all variables
439!
440!-----------------------------------------------------------------------
441
442      ncstat = nf_get_var_int(nc_file_id, nc_srcgrddims_id,  &
443                              grid1_dims)
444      call netcdf_error_handler(ncstat)
445
446      ncstat = nf_get_var_int(nc_file_id, nc_srcgrdimask_id,  &
447                              grid1_mask_int)
448      call netcdf_error_handler(ncstat)
449
450      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlat_id,  &
451                                             grid1_center_lat)
452      call netcdf_error_handler(ncstat)
453
454      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlon_id,  &
455                                             grid1_center_lon)
456      call netcdf_error_handler(ncstat)
457
458      grid1_units = ' '
459      ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcntrlat_id, 'units', &
460                               grid1_units)
461      call netcdf_error_handler(ncstat)
462
463      select case (grid1_units(1:7))
464      case ('degrees')
465        grid1_center_lat = grid1_center_lat*deg2rad
466        grid1_center_lon = grid1_center_lon*deg2rad
467      case ('radians')
468        !*** no conversion necessary
469      case default
470        print *,'unknown units supplied for grid1 center lat/lon: '
471        print *,'proceeding assuming radians'
472      end select
473
474      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlat_id,  &
475                                             grid1_corner_lat)
476      call netcdf_error_handler(ncstat)
477
478      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlon_id,  &
479                                             grid1_corner_lon)
480      call netcdf_error_handler(ncstat)
481
482      grid1_units = ' '
483      ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcrnrlat_id, 'units', &
484                               grid1_units)
485      call netcdf_error_handler(ncstat)
486
487      select case (grid1_units(1:7))
488      case ('degrees')
489        grid1_corner_lat = grid1_corner_lat*deg2rad
490        grid1_corner_lon = grid1_corner_lon*deg2rad
491      case ('radians')
492        !*** no conversion necessary
493      case default
494        print *,'unknown units supplied for grid1 corner lat/lon: '
495        print *,'proceeding assuming radians'
496      end select
497
498      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdarea_id,  &
499                                             grid1_area)
500      call netcdf_error_handler(ncstat)
501
502      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdfrac_id,  &
503                                             grid1_frac)
504      call netcdf_error_handler(ncstat)
505
506      ncstat = nf_get_var_int(nc_file_id, nc_dstgrddims_id,  &
507                              grid2_dims)
508      call netcdf_error_handler(ncstat)
509
510      ncstat = nf_get_var_int(nc_file_id, nc_dstgrdimask_id,  &
511                              grid2_mask_int)
512      call netcdf_error_handler(ncstat)
513
514      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlat_id,  &
515                                             grid2_center_lat)
516      call netcdf_error_handler(ncstat)
517
518      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlon_id,  &
519                                             grid2_center_lon)
520      call netcdf_error_handler(ncstat)
521
522      grid2_units = ' '
523      ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcntrlat_id, 'units', &
524                               grid2_units)
525      call netcdf_error_handler(ncstat)
526
527      select case (grid2_units(1:7))
528      case ('degrees')
529        grid2_center_lat = grid2_center_lat*deg2rad
530        grid2_center_lon = grid2_center_lon*deg2rad
531      case ('radians')
532        !*** no conversion necessary
533      case default
534        print *,'unknown units supplied for grid2 center lat/lon: '
535        print *,'proceeding assuming radians'
536      end select
537
538      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlat_id,  &
539                                             grid2_corner_lat)
540      call netcdf_error_handler(ncstat)
541
542      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlon_id,  &
543                                             grid2_corner_lon)
544      call netcdf_error_handler(ncstat)
545
546      grid2_units = ' '
547      ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcrnrlat_id, 'units', &
548                               grid2_units)
549      call netcdf_error_handler(ncstat)
550
551      select case (grid2_units(1:7))
552      case ('degrees')
553        grid2_corner_lat = grid2_corner_lat*deg2rad
554        grid2_corner_lon = grid2_corner_lon*deg2rad
555      case ('radians')
556        !*** no conversion necessary
557      case default
558        print *,'unknown units supplied for grid2 corner lat/lon: '
559        print *,'proceeding assuming radians'
560      end select
561
562      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdarea_id,  &
563                                             grid2_area)
564      call netcdf_error_handler(ncstat)
565
566      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdfrac_id,  &
567                                             grid2_frac)
568      call netcdf_error_handler(ncstat)
569
570      ncstat = nf_get_var_int(nc_file_id, nc_srcgrdadd_id,  &
571                              grid1_add_map1)
572      call netcdf_error_handler(ncstat)
573
574      ncstat = nf_get_var_int(nc_file_id, nc_dstgrdadd_id,  &
575                              grid2_add_map1)
576      call netcdf_error_handler(ncstat)
577
578      ncstat = nf_get_var_double(nc_file_id, nc_rmpmatrix_id,  &
579                                             wts_map1)
580      call netcdf_error_handler(ncstat)
581
582!-----------------------------------------------------------------------
583!
584!     initialize logical mask
585!
586!-----------------------------------------------------------------------
587
588      where (grid1_mask_int == 1)
589        grid1_mask = .true.
590      elsewhere
591        grid1_mask = .false.
592      endwhere
593      where (grid2_mask_int == 1)
594        grid2_mask = .true.
595      elsewhere
596        grid2_mask = .false.
597      endwhere
598      deallocate(grid1_mask_int, grid2_mask_int)
599
600!-----------------------------------------------------------------------
601!
602!     close input file
603!
604!-----------------------------------------------------------------------
605
606      ncstat = nf_close(nc_file_id)
607      call netcdf_error_handler(ncstat)
608
609!-----------------------------------------------------------------------
610
611      end subroutine read_remap_scrip
612
613!***********************************************************************
614
615      subroutine read_remap_csm
616
617!-----------------------------------------------------------------------
618!
619!     the routine reads a netCDF file to extract remapping info
620!     in NCAR-CSM format
621!
622!-----------------------------------------------------------------------
623!-----------------------------------------------------------------------
624!
625!     local variables
626!
627!-----------------------------------------------------------------------
628
629      character (char_len) :: &
630        grid1_name            & ! grid name for source grid
631      , grid2_name           ! grid name for dest   grid
632
633      integer (kind=int_kind) :: &
634        nc_numwgts1_id     & ! extra netCDF id for num_wgts > 1
635      , nc_rmpmatrix2_id  ! extra netCDF id for high-order remap matrix
636
637      real (kind=dbl_kind), dimension(:),allocatable :: &
638        wts1              ! CSM wants single array for 1st-order wts
639
640      real (kind=dbl_kind), dimension(:,:),allocatable :: &
641        wts2              ! write remaining weights in different array
642
643      integer (kind=int_kind) ::   &
644        n                    ! dummy index
645
646      integer (kind=int_kind), dimension(:), allocatable :: &
647        grid1_mask_int,       & ! integer masks to determine
648        grid2_mask_int       ! cells that participate in map
649
650!-----------------------------------------------------------------------
651!
652!     read some additional global attributes
653!
654!-----------------------------------------------------------------------
655
656      !***
657      !*** source and destination grid names
658      !***
659
660      grid1_name = ' '
661      grid2_name = ' '
662      ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'domain_a', &
663                                grid1_name)
664      call netcdf_error_handler(ncstat)
665
666      ncstat = nf_get_att_text (nc_file_id, NF_GLOBAL, 'domain_b', &
667                                grid2_name)
668      call netcdf_error_handler(ncstat)
669
670      print *,' '
671      print *,'Remapping between:',trim(grid1_name)
672      print *,'and ',trim(grid2_name)
673      print *,' '
674
675!-----------------------------------------------------------------------
676!
677!     read dimension information
678!
679!-----------------------------------------------------------------------
680
681      ncstat = nf_inq_dimid(nc_file_id, 'n_a', nc_srcgrdsize_id)
682      call netcdf_error_handler(ncstat)
683      ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdsize_id, grid1_size)
684      call netcdf_error_handler(ncstat)
685
686      ncstat = nf_inq_dimid(nc_file_id, 'n_b', nc_dstgrdsize_id)
687      call netcdf_error_handler(ncstat)
688      ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdsize_id, grid2_size)
689      call netcdf_error_handler(ncstat)
690
691      ncstat = nf_inq_dimid(nc_file_id, 'nv_a', nc_srcgrdcorn_id)
692      call netcdf_error_handler(ncstat)
693      ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdcorn_id,  &
694                             grid1_corners)
695      call netcdf_error_handler(ncstat)
696
697      ncstat = nf_inq_dimid(nc_file_id, 'nv_b', nc_dstgrdcorn_id)
698      call netcdf_error_handler(ncstat)
699      ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdcorn_id,  &
700                             grid2_corners)
701      call netcdf_error_handler(ncstat)
702
703      ncstat = nf_inq_dimid(nc_file_id, 'src_grid_rank',  &
704                            nc_srcgrdrank_id)
705      call netcdf_error_handler(ncstat)
706      ncstat = nf_inq_dimlen(nc_file_id, nc_srcgrdrank_id,  &
707                             grid1_rank)
708      call netcdf_error_handler(ncstat)
709
710      ncstat = nf_inq_dimid(nc_file_id, 'dst_grid_rank',  &
711                            nc_dstgrdrank_id)
712      call netcdf_error_handler(ncstat)
713      ncstat = nf_inq_dimlen(nc_file_id, nc_dstgrdrank_id,  &
714                             grid2_rank)
715      call netcdf_error_handler(ncstat)
716
717      ncstat = nf_inq_dimid(nc_file_id, 'n_s',  &
718                            nc_numlinks_id)
719      call netcdf_error_handler(ncstat)
720      ncstat = nf_inq_dimlen(nc_file_id, nc_numlinks_id,  &
721                             num_links_map1)
722      call netcdf_error_handler(ncstat)
723
724      ncstat = nf_inq_dimid(nc_file_id, 'num_wgts',  &
725                            nc_numwgts_id)
726      call netcdf_error_handler(ncstat)
727      ncstat = nf_inq_dimlen(nc_file_id, nc_numwgts_id, num_wts)
728      call netcdf_error_handler(ncstat)
729
730      if (num_wts > 1) then
731        ncstat = nf_inq_dimid(nc_file_id, 'num_wgts1',  &
732                              nc_numwgts1_id)
733        call netcdf_error_handler(ncstat)
734      endif
735
736!-----------------------------------------------------------------------
737!
738!     allocate arrays
739!
740!-----------------------------------------------------------------------
741
742      allocate( grid1_dims      (grid1_rank), &
743                grid1_center_lat(grid1_size),  &
744                grid1_center_lon(grid1_size), &
745                grid1_area      (grid1_size), &
746                grid1_frac      (grid1_size), &
747                grid1_mask      (grid1_size), &
748                grid1_mask_int  (grid1_size), &
749                grid1_corner_lat(grid1_corners, grid1_size), &
750                grid1_corner_lon(grid1_corners, grid1_size) )
751
752      allocate( grid2_dims      (grid2_rank), &
753                grid2_center_lat(grid2_size),  &
754                grid2_center_lon(grid2_size), &
755                grid2_area      (grid2_size), &
756                grid2_frac      (grid2_size), &
757                grid2_mask      (grid2_size), &
758                grid2_mask_int  (grid2_size), &
759                grid2_corner_lat(grid2_corners, grid2_size), &
760                grid2_corner_lon(grid2_corners, grid2_size) )
761
762      allocate( grid1_add_map1(num_links_map1), &
763                grid2_add_map1(num_links_map1), &
764                wts_map1(num_wts,num_links_map1), &
765                wts1(num_links_map1), &
766                wts2(num_wts-1,num_links_map1) )
767
768!-----------------------------------------------------------------------
769!
770!     get variable ids
771!
772!-----------------------------------------------------------------------
773
774      ncstat = nf_inq_varid(nc_file_id, 'src_grid_dims',  &
775                            nc_srcgrddims_id)
776      call netcdf_error_handler(ncstat)
777
778      ncstat = nf_inq_varid(nc_file_id, 'mask_a',  &
779                            nc_srcgrdimask_id)
780      call netcdf_error_handler(ncstat)
781
782      ncstat = nf_inq_varid(nc_file_id, 'yc_a', nc_srcgrdcntrlat_id)
783      call netcdf_error_handler(ncstat)
784
785      ncstat = nf_inq_varid(nc_file_id, 'xc_a', nc_srcgrdcntrlon_id)
786      call netcdf_error_handler(ncstat)
787
788      ncstat = nf_inq_varid(nc_file_id, 'yv_a', nc_srcgrdcrnrlat_id)
789      call netcdf_error_handler(ncstat)
790
791      ncstat = nf_inq_varid(nc_file_id, 'xv_a', nc_srcgrdcrnrlon_id)
792      call netcdf_error_handler(ncstat)
793
794      ncstat = nf_inq_varid(nc_file_id, 'area_a', nc_srcgrdarea_id)
795      call netcdf_error_handler(ncstat)
796
797      ncstat = nf_inq_varid(nc_file_id, 'frac_a', nc_srcgrdfrac_id)
798      call netcdf_error_handler(ncstat)
799
800      ncstat = nf_inq_varid(nc_file_id, 'dst_grid_dims',  &
801                            nc_dstgrddims_id)
802      call netcdf_error_handler(ncstat)
803
804      ncstat = nf_inq_varid(nc_file_id, 'mask_b',  &
805                            nc_dstgrdimask_id)
806      call netcdf_error_handler(ncstat)
807
808      ncstat = nf_inq_varid(nc_file_id, 'yc_b', nc_dstgrdcntrlat_id)
809      call netcdf_error_handler(ncstat)
810
811      ncstat = nf_inq_varid(nc_file_id, 'xc_b', nc_dstgrdcntrlon_id)
812      call netcdf_error_handler(ncstat)
813
814      ncstat = nf_inq_varid(nc_file_id, 'yv_b', nc_dstgrdcrnrlat_id)
815      call netcdf_error_handler(ncstat)
816
817      ncstat = nf_inq_varid(nc_file_id, 'xv_b', nc_dstgrdcrnrlon_id)
818      call netcdf_error_handler(ncstat)
819
820      ncstat = nf_inq_varid(nc_file_id, 'area_b', nc_dstgrdarea_id)
821      call netcdf_error_handler(ncstat)
822
823      ncstat = nf_inq_varid(nc_file_id, 'frac_b', nc_dstgrdfrac_id)
824      call netcdf_error_handler(ncstat)
825
826      ncstat = nf_inq_varid(nc_file_id, 'col', nc_srcgrdadd_id)
827      call netcdf_error_handler(ncstat)
828
829      ncstat = nf_inq_varid(nc_file_id, 'row', nc_dstgrdadd_id)
830      call netcdf_error_handler(ncstat)
831
832      ncstat = nf_inq_varid(nc_file_id, 'S', nc_rmpmatrix_id)
833      call netcdf_error_handler(ncstat)
834
835      if (num_wts > 1) then
836        ncstat = nf_inq_varid(nc_file_id, 'S2', nc_rmpmatrix2_id)
837        call netcdf_error_handler(ncstat)
838      endif
839
840!-----------------------------------------------------------------------
841!
842!     read all variables
843!
844!-----------------------------------------------------------------------
845
846      ncstat = nf_get_var_int(nc_file_id, nc_srcgrddims_id,  &
847                              grid1_dims)
848      call netcdf_error_handler(ncstat)
849
850      ncstat = nf_get_var_int(nc_file_id, nc_srcgrdimask_id,  &
851                              grid1_mask_int)
852      call netcdf_error_handler(ncstat)
853
854      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlat_id,  &
855                                             grid1_center_lat)
856      call netcdf_error_handler(ncstat)
857
858      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcntrlon_id,  &
859                                             grid1_center_lon)
860      call netcdf_error_handler(ncstat)
861
862      ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcntrlat_id, 'units', &
863                               grid1_units)
864      call netcdf_error_handler(ncstat)
865
866      select case (grid1_units(1:7))
867      case ('degrees')
868        grid1_center_lat = grid1_center_lat*deg2rad
869        grid1_center_lon = grid1_center_lon*deg2rad
870      case ('radians')
871        !*** no conversion necessary
872      case default
873        print *,'unknown units supplied for grid1 center lat/lon: '
874        print *,'proceeding assuming radians'
875      end select
876
877      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlat_id,  &
878                                             grid1_corner_lat)
879      call netcdf_error_handler(ncstat)
880
881      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdcrnrlon_id,  &
882                                             grid1_corner_lon)
883      call netcdf_error_handler(ncstat)
884
885      ncstat = nf_get_att_text(nc_file_id, nc_srcgrdcrnrlat_id, 'units', &
886                               grid1_units)
887      call netcdf_error_handler(ncstat)
888
889      select case (grid1_units(1:7))
890      case ('degrees')
891        grid1_corner_lat = grid1_corner_lat*deg2rad
892        grid1_corner_lon = grid1_corner_lon*deg2rad
893      case ('radians')
894        !*** no conversion necessary
895      case default
896        print *,'unknown units supplied for grid1 corner lat/lon: '
897        print *,'proceeding assuming radians'
898      end select
899
900      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdarea_id,  &
901                                             grid1_area)
902      call netcdf_error_handler(ncstat)
903
904      ncstat = nf_get_var_double(nc_file_id, nc_srcgrdfrac_id,  &
905                                             grid1_frac)
906      call netcdf_error_handler(ncstat)
907
908      ncstat = nf_get_var_int(nc_file_id, nc_dstgrddims_id,  &
909                              grid2_dims)
910      call netcdf_error_handler(ncstat)
911
912      ncstat = nf_get_var_int(nc_file_id, nc_dstgrdimask_id,  &
913                              grid2_mask_int)
914      call netcdf_error_handler(ncstat)
915
916      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlat_id,  &
917                                             grid2_center_lat)
918      call netcdf_error_handler(ncstat)
919
920      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcntrlon_id,  &
921                                             grid2_center_lon)
922      call netcdf_error_handler(ncstat)
923
924      ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcntrlat_id, 'units', &
925                               grid2_units)
926      call netcdf_error_handler(ncstat)
927
928      select case (grid2_units(1:7))
929      case ('degrees')
930        grid2_center_lat = grid2_center_lat*deg2rad
931        grid2_center_lon = grid2_center_lon*deg2rad
932      case ('radians')
933        !*** no conversion necessary
934      case default
935        print *,'unknown units supplied for grid2 center lat/lon: '
936        print *,'proceeding assuming radians'
937      end select
938
939      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlat_id,  &
940                                             grid2_corner_lat)
941      call netcdf_error_handler(ncstat)
942
943      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdcrnrlon_id,  &
944                                             grid2_corner_lon)
945      call netcdf_error_handler(ncstat)
946
947
948      ncstat = nf_get_att_text(nc_file_id, nc_dstgrdcrnrlat_id, 'units', &
949                               grid2_units)
950      call netcdf_error_handler(ncstat)
951
952      select case (grid2_units(1:7))
953      case ('degrees')
954        grid2_corner_lat = grid2_corner_lat*deg2rad
955        grid2_corner_lon = grid2_corner_lon*deg2rad
956      case ('radians')
957        !*** no conversion necessary
958      case default
959        print *,'unknown units supplied for grid2 corner lat/lon: '
960        print *,'proceeding assuming radians'
961      end select
962
963      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdarea_id,  &
964                                             grid2_area)
965      call netcdf_error_handler(ncstat)
966
967      ncstat = nf_get_var_double(nc_file_id, nc_dstgrdfrac_id,  &
968                                             grid2_frac)
969      call netcdf_error_handler(ncstat)
970
971      ncstat = nf_get_var_int(nc_file_id, nc_srcgrdadd_id,  &
972                              grid1_add_map1)
973      call netcdf_error_handler(ncstat)
974
975      ncstat = nf_get_var_int(nc_file_id, nc_dstgrdadd_id,  &
976                              grid2_add_map1)
977      call netcdf_error_handler(ncstat)
978
979      ncstat = nf_get_var_double(nc_file_id, nc_rmpmatrix_id,  &
980                                             wts1)
981      wts_map1(1,:) = wts1
982      deallocate(wts1)
983
984      if (num_wts > 1) then
985        ncstat = nf_get_var_double(nc_file_id, nc_rmpmatrix2_id,  &
986                                               wts2)
987        wts_map1(2:,:) = wts2
988        deallocate(wts2)
989      endif
990      call netcdf_error_handler(ncstat)
991
992!-----------------------------------------------------------------------
993!
994!     initialize logical mask
995!
996!-----------------------------------------------------------------------
997
998      where (grid1_mask_int == 1)
999        grid1_mask = .true.
1000      elsewhere
1001        grid1_mask = .false.
1002      endwhere
1003      where (grid2_mask_int == 1)
1004        grid2_mask = .true.
1005      elsewhere
1006        grid2_mask = .false.
1007      endwhere
1008      deallocate(grid1_mask_int, grid2_mask_int)
1009
1010!-----------------------------------------------------------------------
1011!
1012!     close input file
1013!
1014!-----------------------------------------------------------------------
1015
1016      ncstat = nf_close(nc_file_id)
1017      call netcdf_error_handler(ncstat)
1018
1019!-----------------------------------------------------------------------
1020
1021      end subroutine read_remap_csm
1022
1023!***********************************************************************
1024
1025      end module remap_read
1026
1027!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Note: See TracBrowser for help on using the repository browser.