/[lmdze]/trunk/Sources/dyn3d/invert_zoom_x.f
ViewVC logotype

Diff of /trunk/Sources/dyn3d/invert_zoom_x.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

trunk/dyn3d/fxhyp_loop_ik.f revision 124 by guez, Thu Feb 5 15:19:37 2015 UTC trunk/Sources/dyn3d/invert_zoom_x.f revision 154 by guez, Tue Jul 7 17:49:23 2015 UTC
# Line 1  Line 1 
1  module fxhyp_loop_ik_m  module invert_zoom_x_m
2    
3    implicit none    implicit none
4    
5    INTEGER, PARAMETER:: nmax = 30000    INTEGER, PARAMETER:: nmax = 30000
6      DOUBLE PRECISION abs_y
7    
8      private abs_y, funcd
9    
10  contains  contains
11    
12    subroutine fxhyp_loop_ik(ik, decalx, xf, xtild, Xprimt, xzoom, xlon, &    subroutine invert_zoom_x(xf, xtild, G, xlon, xprim, xuv)
        xprimm, xuv)  
13    
14      use coefpoly_m, only: coefpoly      use coefpoly_m, only: coefpoly, a1, a2, a3
15      USE dimens_m, ONLY: iim      USE dimens_m, ONLY: iim
16        use dynetat0_m, only: clon
17      use nr_util, only: pi_d, twopi_d      use nr_util, only: pi_d, twopi_d
18      use serre, only: grossismx      use numer_rec_95, only: hunt, rtsafe
19    
20        DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), G(0:) ! (0:nmax)
21    
22      INTEGER, intent(in):: ik      real, intent(out):: xlon(:), xprim(:) ! (iim)
     DOUBLE PRECISION, intent(in):: decalx  
     DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), Xprimt(0:) ! (0:2 * nmax)  
     DOUBLE PRECISION, intent(in):: xzoom  
     real, intent(out):: xlon(:), xprimm(:) ! (iim)  
23    
24      DOUBLE PRECISION, intent(in):: xuv      DOUBLE PRECISION, intent(in):: xuv
25        ! between - 0.25 and 0.5
26      ! 0. si calcul aux points scalaires      ! 0. si calcul aux points scalaires
27      ! 0.5 si calcul aux points U      ! 0.5 si calcul aux points U
28    
29      ! Local:      ! Local:
30      DOUBLE PRECISION xo1, Xfi, xi, a0, a1, a2, a3, Xf1, Xprimin      DOUBLE PRECISION Y
31      integer ii1, ii2, i, it, iter      DOUBLE PRECISION h ! step of the uniform grid
32      DOUBLE PRECISION, parameter:: my_eps = 1e-6      integer i, it
     DOUBLE PRECISION xxprim(iim + 1), xvrai(iim + 1)  
   
     !------------------------------------------------------------------  
   
     IF (ik == 1 .and. grossismx == 1.) THEN  
        ii1 = 2  
        ii2 = iim + 1  
     else  
        ii1=1  
        ii2=iim  
     END IF  
   
     DO i = ii1, ii2  
        Xfi = - pi_d + (i + xuv - decalx) * twopi_d / iim  
   
        it = 2 * nmax  
        do while (xfi < xf(it) .and. it >= 1)  
           it = it - 1  
        end do  
33    
34         ! Calcul de Xf(xi)      DOUBLE PRECISION xvrai(iim), Gvrai(iim)
35        ! intermediary variables because xlon and xprim are simple precision
36    
37         xi = xtild(it)      !------------------------------------------------------------------
38    
39         IF (it == 2 * nmax) THEN      it = 0 ! initial guess
40            it = 2 * nmax -1      h = twopi_d / iim
        END IF  
41    
42         CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), &      DO i = 1, iim
43              xtild(it), xtild(it + 1), a0, a1, a2, a3)         Y = - pi_d + (i + xuv - 0.75d0) * h
44         Xf1 = Xf(it)         ! - pi <= y < pi
45         Xprimin = a1 + xi * (2d0 * a2 + xi * 3d0 * a3)         abs_y = abs(y)
46         xo1 = xi  
47         iter = 1         call hunt(xf, abs_y, it, my_lbound = 0)
48           ! {0 <= it <= nmax - 1}
49         do  
50            xi = xi - (Xf1 - Xfi) / Xprimin         ! Calcul de xvrai(i) et Gvrai(i)
51            IF (ABS(xi - xo1) <= my_eps .or. iter == 300) exit         CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &
52            xo1 = xi              xtild(it + 1))
53            Xf1 = a0 + xi * (a1 + xi * (a2 + xi * a3))         xvrai(i) = rtsafe(funcd, xtild(it), xtild(it + 1), xacc = 1d-6)
54            Xprimin = a1 + xi * (2d0 * a2 + xi * 3d0 * a3)         Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
55         end DO         if (y < 0d0) xvrai(i) = - xvrai(i)
   
        if (ABS(xi - xo1) > my_eps) then  
           ! iter == 300  
           print *, 'Pas de solution.'  
           print *, i, xfi  
           STOP 1  
        end if  
   
        xxprim(i) = twopi_d / (REAL(iim) * Xprimin)  
        xvrai(i) = xi + xzoom  
56      end DO      end DO
57    
     IF (ik == 1 .and. grossismx == 1.) THEN  
        xvrai(1) = xvrai(iim + 1)-twopi_d  
        xxprim(1) = xxprim(iim + 1)  
     END IF  
   
58      DO i = 1, iim -1      DO i = 1, iim -1
59         IF (xvrai(i + 1) < xvrai(i)) THEN         IF (xvrai(i + 1) < xvrai(i)) THEN
60            print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'            print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'
# Line 94  contains Line 62  contains
62         END IF         END IF
63      END DO      END DO
64    
65      DO i = 1, iim      xlon = xvrai + clon
66         xlon(i) = xvrai(i)      xprim = h / Gvrai
67         xprimm(i) = xxprim(i)  
68      END DO    end subroutine invert_zoom_x
69    
70      !**********************************************************************
71    
72      SUBROUTINE funcd(x, fval, fderiv)
73    
74        use coefpoly_m, only: a0, a1, a2, a3
75    
76        DOUBLE PRECISION, INTENT(IN):: x
77        DOUBLE PRECISION, INTENT(OUT):: fval, fderiv
78    
79        fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y
80        fderiv = a1 + x * (2d0 * a2 + x * 3d0 * a3)
81    
82    end subroutine fxhyp_loop_ik    END SUBROUTINE funcd
83    
84  end module fxhyp_loop_ik_m  end module invert_zoom_x_m

Legend:
Removed from v.124  
changed lines
  Added in v.154

  ViewVC Help
Powered by ViewVC 1.1.21