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.
modinterpbasic.F in branches/dev_001_GM/AGRIF/AGRIF_FILES – NEMO

source: branches/dev_001_GM/AGRIF/AGRIF_FILES/modinterpbasic.F @ 4310

Last change on this file since 4310 was 662, checked in by opalod, 17 years ago

RB: update Agrif internal routines with a new update scheme and performance improvment

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 30.9 KB
Line 
1!
2! $Id$
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_Interpbasic
26C
27      Module Agrif_Interpbasic
28C
29CCC   Description:
30CCC   Module containing different procedures of interpolation (linear,lagrange,
31CCC   spline,...) used in the Agrif_Interpolation module.
32C
33C     Modules used:
34      USE Agrif_types
35C
36      IMPLICIT NONE
37C             
38      Real,Dimension(Agrif_MaxRaff) :: tabdiff2, tabdiff3
39      Real,Dimension(:),Allocatable::tabtest4     
40       
41      CONTAINS
42C     Define procedures contained in this module
43C 
44C     ************************************************************************** 
45CCC   Subroutine Linear1d 
46C     ************************************************************************** 
47C 
48      Subroutine Linear1d(x,y,np,nc,
49     &                    s_parent,s_child,ds_parent,ds_child) 
50C
51CCC   Description:
52CCC   Subroutine to do a linear 1D interpolation on a child grid (vector y) from
53CCC   its parent grid (vector x). 
54C
55CC    Method:
56C
57C     Declarations:
58C
59     
60C       
61C     Arguments
62      INTEGER             :: np,nc     
63      REAL,INTENT(IN), DIMENSION(np) :: x     
64      REAL,INTENT(OUT), DIMENSION(nc) :: y 
65      REAL                :: s_parent,s_child,ds_parent,ds_child
66C
67C     Local scalars
68      INTEGER :: i,coeffraf,locind_parent_left
69      REAL    :: ypos,globind_parent_left,globind_parent_right
70      REAL    :: invds, invds2
71      REAL :: ypos2,diff
72C
73C
74
75      coeffraf = nint(ds_parent/ds_child)
76C
77      if (coeffraf == 1) then
78C
79          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
80C       
81          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
82C
83          return
84C
85      endif                         
86C
87      ypos = s_child 
88
89      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent)
90
91        globind_parent_left = s_parent 
92     &                        + (locind_parent_left - 1)*ds_parent
93
94        globind_parent_right = globind_parent_left + ds_parent
95
96C
97      invds = 1./ds_parent
98      invds2 = ds_child/ds_parent
99     
100      ypos2 = ypos*invds
101      globind_parent_right=globind_parent_right*invds
102     
103      do i = 1,nc-1
104C
105        if (ypos2 > globind_parent_right) then
106           locind_parent_left = locind_parent_left + 1.
107           globind_parent_right = globind_parent_right + 1.
108        endif
109       
110        diff=(globind_parent_right - ypos2)
111       
112        y(i) = (diff*x(locind_parent_left)
113     &        + (1.-diff)*x(locind_parent_left+1))
114C
115        ypos2 = ypos2 + invds2
116C
117      enddo
118C
119      ypos = s_child + (nc-1)*ds_child
120      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent)
121C
122      if (locind_parent_left == np) then
123C
124          y(nc) = x(np)
125C
126        else
127C
128          globind_parent_left = s_parent 
129     &                        + (locind_parent_left - 1)*ds_parent 
130C     
131          y(nc) = ((globind_parent_left + ds_parent - ypos)
132     &            *x(locind_parent_left)
133     &          + (ypos - globind_parent_left)
134     &            *x(locind_parent_left+1))*invds
135C
136      endif                                         
137C           
138      Return
139C
140C       
141      End Subroutine Linear1d   
142       
143C
144C
145C
146C     ************************************************************************** 
147CCC   Subroutine Lagrange1d 
148C     **************************************************************************
149C
150      Subroutine Lagrange1d(x,y,np,nc,
151     &                      s_parent,s_child,ds_parent,ds_child)
152C
153CCC   Description:
154CCC   Subroutine to do a lagrange 1D interpolation on a child grid (vector y) 
155CCC   from its parent grid (vector x). 
156C
157CC    Method:
158C
159C     Declarations:
160C
161     
162C             
163C     Arguments
164      INTEGER             :: np,nc     
165      REAL,INTENT(IN), DIMENSION(np) :: x     
166      REAL,INTENT(OUT), DIMENSION(nc) :: y 
167      REAL                :: s_parent,s_child,ds_parent,ds_child 
168C
169C     Local scalars
170      INTEGER :: i,coeffraf,locind_parent_left
171      REAL    :: ypos,globind_parent_left
172      REAL    :: X1,X2,X3 
173      real :: deltax,invdsparent
174      real t1,t2,t3,t4,t5,t6,t7,t8
175C
176C 
177      if (np <= 2) then
178C     
179          Call Linear1D(x,y,np,nc,
180     &                  s_parent,s_child,ds_parent,ds_child)
181C         
182         Return
183C 
184      endif
185C
186      coeffraf = nint(ds_parent/ds_child)
187C
188      if (coeffraf == 1) then
189C
190          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
191C       
192          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
193C
194          return
195C
196      endif
197     
198      invdsparent=1./ds_parent
199C
200      ypos = s_child     
201C
202      do i = 1,nc
203C
204        locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent)
205C
206
207        globind_parent_left = s_parent 
208     &                        + (locind_parent_left - 1)*ds_parent 
209     
210        deltax = invdsparent*(ypos-globind_parent_left)
211        ypos = ypos + ds_child
212         if (abs(deltax).LE.0.0001) then
213           y(i)=x(locind_parent_left)
214           
215           cycle
216         endif
217C           
218C
219        t2 = deltax - 2.
220        t3 = deltax - 1.
221        t4 = deltax + 1.
222       
223        t5 = -(1./6.)*deltax*t2*t3
224        t6 = 0.5*t2*t3*t4
225        t7 = -0.5*deltax*t2*t4
226        t8 = (1./6.)*deltax*t3*t4
227       
228        y(i)=t5*x(locind_parent_left-1)+t6*x(locind_parent_left)
229     &  +t7*x(locind_parent_left+1)+t8*x(locind_parent_left+2)
230C
231C
232      enddo
233C
234      return
235C
236C       
237      End Subroutine Lagrange1d 
238C
239C
240C     ************************************************************************** 
241CCC   Subroutine Constant1d 
242C     ************************************************************************** 
243C 
244      Subroutine constant1d(x,y,np,nc,
245     &                    s_parent,s_child,ds_parent,ds_child) 
246C
247CCC   Description:
248CCC   Subroutine to do a linear 1D interpolation on a child grid (vector y) from
249CCC   its parent grid (vector x). 
250C
251CC    Method:
252C
253C     Declarations:
254C
255     
256C       
257C     Arguments
258      INTEGER             :: np,nc     
259      REAL,INTENT(IN), DIMENSION(np) :: x     
260      REAL,INTENT(OUT), DIMENSION(nc) :: y 
261      REAL                :: s_parent,s_child,ds_parent,ds_child
262C
263C     Local scalars
264      INTEGER :: i,coeffraf,locind_parent
265      REAL    :: ypos
266C
267C
268
269      coeffraf = nint(ds_parent/ds_child)
270C
271      if (coeffraf == 1) then
272C
273          locind_parent = 1 + nint((s_child - s_parent)/ds_parent)
274C       
275          y(1:nc) = x(locind_parent:locind_parent+nc-1)
276C
277          return
278C
279      endif                         
280C
281      ypos = s_child     
282C
283      do i = 1,nc
284C
285        locind_parent = 1 + nint((ypos - s_parent)/ds_parent)
286C
287        y(i) = x(locind_parent)
288C
289        ypos = ypos + ds_child
290C
291      enddo
292C           
293      Return
294C
295C       
296      End Subroutine constant1d   
297C
298C     ************************************************************************** 
299CCC   Subroutine Linear1dconserv
300C     ************************************************************************** 
301C 
302      Subroutine Linear1dconserv(x,y,np,nc,
303     &                    s_parent,s_child,ds_parent,ds_child) 
304C
305CCC   Description:
306CCC   Subroutine to do a linear 1D interpolation on a child grid (vector y) from
307CCC   its parent grid (vector x). 
308C
309CC    Method:
310C
311C     Declarations:
312C
313      Implicit none
314C       
315C     Arguments
316      Integer             :: np,nc     
317      Real, Dimension(np) :: x     
318      Real, Dimension(nc) :: y
319      Real, Dimension(:),Allocatable :: ytemp
320      Real                :: s_parent,s_child,ds_parent,ds_child
321C
322C     Local scalars
323      Integer :: i,coeffraf,locind_parent_left,locind_parent_last
324      Real    :: ypos
325      integer :: i1,i2,ii
326      real :: xpmin,xpmax,slope
327      INTEGER :: diffmod
328      REAL :: xdiffmod
329
330C
331C
332
333      coeffraf = nint(ds_parent/ds_child)
334C
335      If (coeffraf == 1) Then
336C
337          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
338C       
339          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
340C
341          return
342C
343      End If
344C             
345      diffmod = 0
346      IF (mod(coeffraf,2) == 0) diffmod = 1 
347
348      xdiffmod = real(diffmod)/2.
349                         
350      allocate(ytemp(-2*coeffraf:nc+2*coeffraf))
351C
352      ypos = s_child 
353C     
354      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) 
355
356      locind_parent_last = 1 +
357     &  agrif_ceiling((ypos +(nc - 1) *ds_child - s_parent)/ds_parent)
358     
359      xpmin = s_parent + (locind_parent_left-1)*ds_parent 
360      xpmax = s_parent + (locind_parent_last-1)*ds_parent         
361     
362      i1 = 1+agrif_int((xpmin-s_child)/ds_child)
363      i2 = 1+agrif_int((xpmax-s_child)/ds_child)
364
365      i = i1
366     
367      if (locind_parent_left == 1) then
368        slope=
369     &   (x(locind_parent_left+1)-x(locind_parent_left))/(coeffraf)
370      else
371         slope=
372     &   (x(locind_parent_left+1)-x(locind_parent_left-1))/(2.*coeffraf)
373      endif
374     
375        do ii=i-coeffraf/2+diffmod,i+coeffraf/2
376          ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope
377        enddo 
378
379        locind_parent_left = locind_parent_left + 1
380                     
381      do i=i1 +  coeffraf, i2 - coeffraf,coeffraf
382        slope=
383     &   (x(locind_parent_left+1)-x(locind_parent_left-1))/(2.*coeffraf)
384        do ii=i-coeffraf/2+diffmod,i+coeffraf/2
385          ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope
386        enddo
387        locind_parent_left = locind_parent_left + 1
388      enddo
389     
390      i = i2
391     
392      if (locind_parent_left == np) then
393        slope=
394     &   (x(locind_parent_left)-x(locind_parent_left-1))/(coeffraf)
395      else
396         slope=
397     &   (x(locind_parent_left+1)-x(locind_parent_left-1))/(2.*coeffraf)
398      endif
399     
400        do ii=i-coeffraf/2+diffmod,nc
401          ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope
402        enddo     
403C
404      y(1:nc)=ytemp(1:nc)                                   
405C           
406      deallocate(ytemp)
407      Return
408C       
409      End Subroutine Linear1dconserv
410     
411C
412C     ************************************************************************** 
413CCC   Subroutine Linear1dconservlim
414C     ************************************************************************** 
415C 
416      Subroutine Linear1dconservlim(x,y,np,nc,
417     &                    s_parent,s_child,ds_parent,ds_child) 
418C
419CCC   Description:
420CCC   Subroutine to do a linear 1D interpolation on a child grid (vector y) from
421CCC   its parent grid (vector x). 
422C
423CC    Method:
424C
425C     Declarations:
426C
427      Implicit none
428C       
429C     Arguments
430      Integer             :: np,nc     
431      Real, Dimension(np) :: x     
432      Real, Dimension(nc) :: y
433      Real, Dimension(:),Allocatable :: ytemp
434      Real                :: s_parent,s_child,ds_parent,ds_child
435C
436C     Local scalars
437      Integer :: i,coeffraf,locind_parent_left,locind_parent_last
438      Real    :: ypos
439      integer :: i1,i2,ii
440      real :: xpmin,xpmax,slope
441      INTEGER :: diffmod
442      real :: xdiffmod
443C
444C
445
446      coeffraf = nint(ds_parent/ds_child)
447C
448      If (coeffraf == 1) Then
449C
450          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
451C       
452          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
453C
454          return
455C
456      End If
457C
458      IF (coeffraf .NE.3) THEN
459      print *,'LINEARCONSERVLIM not ready for refinement ratio = ',
460     &   coeffraf
461      stop
462      ENDIF   
463     
464      diffmod = 0
465      IF (mod(coeffraf,2) == 0) diffmod = 1       
466
467      xdiffmod = real(diffmod)/2.
468                         
469      allocate(ytemp(-2*coeffraf:nc+2*coeffraf))
470C
471      ypos = s_child 
472C     
473      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) 
474
475      locind_parent_last = 1 +
476     &  agrif_ceiling((ypos +(nc - 1) *ds_child - s_parent)/ds_parent)
477     
478      xpmin = s_parent + (locind_parent_left-1)*ds_parent 
479      xpmax = s_parent + (locind_parent_last-1)*ds_parent         
480     
481      i1 = 1+agrif_int((xpmin-s_child)/ds_child)
482      i2 = 1+agrif_int((xpmax-s_child)/ds_child)
483
484      i = i1
485     
486      if (locind_parent_left == 1) then
487        slope=0.       
488      else
489        slope = vanleer(x(locind_parent_left-1:locind_parent_left+1))
490        slope = slope / coeffraf   
491      endif
492     
493        do ii=i-coeffraf/2+diffmod,i+coeffraf/2
494          ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope
495        enddo 
496
497        locind_parent_left = locind_parent_left + 1
498                     
499      do i=i1 +  coeffraf, i2 - coeffraf,coeffraf     
500        slope = vanleer(x(locind_parent_left-1:locind_parent_left+1))
501        slope = slope / coeffraf
502
503        do ii=i-coeffraf/2+diffmod,i+coeffraf/2
504          ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope
505        enddo
506        locind_parent_left = locind_parent_left + 1
507      enddo
508     
509      i = i2
510     
511      if (locind_parent_left == np) then
512        slope=0.     
513      else
514        slope = vanleer(x(locind_parent_left-1:locind_parent_left+1))
515        slope = slope / coeffraf     
516      endif
517     
518        do ii=i-coeffraf/2+diffmod,nc
519          ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope
520        enddo     
521C
522      y(1:nc)=ytemp(1:nc)                                   
523C           
524      deallocate(ytemp)
525      Return
526C       
527      End Subroutine Linear1dconservlim     
528C         
529
530C     ************************************************************************** 
531CCC   Subroutine ppm1d
532C     ************************************************************************** 
533C 
534      Subroutine ppm1d(x,y,np,nc,
535     &                    s_parent,s_child,ds_parent,ds_child) 
536C
537CCC   Description:
538CCC   Subroutine to do a 1D interpolation and apply monotonicity constraints
539CCC   using piecewise parabolic method 
540CCC   on a child grid (vector y) from its parent grid (vector x).
541CC    Method:
542C
543C     Declarations:
544C
545      Implicit none
546C       
547C     Arguments
548      Integer             :: np,nc     
549      Real, INTENT(IN),Dimension(np) :: x     
550      Real, INTENT(OUT),Dimension(nc) :: y
551C      Real, Dimension(:),Allocatable :: ytemp
552      Real                :: s_parent,s_child,ds_parent,ds_child
553C
554C     Local scalars
555      Integer :: i,coeffraf,locind_parent_left,locind_parent_last
556      Integer :: iparent,ipos,pos,nmin,nmax
557      Real    :: ypos
558      integer :: i1,jj
559      Real :: xpmin,a
560C     
561      Real :: xrmin,xrmax,am3,s2,s1 
562      Real, Dimension(np) :: xl,delta,a6,slope
563C      Real, Dimension(:),Allocatable  :: diff,diff2,diff3   
564      INTEGER :: diffmod
565      REAL :: invcoeffraf
566C     
567      coeffraf = nint(ds_parent/ds_child)
568C
569      If (coeffraf == 1) Then
570          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
571          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
572          return
573      End If
574      invcoeffraf = ds_child/ds_parent
575C     
576
577      IF( .NOT. allocated(tabtest4) ) THEN   
578      Allocate(tabtest4(-2*coeffraf:nc+2*coeffraf))
579      ELSE
580         IF (size(tabtest4) .LT. nc+4*coeffraf+1)THEN
581         deallocate( tabtest4 )
582         Allocate(tabtest4(-2*coeffraf:nc+2*coeffraf))
583         ENDIF
584      ENDIF 
585      ypos = s_child 
586C
587      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) 
588      locind_parent_last = 1 +
589     &      agrif_ceiling((ypos +(nc - 1) 
590     &      *ds_child - s_parent)/ds_parent) 
591C
592      xpmin = s_parent + (locind_parent_left-1)*ds_parent       
593      i1 = 1+agrif_int((xpmin-s_child)/ds_child)       
594C     
595C
596     
597      Do i=1,coeffraf
598        tabdiff2(i)=(real(i)-0.5)*invcoeffraf
599      EndDo
600
601      a = invcoeffraf**2 
602      tabdiff3(1) = (1./3.)*a
603      a=2.*a
604      Do i=2,coeffraf
605        tabdiff3(i) = tabdiff3(i-1)+(real(i)-1)*a
606      EndDo
607C
608      if( locind_parent_last+2 <= np ) then
609           nmax = locind_parent_last+2   
610      else if( locind_parent_last+1 <= np ) then
611           nmax = locind_parent_last+1
612      else
613           nmax = locind_parent_last 
614      endif     
615C     
616      if(locind_parent_left-1 >= 1) then
617          nmin = locind_parent_left-1
618      else
619          nmin = locind_parent_left
620      endif   
621C 
622C
623      Do i = nmin,nmax
624         slope(i) = x(i) - x(i-1)
625      Enddo
626
627      Do i = nmin+1,nmax-1
628         xl(i)= 0.5*(x(i-1)+x(i))
629     &      -0.08333333333333*(slope(i+1)-slope(i-1)) 
630      Enddo
631C
632C apply parabolic monotonicity
633       Do i = locind_parent_left,locind_parent_last
634          delta(i) = xl(i+1) - xl(i)
635          a6(i) = 6.*x(i)-3.*(xl(i) +xl(i+1))
636C
637       End do   
638C
639        diffmod = 0
640       IF (mod(coeffraf,2) == 0) diffmod = 1
641C
642        ipos = i1
643C               
644        Do iparent = locind_parent_left,locind_parent_last       
645             pos=1
646             Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2
647C
648               tabtest4(jj) = xl(iparent)   
649     &             + tabdiff2(pos) *  (delta(iparent)+a6(iparent))
650     &             - tabdiff3(pos) *  a6(iparent)
651               pos = pos+1 
652             End do
653             ipos = ipos + coeffraf
654C
655        End do     
656C
657C
658        y(1:nc)=tabtest4(1:nc)
659       
660      Return
661      End Subroutine ppm1d
662     
663C     ************************************************************************** 
664CCC   Subroutine weno1d
665C     ************************************************************************** 
666C 
667      Subroutine weno1dnew(x,y,np,nc,
668     &                    s_parent,s_child,ds_parent,ds_child) 
669C
670CCC   Description:
671CCC   Subroutine to do a 1D interpolation and apply monotonicity constraints
672CCC   using piecewise parabolic method 
673CCC   on a child grid (vector y) from its parent grid (vector x).
674CC    Method:
675C
676C     Declarations:
677C
678      Implicit none
679C       
680C     Arguments
681      Integer             :: np,nc
682      Real, Dimension(np) :: x
683      Real, Dimension(nc) :: y
684      Real, Dimension(:),Allocatable :: ytemp
685      Real                :: s_parent,s_child,ds_parent,ds_child
686C
687C     Local scalars
688      Integer :: i,coeffraf,locind_parent_left,locind_parent_last
689      Integer :: iparent,ipos,pos,nmin,nmax
690      Real    :: ypos
691      integer :: i1,jj
692      Real :: xpmin,cavg,a,b
693C     
694      Real :: xrmin,xrmax,am3,s2,s1
695      Real, Dimension(np) :: xr,xl,delta,a6,slope,slope2,smooth
696      Real, Dimension(:),Allocatable  :: diff,diff2,diff3
697      INTEGER :: diffmod
698      REAL :: invcoeffraf
699      integer :: s,l,k
700      integer :: etan, etap
701      real :: delta0, delta1, delta2
702      real :: epsilon
703      parameter (epsilon = 1.D-8)
704      real, dimension(:,:), allocatable :: ak, ck
705C     
706      coeffraf = nint(ds_parent/ds_child)
707C
708      If (coeffraf == 1) Then
709          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
710          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
711          return
712      End If
713      invcoeffraf = ds_child/ds_parent
714      Allocate(ak(0:1,coeffraf))
715      Allocate(ck(0:1,coeffraf))
716           
717C     
718      Allocate(ytemp(-2*coeffraf:nc+2*coeffraf))
719      ypos = s_child
720C
721      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent)
722      locind_parent_last = 1 +
723     &      agrif_ceiling((ypos +(nc - 1)
724     &      *ds_child - s_parent)/ds_parent)
725C
726      xpmin = s_parent + (locind_parent_left-1)*ds_parent
727      i1 = 1+agrif_int((xpmin-s_child)/ds_child)
728C     
729      Allocate( diff(coeffraf),diff2(coeffraf),diff3(coeffraf) )
730C     
731      diff(1)=0.5*invcoeffraf
732      do i=2,coeffraf
733         diff(i) = diff(i-1)+invcoeffraf
734      enddo
735     
736      ak = 0.
737      ck = 0.
738     
739      do i=1,coeffraf
740         do k=0,1
741         do s=0,2
742         do l=0,2
743           if (l /= s) then
744            ak(k,i) = ak(k,i)+(diff(i)-(k-l+1.))
745           endif
746         enddo
747         enddo
748         enddo
749               
750         etap = 0
751         etan = 0
752         do k=0,1
753          if (ak(k,i) > 0) then
754            etap = etap+1
755          else if (ak(k,i) < 0) then
756            etan = etan + 1
757          endif
758         enddo
759               
760         do k=0,1
761           if (ak(k,i) == 0) then
762            Ck(k,i) = 1.
763           else if (ak(k,i) > 0) then
764            Ck(k,i) = 1./(etap * ak(k,i))
765           else
766            Ck(k,i) = -1./(etan * ak(k,i))
767           endif
768         enddo
769      enddo
770                     
771C     
772      a = 0.
773      b = invcoeffraf
774      Do i=1,coeffraf
775         diff2(i) = 0.5*(b*b - a*a) 
776         diff3(i) = (1./3.)*(b*b*b - a*a*a)
777         a = a + invcoeffraf
778         b = b + invcoeffraf
779      End do
780C
781      if( locind_parent_last+2 <= np ) then
782           nmax = locind_parent_last+2   
783      elseif( locind_parent_last+1 <= np ) then
784           nmax = locind_parent_last+1
785      else
786           nmax = locind_parent_last 
787      endif     
788C     
789      if(locind_parent_left-2 >= 1) then
790          nmin = locind_parent_left-2
791      elseif(locind_parent_left-1 >= 1) then
792          nmin = locind_parent_left-1
793      else
794          nmin = locind_parent_left
795      endif   
796C 
797      Do i = nmin+1,nmax
798         slope(i) = (x(i) - x(i-1))
799      Enddo 
800      DO i=nmin+2,nmax
801        smooth(i) = 0.5*(slope(i)**2+slope(i-1)**2)
802     &       +(slope(i)-slope(i-1))**2
803      enddo
804C
805        diffmod = 0
806        IF (mod(coeffraf,2) == 0) diffmod = 1           
807C
808        ipos = i1
809C               
810        Do iparent = locind_parent_left,locind_parent_last       
811             pos=1
812             
813            delta0=1./(epsilon+smooth(iparent  ))**3
814            delta1=1./(epsilon+smooth(iparent+1))**3
815            delta2=1./(epsilon+smooth(iparent+2))**3   
816                     
817             Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2
818C
819               pos = pos+1 
820             End do
821             ipos = ipos + coeffraf
822C
823        End do     
824C
825C
826        y(1:nc)=ytemp(1:nc)                                 
827        deallocate(ytemp)               
828        deallocate(diff, diff2, diff3)
829       
830        deallocate(ak,ck)
831       
832      Return
833      End Subroutine weno1dnew
834     
835C     ************************************************************************** 
836CCC   Subroutine weno1d
837C     ************************************************************************** 
838C 
839      Subroutine weno1d(x,y,np,nc,
840     &                    s_parent,s_child,ds_parent,ds_child) 
841C
842CCC   Description:
843CCC   Subroutine to do a 1D interpolation and apply monotonicity constraints
844CCC   using piecewise parabolic method 
845CCC   on a child grid (vector y) from its parent grid (vector x).
846CC    Method:
847C
848C     Declarations:
849C
850      Implicit none
851C       
852C     Arguments
853      Integer             :: np,nc
854      Real, Dimension(np) :: x
855      Real, Dimension(nc) :: y
856      Real, Dimension(:),Allocatable :: ytemp
857      Real                :: s_parent,s_child,ds_parent,ds_child
858C
859C     Local scalars
860      Integer :: i,coeffraf,locind_parent_left,locind_parent_last
861      Integer :: iparent,ipos,pos,nmin,nmax
862      Real    :: ypos
863      integer :: i1,jj
864      Real :: xpmin,cavg,a,b
865C     
866      Real :: xrmin,xrmax,am3,s2,s1
867      Real, Dimension(np) :: xr,xl,delta,a6,slope,slope2
868      Real, Dimension(:),Allocatable  :: diff,diff2,diff3
869      INTEGER :: diffmod
870      REAL :: invcoeffraf
871      integer :: s,l,k
872      integer :: etan, etap
873      real :: delta0, delta1,sumdelta
874      real :: epsilon
875      parameter (epsilon = 1.D-8)
876C     
877      coeffraf = nint(ds_parent/ds_child)
878C
879      If (coeffraf == 1) Then
880          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
881          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
882          return
883      End If
884      invcoeffraf = ds_child/ds_parent
885           
886C     
887      Allocate(ytemp(-2*coeffraf:nc+2*coeffraf))
888      ypos = s_child
889C
890      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent)
891      locind_parent_last = 1 +
892     &      agrif_ceiling((ypos +(nc - 1)
893     &      *ds_child - s_parent)/ds_parent)
894C
895      xpmin = s_parent + (locind_parent_left-1)*ds_parent
896      i1 = 1+agrif_int((xpmin-s_child)/ds_child)
897C     
898      Allocate( diff(coeffraf))
899C     
900      diff(1)=0.5*invcoeffraf
901      do i=2,coeffraf
902         diff(i) = diff(i-1)+invcoeffraf
903      enddo
904C
905      if( locind_parent_last+2 <= np ) then
906           nmax = locind_parent_last+2   
907      else if( locind_parent_last+1 <= np ) then
908           nmax = locind_parent_last+1
909      else
910           nmax = locind_parent_last 
911      endif     
912C     
913      if(locind_parent_left-1 >= 1) then
914          nmin = locind_parent_left-1
915      else
916          nmin = locind_parent_left
917      endif   
918C 
919      Do i = nmin+1,nmax
920         slope(i) = (x(i) - x(i-1))
921      Enddo 
922C
923        diffmod = 0
924        IF (mod(coeffraf,2) == 0) diffmod = 1           
925C
926        ipos = i1
927C               
928        Do iparent = locind_parent_left,locind_parent_last       
929             pos=1
930            delta0=1./(epsilon+slope(iparent  )**2)**2
931            delta1=1./(epsilon+slope(iparent+1)**2)**2
932            sumdelta = 1./(delta0+delta1)             
933             Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2
934C                             
935          ytemp(jj) = x(iparent)+(diff(pos)-0.5)*(
936     &                delta0*slope(iparent)+
937     &                delta1*slope(iparent+1))*sumdelta
938               pos = pos+1 
939             End do
940             ipos = ipos + coeffraf
941C
942        End do     
943C
944C
945        y(1:nc)=ytemp(1:nc)                                 
946        deallocate(ytemp)               
947        deallocate(diff)
948       
949      Return
950      End Subroutine weno1d
951                 
952C                             
953C     ************************************************************************** 
954CCC   Subroutine eno1d
955C     ************************************************************************** 
956C 
957      Subroutine eno1d(x,y,np,nc,
958     &                    s_parent,s_child,ds_parent,ds_child) 
959C
960CCC   Description:
961CCC   ---- p 163-164 Computational gasdynamics ----
962CCC   Subroutine to do a 1D interpolation 
963CCC   using piecewise polynomial ENO reconstruction technique 
964CCC   on a child grid (vector y) from its parent grid (vector x).
965CC    Method:
966C
967C     Declarations:
968C
969      Implicit none
970C       
971C     Arguments
972      Integer             :: np,nc     
973      Real, Dimension(np) :: x     
974      Real, Dimension(nc) :: y
975      Real, Dimension(:),Allocatable :: ytemp
976      Real                :: s_parent,s_child,ds_parent,ds_child
977C
978C     Local scalars
979      Integer :: i,coeffraf,locind_parent_left,locind_parent_last
980      Integer :: ipos,pos
981      Real    :: ypos,xi
982      integer :: i1,jj
983      Real :: xpmin,cavg
984C 
985      Real, Dimension(3,np) :: dd,c
986      Integer :: left
987C           
988      Real, DImension(1:np+1) :: xhalf
989      Real, Dimension(:,:),Allocatable  :: Xbar 
990      INTEGER :: diffmod     
991C
992      coeffraf = nint(ds_parent/ds_child)
993C     
994      If (coeffraf == 1) Then
995          locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent)
996          y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1)
997          return
998      End If     
999     
1000      diffmod = 0
1001      IF (mod(coeffraf,2) == 0) diffmod = 1 
1002C     
1003      Allocate(ytemp(-2*coeffraf:nc+2*coeffraf))
1004      ypos = s_child 
1005      locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) 
1006      locind_parent_last = 1 +
1007     &      agrif_ceiling((ypos +(nc - 1) *ds_child - 
1008     &      s_parent)/ds_parent)       
1009      xpmin = s_parent + (locind_parent_left-1)*ds_parent     
1010      i1 = 1+agrif_int((xpmin-s_child)/ds_child)           
1011C     
1012      xhalf(np+1) = np + 0.5
1013      Do i = 1,np
1014          xhalf(i) = i - 0.5
1015      Enddo
1016C
1017C compute divided differences
1018C
1019      dd(1,1:np) = x(1:np)
1020      dd(2,1:np-1) = 0.5*( dd(1,2:np) - dd(1,1:np-1) )
1021      dd(3,1:np-2) = (1./3.)*( dd(2,2:np-1) - dd(2,1:np-2) )
1022C
1023      Allocate( Xbar( coeffraf,2 ) )
1024      xi = 0.5
1025      Do i = 1,coeffraf
1026        Xbar(i,1) = (i-1)*ds_child/ds_parent - xi
1027        Xbar(i,2) = i*ds_child/ds_parent - xi
1028      Enddo
1029C
1030      ipos = i1
1031C           
1032      DO i = locind_parent_left,locind_parent_last           
1033         left = i           
1034         do jj = 2,3
1035             If(abs(dd(jj,left)) .gt. abs(dd(jj,left-1)))
1036     &               left = left-1           
1037         enddo
1038C           
1039C  convert to Taylor series form
1040C
1041         Call Taylor(i,xhalf(left:left+2),dd(1:3,left),c(1:3,i))
1042      ENDDO     
1043C
1044C evaluate the reconstruction on each cell
1045C
1046       DO i = locind_parent_left,locind_parent_last 
1047C
1048         cavg = 0.
1049         pos = 1.           
1050C         
1051         Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2
1052           ytemp(jj) =(c(1,i)*(Xbar(pos,2)-Xbar(pos,1))
1053     &                +c(2,i)*(Xbar(pos,2)*Xbar(pos,2)-
1054     &                         Xbar(pos,1)*Xbar(pos,1))
1055     &                +c(3,i)*(Xbar(pos,2)*Xbar(pos,2)*Xbar(pos,2)-
1056     &                         Xbar(pos,1)*Xbar(pos,1)*Xbar(pos,1)))
1057     &                         *coeffraf
1058           cavg = cavg + ytemp(jj)
1059           pos = pos+1
1060         Enddo               
1061         ipos = ipos + coeffraf                 
1062      ENDDO
1063C
1064      y(1:nc)=ytemp(1:nc)                                           
1065      deallocate(ytemp,Xbar)                 
1066C     
1067      Return       
1068      End Subroutine eno1d
1069C 
1070C     
1071C     ************************************************************************** 
1072CCC   Subroutine taylor
1073C     ************************************************************************** 
1074C       
1075      subroutine taylor(ind,xhalf,dd,c)     
1076C     
1077      Integer :: ind
1078      real,dimension(3) :: dd,c     
1079      real,dimension(0:3,0:3) :: d 
1080      real,dimension(3) :: xhalf   
1081      integer ::i,j
1082C     
1083C
1084      d(0,0:3)=1.
1085      do i = 1,3
1086         d(i,0)=(ind-xhalf(i))*d(i-1,0)
1087      enddo 
1088C     
1089      do i = 1,3
1090         do j = 1,3-i
1091           d(i,j) = d(i,j-1) + (ind-xhalf(i+j))*d(i-1,j)
1092         enddo
1093      enddo
1094C         
1095      do j = 1,3       
1096         c(j) = 0.
1097         do i=0,3-j
1098            c(j) = c(j) + d(i,j)*dd(i+j)         
1099         enddo
1100      enddo
1101C     
1102      end subroutine taylor   
1103     
1104     
1105      REAL FUNCTION vanleer(tab)
1106      REAL, DIMENSION(3) :: tab
1107       real res1
1108       real p1,p2,p3
1109       
1110       p1=(tab(3)-tab(1))/2.
1111       p2=2.*(tab(2)-tab(1))
1112       p3=2.*(tab(3)-tab(2))
1113       
1114       if ((p1>0.).AND.(p2>0.).AND.(p3>0)) then
1115          res1=minval((/p1,p2,p3/))
1116       elseif ((p1<0.).AND.(p2<0.).AND.(p3<0)) then
1117          res1=maxval((/p1,p2,p3/))
1118       else
1119          res1=0.
1120       endif
1121         
1122          vanleer = res1   
1123     
1124     
1125      END FUNCTION vanleer 
1126
1127C     
1128      End Module Agrif_Interpbasic
Note: See TracBrowser for help on using the repository browser.