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.
modupdatebasic.F in branches/UKMO/dev_r5518_fa_am_dt_deltadelta_toa/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/dev_r5518_fa_am_dt_deltadelta_toa/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modupdatebasic.F @ 7051

Last change on this file since 7051 was 7051, checked in by kuniko, 7 years ago

Removed $Id ...

File size: 16.7 KB
Line 
1!
2! $Id: modupdatebasic.F 2715 2011-03-30 15:58:35Z rblod $
3!
4C     AGRIF (Adaptive Grid Refinement In Fortran)
5C
6C     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7C                        Christophe Vouland (Christophe.Vouland@imag.fr)   
8C
9C     This program is free software; you can redistribute it and/or modify
10C     it under the terms of the GNU General Public License as published by
11C     the Free Software Foundation; either version 2 of the License, or
12C     (at your option) any later version.
13C
14C     This program is distributed in the hope that it will be useful,
15C     but WITHOUT ANY WARRANTY; without even the implied warranty of
16C     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17C     GNU General Public License for more details.
18C
19C     You should have received a copy of the GNU General Public License
20C     along with this program; if not, write to the Free Software
21C     Foundation, Inc., 59 Temple Place -  Suite 330, Boston, MA 02111-1307, USA.
22C
23C
24C
25CCC   Module Agrif_Updatebasic
26C     
27C
28      Module Agrif_Updatebasic
29C
30CCC   Description:
31CCC   Module containing different procedures of update (copy,average,
32CCC   full_weighting) used in the Agrif_Update module.
33C
34C     Modules used:
35C
36      USE Agrif_types
37     
38      IMPLICIT NONE
39     
40      integer,dimension(:,:),allocatable :: indchildcopy
41      integer,dimension(:,:),allocatable :: indchildaverage
42C             
43
44      CONTAINS
45C     Define procedures contained in this module
46C
47C
48C
49C     ************************************************************************** 
50CCC   Subroutine Copy1d 
51C     ************************************************************************** 
52C
53      Subroutine agrif_copy1d(x,y,np,nc,
54     &                  s_parent,s_child,ds_parent,ds_child) 
55C
56CCC   Description:
57CCC   Subroutine to do a copy on a parent grid (vector x) from its child grid 
58CCC   (vector y). 
59C
60CC    Method:
61C
62C     Declarations:
63C
64     
65C       
66C     Arguments
67      INTEGER             :: np,nc     
68      REAL, DIMENSION(np) :: x     
69      REAL, DIMENSION(nc) :: y 
70      REAL                :: s_parent,s_child
71      REAL                :: ds_parent,ds_child
72C
73C     Local variables
74      INTEGER :: i,locind_child_left,coeffraf
75C 
76C
77      coeffraf = nint(ds_parent/ds_child)
78C
79      if (coeffraf == 1) then
80C
81          locind_child_left = 1 + nint((s_parent - s_child)/ds_child)
82C       
83!CDIR ALTCODE
84          x(1:np) = y(locind_child_left:locind_child_left+np-1)
85C
86          return
87C
88      endif
89C
90     
91      locind_child_left = 1 + nint((s_parent - s_child)/ds_child)
92     
93!CDIR ALTCODE
94      do i = 1,np
95C     
96         x(i) = y(locind_child_left)
97C
98         locind_child_left = locind_child_left + coeffraf
99C         
100      enddo   
101       
102C
103      Return
104C
105C
106      End Subroutine agrif_copy1d
107
108C     ************************************************************************** 
109CCC   Subroutine Copy1dprecompute 
110C     ************************************************************************** 
111C
112      Subroutine copy1dprecompute2d(nc2,np,nc,
113     &                  s_parent,s_child,ds_parent,ds_child,dir)
114C
115CCC   Description:
116CCC   Subroutine to precompute index for a copy on a parent grid (vector x) from
117CCC   its child grid (vector y).
118C
119CC    Method:
120C
121C     Declarations:
122C
123
124C       
125C     Arguments
126      INTEGER             :: nc2,np,nc
127      INTEGER             :: dir   
128      REAL                :: s_parent,s_child
129      REAL                :: ds_parent,ds_child
130C
131C     Local variables
132      INTEGER,DIMENSION(:,:),ALLOCATABLE :: indchildcopy_tmp
133      INTEGER :: i,locind_child_left,coeffraf
134C 
135C
136      coeffraf = nint(ds_parent/ds_child)
137C
138      locind_child_left = 1 + nint((s_parent - s_child)/ds_child)
139
140      if (.not.allocated(indchildcopy)) then
141      allocate(indchildcopy(np*nc2,3))
142      else
143      if (size(indchildcopy,1)<np*nc2) then
144      allocate( indchildcopy_tmp(
145     &         size(indchildcopy,1),size(indchildcopy,2)))
146      indchildcopy_tmp = indchildcopy
147      deallocate(indchildcopy)
148      allocate(indchildcopy(np*nc2,3))
149      indchildcopy(1:size(indchildcopy_tmp,1),
150     &              1:size(indchildcopy_tmp,2)) = indchildcopy_tmp
151      deallocate(indchildcopy_tmp)
152      endif
153      endif
154
155
156      do i = 1,np
157C     
158         indchildcopy(i,dir) = locind_child_left
159         locind_child_left = locind_child_left + coeffraf
160C         
161      enddo
162
163      do i =2, nc2
164        indchildcopy(1+(i-1)*np:i*np,dir)=
165     &       indchildcopy(1+(i-2)*np:(i-1)*np,dir) + nc 
166      enddo
167
168C
169      Return
170C
171C
172      End Subroutine copy1dprecompute2d
173
174C
175C
176C     ************************************************************************** 
177CCC   Subroutine Copy1dprecompute 
178C     ************************************************************************** 
179C
180      Subroutine copy1dprecompute(np,nc,
181     &                  s_parent,s_child,ds_parent,ds_child,dir)
182C
183CCC   Description:
184CCC   Subroutine to precompute index for a copy on a parent grid (vector x) from
185CCC   its child grid (vector y).
186C
187CC    Method:
188C
189C     Declarations:
190C
191
192C       
193C     Arguments
194      INTEGER             :: np,nc
195      INTEGER             :: dir
196      REAL                :: s_parent,s_child
197      REAL                :: ds_parent,ds_child
198C
199C     Local variables
200      INTEGER :: i,locind_child_left,coeffraf
201C 
202C
203      coeffraf = nint(ds_parent/ds_child)
204C
205      locind_child_left = 1 + nint((s_parent - s_child)/ds_child)
206
207      if (.not.allocated(indchildcopy)) then
208      allocate(indchildcopy(np,1))
209      else
210      if (size(indchildcopy)<np) then
211      deallocate(indchildcopy)
212      allocate(indchildcopy(np,1))
213      endif
214      endif
215
216!CDIR ALTCODE
217      do i = 1,np
218C     
219         indchildcopy(i,1) = locind_child_left
220         locind_child_left = locind_child_left + coeffraf
221C         
222      enddo
223
224C
225      Return
226C
227C
228      End Subroutine copy1dprecompute
229C
230C
231C     ************************************************************************** 
232CCC   Subroutine Copy1daftercompute 
233C     ************************************************************************** 
234C
235      Subroutine copy1daftercompute(x,y,np,nc,
236     &                  s_parent,s_child,ds_parent,ds_child,dir)
237C
238CCC   Description:
239CCC   Subroutine to do a copy on a parent grid (vector x) from its child grid 
240CCC   (vector y) using precomputed index. 
241C
242CC    Method:
243C
244C     Declarations:
245C
246 
247C       
248C     Arguments
249      INTEGER             :: np,nc
250      INTEGER             :: dir
251      REAL, DIMENSION(np) :: x
252      REAL, DIMENSION(nc) :: y
253      REAL                :: s_parent,s_child
254      REAL                :: ds_parent,ds_child
255C
256C     Local variables
257      INTEGER :: i
258C 
259C
260     
261!CDIR ALTCODE
262      do i = 1,np
263C     
264         x(i) = y(indchildcopy(i,dir))
265C
266      enddo
267 
268C
269      Return
270C
271C
272      End Subroutine copy1daftercompute
273
274
275C
276C     ************************************************************************** 
277CCC   Subroutine Average1d 
278C     ************************************************************************** 
279C   
280      Subroutine average1d(x,y,np,nc,
281     &                     s_parent,s_child,ds_parent,ds_child) 
282C
283CCC   Description:
284CCC   Subroutine to do an update by average on a parent grid (vector x)from its 
285CCC   child grid (vector y).
286C
287C     Arguments
288      INTEGER             :: np,nc     
289      REAL, DIMENSION(np) :: x     
290      REAL, DIMENSION(nc) :: y 
291      REAL                :: s_parent,s_child
292      REAL                :: ds_parent,ds_child
293C
294C     Local variables
295      INTEGER :: i,locind_child_left,coeffraf,ii
296      REAL    :: xpos, invcoeffraf 
297      INTEGER :: nbnonnuls
298      INTEGER :: diffmod
299C 
300C
301      coeffraf = nint(ds_parent/ds_child)
302      invcoeffraf = 1./coeffraf
303C
304      if (coeffraf == 1) then
305C
306          locind_child_left = 1 + nint((s_parent - s_child)/ds_child)
307C       
308          x(1:np) = y(locind_child_left:locind_child_left+np-1)
309C
310          return
311C
312      endif
313C
314      xpos = s_parent     
315     
316      x = 0.
317C
318      diffmod = 0
319     
320      IF ( mod(coeffraf,2) == 0 ) diffmod = 1
321
322        locind_child_left = 1 + agrif_int((xpos - s_child)/ds_child)
323C
324      IF (Agrif_UseSpecialValueInUpdate) THEN
325
326      do i = 1,np
327          nbnonnuls = 0
328!CDIR NOVECTOR
329          Do ii = -coeffraf/2+locind_child_left+diffmod,
330     &                coeffraf/2+locind_child_left
331C 
332               IF (y(ii) .NE. Agrif_SpecialValueFineGrid) THEN
333               nbnonnuls = nbnonnuls + 1
334               x(i) = x(i) + y(ii)
335               ENDIF
336          End Do
337               IF (nbnonnuls .NE. 0) THEN
338                  x(i) = x(i)/nbnonnuls
339               ELSE
340                  x(i) = Agrif_SpecialValueFineGrid
341               ENDIF
342        locind_child_left = locind_child_left + coeffraf
343      enddo
344
345      ELSE 
346
347!CDIR ALTCODE
348      do i = 1,np
349!CDIR NOVECTOR
350          Do ii = -coeffraf/2+locind_child_left+diffmod,
351     &                coeffraf/2+locind_child_left
352C 
353               x(i) = x(i) + y(ii)
354          End Do
355                 x(i) = x(i)*invcoeffraf
356        locind_child_left = locind_child_left + coeffraf
357      enddo
358
359      ENDIF
360
361      Return
362C           
363C     
364      End Subroutine average1d
365     
366      Subroutine average1dprecompute2d(nc2,np,nc,
367     &                     s_parent,s_child,ds_parent,ds_child,dir)
368C
369CCC   Description:
370CCC   Subroutine to do an update by average on a parent grid (vector x)from its 
371CCC   child grid (vector y).
372C
373C     Arguments
374      INTEGER             :: nc2,np,nc,dir
375      REAL                :: s_parent,s_child
376      REAL                :: ds_parent,ds_child
377C
378C     Local variables
379      INTEGER :: i,locind_child_left,coeffraf,ii
380      INTEGER,DIMENSION(:,:),ALLOCATABLE :: indchildaverage_tmp
381      REAL    :: xpos, invcoeffraf
382      INTEGER :: nbnonnuls
383      INTEGER :: diffmod
384C 
385C
386      coeffraf = nint(ds_parent/ds_child)
387C
388      xpos = s_parent
389C
390      diffmod = 0
391
392      IF ( mod(coeffraf,2) == 0 ) diffmod = 1
393
394      locind_child_left = 1 + agrif_int((xpos - s_child)/ds_child)
395
396      if (.not.allocated(indchildaverage)) then
397      allocate(indchildaverage(np*nc2,3))
398      else
399      if (size(indchildaverage,1)<np*nc2) then
400      allocate( indchildaverage_tmp(
401     &         size(indchildaverage,1),size(indchildaverage,2)))
402      indchildaverage_tmp = indchildaverage
403      deallocate(indchildaverage)
404      allocate(indchildaverage(np*nc2,3))
405      indchildaverage(1:size(indchildaverage_tmp,1),
406     &              1:size(indchildaverage_tmp,2)) = indchildaverage_tmp
407      deallocate(indchildaverage_tmp)
408      endif
409      endif
410
411      do i = 1,np
412        indchildaverage(i,dir)= -coeffraf/2+locind_child_left
413     &                                    +diffmod
414        locind_child_left = locind_child_left + coeffraf
415      enddo
416      do i = 2,nc2
417        indchildaverage(1+(i-1)*np:i*np,dir)= 
418     &                 indchildaverage(1+(i-2)*np:(i-1)*np,dir) + nc 
419      enddo
420
421      Return
422C           
423C     
424      End Subroutine average1dprecompute2d
425     
426
427      Subroutine average1dprecompute(np,nc,
428     &                     s_parent,s_child,ds_parent,ds_child) 
429C
430CCC   Description:
431CCC   Subroutine to do an update by average on a parent grid (vector x)from its 
432CCC   child grid (vector y).
433C
434C     Arguments
435      INTEGER             :: np,nc 
436      REAL                :: s_parent,s_child
437      REAL                :: ds_parent,ds_child
438C
439C     Local variables
440      INTEGER :: i,locind_child_left,coeffraf,ii
441      REAL    :: xpos, invcoeffraf 
442      INTEGER :: nbnonnuls
443      INTEGER :: diffmod
444C 
445C
446      coeffraf = nint(ds_parent/ds_child)
447
448C
449      if (coeffraf == 1) then
450C
451          return
452C
453      endif
454           
455C
456      xpos = s_parent
457C
458      diffmod = 0
459     
460      IF ( mod(coeffraf,2) == 0 ) diffmod = 1
461
462        locind_child_left = 1 + agrif_int((xpos - s_child)/ds_child)
463
464      if (.not.allocated(indchildaverage)) then
465      allocate(indchildaverage(np,1))
466      else
467      if (size(indchildaverage,1)<np) then
468      deallocate(indchildaverage)
469      allocate(indchildaverage(np,1))
470      endif
471      endif
472     
473!CDIR ALTCODE
474      do i = 1,np
475C 
476        indchildaverage(i,1)=-coeffraf/2+locind_child_left+diffmod
477        locind_child_left = locind_child_left + coeffraf
478      enddo
479
480      Return
481C           
482C     
483      End Subroutine average1dprecompute
484     
485      Subroutine average1daftercompute(x,y,np,nc,
486     &                     s_parent,s_child,ds_parent,ds_child,dir) 
487C
488CCC   Description:
489CCC   Subroutine to do an update by average on a parent grid (vector x)from its 
490CCC   child grid (vector y).
491C
492C     Arguments
493      INTEGER             :: np,nc,dir     
494      REAL, DIMENSION(np) :: x     
495      REAL, DIMENSION(nc) :: y 
496      REAL                :: s_parent,s_child
497      REAL                :: ds_parent,ds_child
498C
499C     Local variables
500      INTEGER :: i,locind_child_left,coeffraf,ii,j
501      REAL    :: xpos, invcoeffraf 
502      REAL, PARAMETER :: one_third=1./3.
503      INTEGER, DIMENSION(np) :: nbnonnuls
504      INTEGER :: diffmod
505      REAL, DIMENSION(0:5) :: invcoeff=(/1.,1.,0.5,one_third,0.25,0.2/)
506C 
507C
508      coeffraf = nint(ds_parent/ds_child)
509      invcoeffraf = 1./coeffraf     
510C     
511C
512      IF (Agrif_UseSpecialValueInUpdate) THEN
513
514      nbnonnuls = 0
515      do  j =1, coeffraf
516          do i=1, np
517               IF (y(indchildaverage(i,dir) + j -1) .NE. 
518     &               Agrif_SpecialValueFineGrid) THEN
519               nbnonnuls(i) = nbnonnuls(i) + 1
520               x(i) = x(i) +  y(indchildaverage(i,dir) + j -1 )
521               ENDIF
522          End Do
523      enddo
524      do i=1,np
525      x(i)= x(i)*invcoeff(nbnonnuls(i)) 
526      enddo
527
528      ELSE
529!CDIR NOLOOPCHG
530      do  j =1, coeffraf
531!CDIR VECTOR
532          do i=1, np
533             x(i) = x(i) + y(indchildaverage(i,dir) + j -1 )
534          enddo
535      enddo
536       x = x * invcoeffraf
537      ENDIF
538
539      Return
540C           
541C     
542      End Subroutine average1daftercompute           
543C
544C
545C
546C     ************************************************************************** 
547CCC   Subroutine Full_weighting1d 
548C     ************************************************************************** 
549C
550      Subroutine full_weighting1D(x,y,np,nc,
551     &                            s_parent,s_child,ds_parent,ds_child,
552     &                            coeffraf,locind_child_left) 
553C
554CCC   Description:
555CCC   Subroutine to do an update by full_weighting on a parent grid (vector x) 
556CCC   from its child grid (vector y).
557C 
558C     Arguments
559      INTEGER             :: np,nc     
560      REAL, DIMENSION(np) :: x     
561      REAL, DIMENSION(nc) :: y 
562      REAL                :: s_parent,s_child
563      REAL                :: ds_parent,ds_child
564C
565C     Local variables
566      INTEGER :: i,locind_child_left,coeffraf
567      REAL    :: xpos,sumweight,weight
568      INTEGER :: ii,diffmod
569      REAL :: xposfin
570      INTEGER :: it1,it2
571      INTEGER :: i1,i2
572      REAL :: invsumweight
573      REAL :: weights(-(coeffraf):coeffraf)
574     
575C
576C
577      if (coeffraf == 1) then
578C       
579          x(1:np) = y(locind_child_left:locind_child_left+np-1)
580C
581          return
582C
583      endif
584C
585      xpos = s_parent     
586     
587       x = 0.
588
589       xposfin = s_child + ds_child * (locind_child_left - 1)
590       IF (abs(xposfin - xpos).LT.0.001) THEN
591          diffmod = 0
592       ELSE
593          diffmod = 1
594       ENDIF
595C
596             
597      if (diffmod == 1) THEN
598        invsumweight=1./(2.*coeffraf**2)
599        do i=-coeffraf,-1
600          weights(i) = invsumweight*(2*(coeffraf+i)+1)
601        enddo
602        do i=0,coeffraf-1
603          weights(i)=weights(-(i+1))
604        enddo
605        it1 = -coeffraf
606        i1 = -(coeffraf-1)+locind_child_left
607        i2 = 2*coeffraf - 1
608      else
609      invsumweight=1./coeffraf**2
610      do i=-(coeffraf-1),0
611        weights(i) = invsumweight*(coeffraf + i)
612      enddo
613      do i=1,coeffraf-1
614        weights(i) = invsumweight*(coeffraf - i)
615      enddo
616        it1 = -(coeffraf-1)
617        i1 = -(coeffraf-1)+locind_child_left
618        i2 = 2*coeffraf - 2
619      endif
620
621      sumweight = 0                   
622      do i = 1,np
623C
624          it2 = it1
625          Do ii = i1,i1+i2
626C
627           IF (Agrif_UseSpecialValueInUpdate) THEN
628            IF (y(ii) .NE. Agrif_SpecialValueFineGrid) THEN
629               x(i) = x(i) + weights(it2)*y(ii)
630               sumweight = sumweight+weights(it2)
631            ENDIF
632           ELSE           
633               x(i) = x(i) + weights(it2)*y(ii)
634           ENDIF           
635           
636          it2 = it2+1
637          End Do
638
639           IF (Agrif_UseSpecialValueInUpdate) THEN         
640                 IF (sumweight .NE. 0.) THEN
641                    x(i) = x(i)/sumweight
642                    sumweight = 0
643                 ELSE
644                    x(i) = Agrif_SpecialValueFineGrid
645                 ENDIF
646           ENDIF
647       
648        i1 = i1 + coeffraf
649C
650      enddo   
651C
652
653      Return
654C           
655C
656      End Subroutine full_weighting1D 
657
658C
659      End module AGRIF_updatebasic
Note: See TracBrowser for help on using the repository browser.