source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/scrip/src/remap_vars.F @ 4775

Last change on this file since 4775 was 4775, checked in by aclsce, 4 years ago
  • Imported oasis3-mct from Cerfacs svn server (not suppotred anymore).

The version has been extracted from https://oasis3mct.cerfacs.fr/svn/branches/OASIS3-MCT_2.0_branch/oasis3-mct@1818

File size: 15.1 KB
Line 
1C****
2C                    ************************
3C                    *     OASIS MODULE     *
4C                    *     ------------     *
5C                    ************************
6C****
7C***********************************************************************
8C     This module belongs to the SCRIP library. It is modified to run
9C     within OASIS. 
10C     Modifications:
11C       - introduction of logical flags to allow multiple calls of SCRIP
12C       - deallocation of arrays not needed any more
13C       - added CASE for GAUSWGT
14C
15C     Modified by            V. Gayler,  M&D                  20.09.2001
16C     Modified by            D. Declat,  CERFACS              27.06.2002
17C***********************************************************************
18!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
19!
20!     this module contains necessary variables for remapping between
21!     two grids.  also routines for resizing and initializing these
22!     variables.
23!
24!-----------------------------------------------------------------------
25!
26!     CVS:$Id: remap_vars.F 2826 2010-12-10 11:14:21Z valcke $
27!
28!     Copyright (c) 1997, 1998 the Regents of the University of
29!       California.
30!
31!     This software and ancillary information (herein called software)
32!     called SCRIP is made available under the terms described here. 
33!     The software has been approved for release with associated
34!     LA-CC Number 98-45.
35!
36!     Unless otherwise indicated, this software has been authored
37!     by an employee or employees of the University of California,
38!     operator of the Los Alamos National Laboratory under Contract
39!     No. W-7405-ENG-36 with the U.S. Department of Energy.  The U.S.
40!     Government has rights to use, reproduce, and distribute this
41!     software.  The public may copy and use this software without
42!     charge, provided that this Notice and any statement of authorship
43!     are reproduced on all copies.  Neither the Government nor the
44!     University makes any warranty, express or implied, or assumes
45!     any liability or responsibility for the use of this software.
46!
47!     If software is modified to produce derivative works, such modified
48!     software should be clearly marked, so as not to confuse it with
49!     the version available from Los Alamos National Laboratory.
50!
51!***********************************************************************
52
53      module remap_vars
54
55      use kinds_mod
56      use constants
57      use grids
58      USE mod_oasis_flush
59
60      implicit none
61
62!-----------------------------------------------------------------------
63!
64!     module variables
65!
66!-----------------------------------------------------------------------
67
68      integer (kind=int_kind), parameter ::
69     &      norm_opt_none    = 1
70     &,     norm_opt_dstarea = 2
71     &,     norm_opt_frcarea = 3
72     &,     norm_opt_nonorm  = 4
73
74      integer (kind=int_kind), parameter ::
75     &      map_type_conserv  = 1
76     &,     map_type_bilinear = 2
77     &,     map_type_bicubic  = 3
78     &,     map_type_distwgt  = 4
79     &,     map_type_gauswgt  = 5
80
81      integer (kind=int_kind), save :: 
82     &      max_links_map1  ! current size of link arrays
83     &,     num_links_map1  ! actual number of links for remapping
84     &,     max_links_map2  ! current size of link arrays
85     &,     num_links_map2  ! actual number of links for remapping
86     &,     num_maps        ! num of remappings for this grid pair
87     &,     num_wts         ! num of weights used in remapping
88     &,     map_type        ! identifier for remapping method
89     &,     norm_opt        ! option for normalization (conserv only)
90     &,     resize_increment ! default amount to increase array size
91
92      integer (kind=int_kind), dimension(:), allocatable, save ::
93     &      grid1_add_map1, ! grid1 address for each link in mapping 1
94     &      grid2_add_map1, ! grid2 address for each link in mapping 1
95     &      grid1_add_map2, ! grid1 address for each link in mapping 2
96     &      grid2_add_map2  ! grid2 address for each link in mapping 2
97#ifdef TREAT_OVERLAY
98      INTEGER (kind=int_kind), dimension(:), allocatable, save ::
99     &      grid1_add_repl1 ! grid1 address to use after overlap calculation
100#endif TREAT_OVERLAY
101      real (kind=dbl_kind), dimension(:,:), allocatable, save ::
102     &      wts_map1, ! map weights for each link (num_wts,max_links)
103     &      wts_map2  ! map weights for each link (num_wts,max_links)
104
105      logical (kind=log_kind), save :: lfracnnei = .false.
106 
107      logical (kind=log_kind), save :: first_conserv = .true. ! flag to
108                                ! indicate, whether scrip is called from
109                                ! oasis for the first time
110      logical (kind=log_kind), save :: first_call = .true. ! flag used in
111                                ! remap_conserve (store_link_cnsrv)
112
113!***********************************************************************
114
115      contains
116
117!***********************************************************************
118
119      subroutine init_remap_vars (id_scripvoi)
120
121!-----------------------------------------------------------------------
122!
123!     this routine initializes some variables and provides an initial
124!     allocation of arrays (fairly large so frequent resizing
125!     unnecessary).
126!
127!-----------------------------------------------------------------------
128!
129!     input variables
130!
131!-----------------------------------------------------------------------
132
133      INTEGER (kind=int_kind)::
134     &     id_scripvoi          ! number of neighbours for DISTWGT and GAUSWGT
135!-----------------------------------------------------------------------
136!
137      IF (nlogprt .GE. 2) THEN
138         WRITE (UNIT = nulou,FMT = *)' '
139         WRITE (UNIT = nulou,FMT = *)'Entering routine init_remap_vars'
140         CALL OASIS_FLUSH_SCRIP(nulou)
141      ENDIF
142!
143!-----------------------------------------------------------------------
144!
145!     determine the number of weights
146!
147!-----------------------------------------------------------------------
148
149      select case (map_type)
150      case(map_type_conserv)
151        num_wts = 3
152      case(map_type_bilinear)
153        num_wts = 1
154      case(map_type_bicubic)
155          IF (restrict_type == 'REDUCED') THEN
156              num_wts = 1
157          ELSE
158              num_wts = 4
159          ENDIF
160      case(map_type_distwgt)
161        num_wts = 1
162      case(map_type_gauswgt)
163        num_wts = 1
164      end select
165
166!-----------------------------------------------------------------------
167!
168!     initialize num_links and set max_links to four times the largest
169!     of the destination grid sizes initially (can be changed later).
170!     set a default resize increment to increase the size of link
171!     arrays if the number of links exceeds the initial size
172!   
173!-----------------------------------------------------------------------
174
175      num_links_map1 = 0
176      select case (map_type)
177      case(map_type_conserv)
178          max_links_map1 = 4*grid2_size
179      case(map_type_bilinear)
180          max_links_map1 = 4*grid2_size
181      case(map_type_bicubic)
182          IF (restrict_type == 'REDUCED') THEN
183              max_links_map1 = 16*grid2_size
184          ELSE
185              max_links_map1 = 4*grid2_size
186          ENDIF
187      case(map_type_distwgt)
188          max_links_map1 = id_scripvoi*grid2_size
189      case(map_type_gauswgt)
190          max_links_map1 = id_scripvoi*grid2_size
191      END select
192
193      if (num_maps > 1) then
194        num_links_map2 = 0
195        max_links_map1 = max(4*grid1_size,4*grid2_size)
196        max_links_map2 = max_links_map1
197      endif
198
199      resize_increment = 0.1*max(grid1_size,grid2_size)
200
201!-----------------------------------------------------------------------
202!
203!     allocate address and weight arrays for mapping 1
204!   
205!-----------------------------------------------------------------------
206
207      allocate (grid1_add_map1(max_links_map1),
208     &          grid2_add_map1(max_links_map1),
209     &          wts_map1(num_wts, max_links_map1))
210#ifdef TREAT_OVERLAY
211      allocate (grid1_add_repl1(grid1_size))
212#endif TREAT_OVERLAY
213
214!-----------------------------------------------------------------------
215!
216!     allocate address and weight arrays for mapping 2 if necessary
217!   
218!-----------------------------------------------------------------------
219
220      if (num_maps > 1) then
221        allocate (grid1_add_map2(max_links_map2),
222     &            grid2_add_map2(max_links_map2),
223     &            wts_map2(num_wts, max_links_map2))
224      endif
225!-----------------------------------------------------------------------
226!
227!     initialize flag for routine store_link_cnsrv
228!   
229!-----------------------------------------------------------------------
230      first_call = .true.
231
232!-----------------------------------------------------------------------
233!
234      IF (nlogprt .GE. 2) THEN
235         WRITE (UNIT = nulou,FMT = *)' '
236         WRITE (UNIT = nulou,FMT = *)'Leaving routine init_remap_vars'
237         CALL OASIS_FLUSH_SCRIP(nulou)
238      ENDIF
239!
240      end subroutine init_remap_vars
241
242!***********************************************************************
243
244      subroutine resize_remap_vars(nmap, increment)
245
246!-----------------------------------------------------------------------
247!
248!     this routine resizes remapping arrays by increasing(decreasing)
249!     the max_links by increment
250!
251!-----------------------------------------------------------------------
252
253!-----------------------------------------------------------------------
254!
255!     input variables
256!
257!-----------------------------------------------------------------------
258
259      integer (kind=int_kind), intent(in) ::
260     &     nmap,      ! identifies which mapping array to resize
261     &     increment  ! the number of links to add(subtract) to arrays
262
263!-----------------------------------------------------------------------
264!
265!     local variables
266!
267!-----------------------------------------------------------------------
268
269      integer (kind=int_kind) ::
270     &   ierr,     ! error flag
271     &   mxlinks   ! size of link arrays
272
273      integer (kind=int_kind), dimension(:), allocatable ::
274     &   add1_tmp, ! temp array for resizing address arrays
275     &   add2_tmp  ! temp array for resizing address arrays
276
277      real (kind=dbl_kind), dimension(:,:), allocatable ::
278     &   wts_tmp   ! temp array for resizing weight arrays
279!
280      IF (nlogprt .GE. 2) THEN
281         WRITE (UNIT = nulou,FMT = *)' '
282         WRITE (UNIT = nulou,FMT = *)
283     &       'Entering routine resize_remap_vars'
284         CALL OASIS_FLUSH_SCRIP(nulou)
285      ENDIF
286!
287!-----------------------------------------------------------------------
288!
289!     resize map 1 arrays if required.
290!
291!-----------------------------------------------------------------------
292
293      select case (nmap)
294      case(1)
295
296        !***
297        !*** allocate temporaries to hold original values
298        !***
299
300        mxlinks = size(grid1_add_map1)
301        allocate (add1_tmp(mxlinks), add2_tmp(mxlinks), 
302     &            wts_tmp(num_wts,mxlinks))
303
304        add1_tmp = grid1_add_map1
305        add2_tmp = grid2_add_map1
306        wts_tmp  = wts_map1
307       
308        !***
309        !*** deallocate originals and increment max_links then
310        !*** reallocate arrays at new size
311        !***
312
313        deallocate (grid1_add_map1, grid2_add_map1, wts_map1)
314        max_links_map1 = mxlinks + increment
315        allocate (grid1_add_map1(max_links_map1),
316     &            grid2_add_map1(max_links_map1),
317     &            wts_map1(num_wts,max_links_map1))
318
319        !***
320        !*** restore original values from temp arrays and
321        !*** deallocate temps
322        !***
323
324        mxlinks = min(mxlinks, max_links_map1)
325        grid1_add_map1(1:mxlinks) = add1_tmp (1:mxlinks)
326        grid2_add_map1(1:mxlinks) = add2_tmp (1:mxlinks)
327        wts_map1    (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
328        deallocate(add1_tmp, add2_tmp, wts_tmp)
329
330!-----------------------------------------------------------------------
331!
332!     resize map 2 arrays if required.
333!
334!-----------------------------------------------------------------------
335
336      case(2)
337
338        !***
339        !*** allocate temporaries to hold original values
340        !***
341
342        mxlinks = size(grid1_add_map2)
343        allocate (add1_tmp(mxlinks), add2_tmp(mxlinks), 
344     &            wts_tmp(num_wts,mxlinks),stat=ierr)
345        if (ierr .ne. 0) THEN
346            WRITE(nulou,*) '   '
347            WRITE(nulou,*)'error allocating temps in resize: ',ierr
348            CALL OASIS_FLUSH_SCRIP(nulou)
349            stop
350        endif
351
352        add1_tmp = grid1_add_map2
353        add2_tmp = grid2_add_map2
354        wts_tmp  = wts_map2
355       
356        !***
357        !*** deallocate originals and increment max_links then
358        !*** reallocate arrays at new size
359        !***
360
361        deallocate (grid1_add_map2, grid2_add_map2, wts_map2)
362        max_links_map2 = mxlinks + increment
363        allocate (grid1_add_map2(max_links_map2),
364     &            grid2_add_map2(max_links_map2),
365     &            wts_map2(num_wts,max_links_map2),stat=ierr)
366        if (ierr .ne. 0) then
367            WRITE(nulou,*) '   '
368            WRITE(nulou,*)'error allocating new arrays in resize: ',ierr
369            CALL OASIS_FLUSH_SCRIP(nulou)
370            stop
371        endif
372
373
374        !***
375        !*** restore original values from temp arrays and
376        !*** deallocate temps
377        !***
378
379        mxlinks = min(mxlinks, max_links_map2)
380        grid1_add_map2(1:mxlinks) = add1_tmp (1:mxlinks)
381        grid2_add_map2(1:mxlinks) = add2_tmp (1:mxlinks)
382        wts_map2    (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
383        deallocate(add1_tmp, add2_tmp, wts_tmp)
384
385      end select
386!
387      IF (nlogprt .GE. 2) THEN
388         WRITE (UNIT = nulou,FMT = *)' '
389         WRITE (UNIT = nulou,FMT = *)
390     &       'Leaving routine resize_remap_vars'
391         CALL OASIS_FLUSH_SCRIP(nulou)
392      ENDIF
393!
394!-----------------------------------------------------------------------
395
396      end subroutine resize_remap_vars
397!-----------------------------------------------------------------------
398!-----------------------------------------------------------------------
399      subroutine free_remap_vars
400!-----------------------------------------------------------------------
401!
402!     subroutine to deallocate allocated arrays
403!
404!-----------------------------------------------------------------------
405!
406      IF (nlogprt .GE. 2) THEN
407         WRITE (UNIT = nulou,FMT = *)' '
408         WRITE (UNIT = nulou,FMT = *)'Entering routine free_remap_vars'
409         CALL OASIS_FLUSH_SCRIP(nulou)
410      ENDIF
411!
412      deallocate (grid1_add_map1, grid2_add_map1, wts_map1)
413#ifdef TREAT_OVERLAY
414      deallocate (grid1_add_repl1)
415#endif TREAT_OVERLAY
416
417      if (num_maps > 1) then
418        deallocate (grid1_add_map2, grid2_add_map2, wts_map2)
419      endif
420      if (map_type == map_type_conserv) then
421        first_call = .true.
422        first_conserv = .false.
423      endif
424!
425      IF (nlogprt .GE. 2) THEN
426         WRITE (UNIT = nulou,FMT = *)' '
427         WRITE (UNIT = nulou,FMT = *)'Leaving routine free_remap_vars'
428         CALL OASIS_FLUSH_SCRIP(nulou)
429      ENDIF
430!
431!-----------------------------------------------------------------------
432
433      end subroutine free_remap_vars
434
435!***********************************************************************
436
437      end module remap_vars
438
439!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Note: See TracBrowser for help on using the repository browser.