/[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

revision 145 by guez, Tue Jun 16 15:23:29 2015 UTC revision 167 by guez, Mon Aug 24 16:30:33 2015 UTC
# Line 3  module invert_zoom_x_m Line 3  module invert_zoom_x_m
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 invert_zoom_x(xf, xtild, Xprimt, xlon, xprimm, xuv)    subroutine invert_zoom_x(beta, xf, xtild, G, xlon, xprim, 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      use dynetat0_m, only: clon, grossismx
17      use nr_util, only: pi_d, twopi_d      use nr_util, only: pi_d, twopi_d
18      use numer_rec_95, only: hunt      use numer_rec_95, only: hunt, rtsafe
19    
20        DOUBLE PRECISION, intent(in):: beta, Xf(0:), xtild(0:), G(0:) ! (0:nmax)
21    
22      DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), Xprimt(0:) ! (0:2 * nmax)      real, intent(out):: xlon(:), xprim(:) ! (iim)
     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      ! between - 0.25 and 0.5
# Line 23  contains Line 27  contains
27      ! 0.5 si calcul aux points U      ! 0.5 si calcul aux points U
28    
29      ! Local:      ! Local:
30      DOUBLE PRECISION xo1, Xfi, a0, a1, a2, a3, Xf1, Xprimin      DOUBLE PRECISION Y
31      integer i, it, iter      DOUBLE PRECISION h ! step of the uniform grid
32      DOUBLE PRECISION, parameter:: my_eps = 1d-6      integer i, it
33    
34      DOUBLE PRECISION xxprim(iim), xvrai(iim)      DOUBLE PRECISION xvrai(iim), Gvrai(iim)
35      ! intermediary variables because xlon and xprimm are simple precision      ! intermediary variables because xlon and xprim are simple precision
36    
37      !------------------------------------------------------------------      !------------------------------------------------------------------
38    
39        print *, "Call sequence information: invert_zoom_x"
40      it = 0 ! initial guess      it = 0 ! initial guess
41        h = twopi_d / iim
42    
43      DO i = 1, iim      DO i = 1, iim
44         Xfi = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim         Y = - pi_d + (i + xuv - 0.75d0) * h
45         ! - pi <= xfi < pi         ! - pi <= y < pi
46           abs_y = abs(y)
47         call hunt(xf, xfi, it, my_lbound = 0)  
48         it = max(0, it)         ! Distinguish boundaries in order to avoid roundoff error.
49         ! In principle, xfi >= xf(0), max with 0 just in case of         ! funcd should be exactly equal to 0 at xtild(it) or xtild(it +
50         ! roundoff error. {0 <= it <= 2 * nmax - 1}         ! 1) and could be very small with the wrong sign so rtsafe
51           ! would fail.
52         ! Calcul de Xf(xvrai(i))         if (abs_y == 0d0) then
53              xvrai(i) = 0d0
54         CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), &            gvrai(i) = grossismx
55              xtild(it), xtild(it + 1), a0, a1, a2, a3)         else if (abs_y == pi_d) then
56         xvrai(i) = xtild(it)            xvrai(i) = pi_d
57         Xf1 = Xf(it)            gvrai(i) = 2d0 * beta - grossismx
58         Xprimin = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)         else
59         xo1 = xvrai(i)            call hunt(xf, abs_y, it, my_lbound = 0)
60         iter = 1            ! {0 <= it <= nmax - 1}
61    
62         do            ! Calcul de xvrai(i) et Gvrai(i)
63            xvrai(i) = xvrai(i) - (Xf1 - Xfi) / Xprimin            CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &
64            IF (ABS(xvrai(i) - xo1) <= my_eps .or. iter == 300) exit                 xtild(it + 1))
65            xo1 = xvrai(i)            xvrai(i) = rtsafe(funcd, xtild(it), xtild(it + 1), xacc = 1d-6)
66            Xf1 = a0 + xvrai(i) * (a1 + xvrai(i) * (a2 + xvrai(i) * a3))            Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
           Xprimin = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)  
        end DO  
   
        if (ABS(xvrai(i) - xo1) > my_eps) then  
           ! iter == 300  
           print *, 'Pas de solution.'  
           print *, i, xfi  
           STOP 1  
67         end if         end if
68    
69         xxprim(i) = twopi_d / (iim * Xprimin)         if (y < 0d0) xvrai(i) = - xvrai(i)
70      end DO      end DO
71    
72      DO i = 1, iim -1      DO i = 1, iim -1
# Line 79  contains Line 77  contains
77      END DO      END DO
78    
79      xlon = xvrai + clon      xlon = xvrai + clon
80      xprimm = xxprim      xprim = h / Gvrai
81    
82    end subroutine invert_zoom_x    end subroutine invert_zoom_x
83    
84      !**********************************************************************
85    
86      SUBROUTINE funcd(x, fval, fderiv)
87    
88        use coefpoly_m, only: a0, a1, a2, a3
89    
90        DOUBLE PRECISION, INTENT(IN):: x
91        DOUBLE PRECISION, INTENT(OUT):: fval, fderiv
92    
93        fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y
94        fderiv = a1 + x * (2d0 * a2 + x * 3d0 * a3)
95    
96      END SUBROUTINE funcd
97    
98  end module invert_zoom_x_m  end module invert_zoom_x_m

Legend:
Removed from v.145  
changed lines
  Added in v.167

  ViewVC Help
Powered by ViewVC 1.1.21