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.
scrip.F90 in utils/tools/WEIGHTS/src – NEMO

source: utils/tools/WEIGHTS/src/scrip.F90 @ 13204

Last change on this file since 13204 was 13204, checked in by smasson, 4 years ago

tools: update with tools_dev_r12970_AGRIF_CMEMS

File size: 7.9 KB
Line 
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!
3!     This routine is the driver for computing the addresses and weights
4!     for interpolating between two grids on a sphere.
5!
6!     Modified slightly to get name of namelist file from command line - sga 2/12/05
7!
8!-----------------------------------------------------------------------
9!
10!     CVS:$Id: scrip.f,v 1.6 2001/08/21 21:06:44 pwjones Exp $
11!
12!     Copyright (c) 1997, 1998 the Regents of the University of
13!       California.
14!
15!     This software and ancillary information (herein called software)
16!     called SCRIP is made available under the terms described here. 
17!     The software has been approved for release with associated
18!     LA-CC Number 98-45.
19!
20!     Unless otherwise indicated, this software has been authored
21!     by an employee or employees of the University of California,
22!     operator of the Los Alamos National Laboratory under Contract
23!     No. W-7405-ENG-36 with the U.S. Department of Energy.  The U.S.
24!     Government has rights to use, reproduce, and distribute this
25!     software.  The public may copy and use this software without
26!     charge, provided that this Notice and any statement of authorship
27!     are reproduced on all copies.  Neither the Government nor the
28!     University makes any warranty, express or implied, or assumes
29!     any liability or responsibility for the use of this software.
30!
31!     If software is modified to produce derivative works, such modified
32!     software should be clearly marked, so as not to confuse it with
33!     the version available from Los Alamos National Laboratory.
34!
35!***********************************************************************
36
37      program scrip
38
39!-----------------------------------------------------------------------
40
41      use kinds_mod                  ! module defining data types
42      use constants                  ! module for common constants
43      use iounits                    ! I/O unit manager
44      use timers                     ! CPU timers
45      use grids                      ! module with grid information
46      use remap_vars                 ! common remapping variables
47      use remap_conservative         ! routines for conservative remap
48      use remap_distance_weight      ! routines for dist-weight remap
49      use remap_bilinear             ! routines for bilinear interp
50      use remap_bicubic              ! routines for bicubic  interp
51      use remap_write                ! routines for remap output
52
53      implicit none
54
55!-----------------------------------------------------------------------
56!
57!     input namelist variables
58!
59!-----------------------------------------------------------------------
60
61      character (char_len) :: &
62         grid1_file,       &  ! filename of grid file containing grid1
63         grid2_file,       &  ! filename of grid file containing grid2
64         interp_file1,     &  ! filename for output remap data (map1)
65         interp_file2,     &  ! filename for output remap data (map2)
66         map1_name,        &  ! name for mapping from grid1 to grid2
67         map2_name,        &  ! name for mapping from grid2 to grid1
68         map_method,       &  ! choice for mapping method
69         normalize_opt,    &  ! option for normalizing weights
70         output_opt           ! option for output conventions
71
72      integer (kind=int_kind) :: &
73         nmap                 ! number of mappings to compute (1 or 2)
74
75      namelist /remap_inputs/ grid1_file, grid2_file,  &
76         interp_file1, interp_file2,                   &
77         map1_name, map2_name, num_maps,               &
78         luse_grid1_area, luse_grid2_area,             &
79         map_method, normalize_opt, output_opt,        &
80         restrict_type, num_srch_bins
81
82!-----------------------------------------------------------------------
83!
84!     local variables
85!
86!-----------------------------------------------------------------------
87
88      integer (kind=int_kind) :: n,   &  ! dummy counter
89                                 iunit   ! unit number for namelist file
90
91      character (char_len) :: nm_in
92
93  if (COMMAND_ARGUMENT_COUNT() == 1) then
94    CALL GET_COMMAND_ARGUMENT(1, nm_in)
95  else
96  write(6,*) 'enter name of namelist file'
97  read(5,*) nm_in
98  endif
99
100!-----------------------------------------------------------------------
101!
102!     initialize timers
103!
104!-----------------------------------------------------------------------
105
106      call timers_init
107      do n=1,max_timers
108        call timer_clear(n)
109      end do
110
111!-----------------------------------------------------------------------
112!
113!     read input namelist
114!
115!-----------------------------------------------------------------------
116
117      grid1_file    = 'unknown'
118      grid2_file    = 'unknown'
119      interp_file1  = 'unknown'
120      interp_file2  = 'unknown'
121      map1_name     = 'unknown'
122      map2_name     = 'unknown'
123      luse_grid1_area = .false.
124      luse_grid2_area = .false.
125      num_maps      = 2
126      map_type      = 1
127      normalize_opt = 'fracarea'
128      output_opt    = 'scrip'
129      restrict_type = 'latitude'
130      num_srch_bins = 900
131
132      call get_unit(iunit)
133      open(iunit, file=nm_in, status='old', form='formatted')
134      read(iunit, nml=remap_inputs)
135      call release_unit(iunit)
136
137      select case(map_method)
138      case ('conservative')
139        map_type = map_type_conserv
140        luse_grid_centers = .false.
141      case ('bilinear')
142        map_type = map_type_bilinear
143        luse_grid_centers = .true.
144      case ('bicubic')
145        map_type = map_type_bicubic
146        luse_grid_centers = .true.
147      case ('distwgt')
148        map_type = map_type_distwgt
149        luse_grid_centers = .true.
150      case default
151        stop 'unknown mapping method'
152      end select
153
154      select case(normalize_opt(1:4))
155      case ('none')
156        norm_opt = norm_opt_none
157      case ('frac')
158        norm_opt = norm_opt_frcarea
159      case ('dest')
160        norm_opt = norm_opt_dstarea
161      case default
162        stop 'unknown normalization option'
163      end select
164
165!-----------------------------------------------------------------------
166!
167!     initialize grid information for both grids
168!
169!-----------------------------------------------------------------------
170
171      call grid_init(grid1_file, grid2_file)
172
173      write(stdout, *) ' Computing remappings between: ',grid1_name
174      write(stdout, *) '                          and  ',grid2_name
175
176!-----------------------------------------------------------------------
177!
178!     initialize some remapping variables.
179!
180!-----------------------------------------------------------------------
181
182      call init_remap_vars
183
184!-----------------------------------------------------------------------
185!
186!     call appropriate interpolation setup routine based on type of
187!     remapping requested.
188!
189!-----------------------------------------------------------------------
190
191      select case(map_type)
192      case(map_type_conserv)
193        call remap_conserv
194      case(map_type_bilinear)
195        call remap_bilin
196      case(map_type_distwgt)
197        call remap_distwgt
198      case(map_type_bicubic)
199        call remap_bicub
200      case default
201        stop 'Invalid Map Type'
202      end select
203
204!-----------------------------------------------------------------------
205!
206!     reduce size of remapping arrays and then write remapping info
207!     to a file.
208!
209!-----------------------------------------------------------------------
210
211      if (num_links_map1 /= max_links_map1) then
212        call resize_remap_vars(1, num_links_map1-max_links_map1)
213      endif
214      if ((num_maps > 1) .and. (num_links_map2 /= max_links_map2)) then
215        call resize_remap_vars(2, num_links_map2-max_links_map2)
216      endif
217
218      call write_remap(map1_name, map2_name, &
219                       interp_file1, interp_file2, output_opt)
220
221!-----------------------------------------------------------------------
222
223      end program scrip
224
225!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Note: See TracBrowser for help on using the repository browser.