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

Diff of /trunk/dyn3d/invert_zoom_x.f

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

revision 148 by guez, Wed Jun 17 16:40:24 2015 UTC revision 149 by guez, Thu Jun 18 12:23:44 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, G, xlon, xprim, xuv)    subroutine invert_zoom_x(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
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):: Xf(0:), xtild(0:), G(0:) ! (0:nmax)      DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), G(0:) ! (0:nmax)
21    
# Line 24  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 Y, abs_y, a0, a1, a2, a3      DOUBLE PRECISION Y
31      integer i, it, iter      integer i, it
     real, parameter:: my_eps = 1e-6  
32    
     real xo1, Xf1  
33      real xvrai(iim), Gvrai(iim)      real xvrai(iim), Gvrai(iim)
34      ! intermediary variables because xlon and xprim are simple precision      ! intermediary variables because xlon and xprim are simple precision
35    
# Line 45  contains Line 46  contains
46         ! {0 <= it <= nmax - 1}         ! {0 <= it <= nmax - 1}
47    
48         ! Calcul de xvrai(i) et Gvrai(i)         ! Calcul de xvrai(i) et Gvrai(i)
   
49         CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &         CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &
50              xtild(it + 1), a0, a1, a2, a3)              xtild(it + 1))
51         xvrai(i) = xtild(it)         xvrai(i) = rtsafe(funcd, real(xtild(it)), real(xtild(it + 1)), &
52         Xf1 = Xf(it)              xacc = 1e-6)
53         Gvrai(i) = G(it)         Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
        xo1 = xvrai(i)  
        iter = 1  
   
        do  
           xvrai(i) = xvrai(i) - (Xf1 - abs_y) / Gvrai(i)  
           IF (ABS(xvrai(i) - xo1) <= my_eps .or. iter == 300) exit  
           xo1 = xvrai(i)  
           Xf1 = a0 + xvrai(i) * (a1 + xvrai(i) * (a2 + xvrai(i) * a3))  
           Gvrai(i) = 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, abs_y  
           STOP 1  
        end if  
   
54         if (y < 0d0) xvrai(i) = - xvrai(i)         if (y < 0d0) xvrai(i) = - xvrai(i)
55      end DO      end DO
56    
# Line 84  contains Line 66  contains
66    
67    end subroutine invert_zoom_x    end subroutine invert_zoom_x
68    
69      !**********************************************************************
70    
71      SUBROUTINE funcd(x, fval, fderiv)
72    
73        use coefpoly_m, only: a0, a1, a2, a3
74    
75        REAL, INTENT(IN):: x
76        REAL, INTENT(OUT):: fval, fderiv
77    
78        fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y
79        fderiv = a1 + x * (2. * a2 + x * 3. * a3)
80    
81      END SUBROUTINE funcd
82    
83  end module invert_zoom_x_m  end module invert_zoom_x_m

Legend:
Removed from v.148  
changed lines
  Added in v.149

  ViewVC Help
Powered by ViewVC 1.1.21