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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 148 - (show annotations)
Wed Jun 17 16:40:24 2015 UTC (9 years ago) by guez
File size: 2210 byte(s)
Renamed variables in fxhyp and invert_zoom_x to be closer to external
documentation. Changed from double precision to real in invert_zoom_x:
I do not see the need for double precision and double precision is
annoying because I want to use rtsafe without overloading it in
Numer_Rec_95.

Results are changed at the numerical noise level.

1 module invert_zoom_x_m
2
3 implicit none
4
5 INTEGER, PARAMETER:: nmax = 30000
6
7 contains
8
9 subroutine invert_zoom_x(xf, xtild, G, xlon, xprim, xuv)
10
11 use coefpoly_m, only: coefpoly
12 USE dimens_m, ONLY: iim
13 use dynetat0_m, only: clon
14 use nr_util, only: pi_d, twopi_d
15 use numer_rec_95, only: hunt
16
17 DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), G(0:) ! (0:nmax)
18
19 real, intent(out):: xlon(:), xprim(:) ! (iim)
20
21 DOUBLE PRECISION, intent(in):: xuv
22 ! between - 0.25 and 0.5
23 ! 0. si calcul aux points scalaires
24 ! 0.5 si calcul aux points U
25
26 ! Local:
27 DOUBLE PRECISION Y, abs_y, a0, a1, a2, a3
28 integer i, it, iter
29 real, parameter:: my_eps = 1e-6
30
31 real xo1, Xf1
32 real xvrai(iim), Gvrai(iim)
33 ! intermediary variables because xlon and xprim are simple precision
34
35 !------------------------------------------------------------------
36
37 it = 0 ! initial guess
38
39 DO i = 1, iim
40 Y = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim
41 ! - pi <= y < pi
42 abs_y = abs(y)
43
44 call hunt(xf, abs_y, it, my_lbound = 0)
45 ! {0 <= it <= nmax - 1}
46
47 ! Calcul de xvrai(i) et Gvrai(i)
48
49 CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &
50 xtild(it + 1), a0, a1, a2, a3)
51 xvrai(i) = xtild(it)
52 Xf1 = Xf(it)
53 Gvrai(i) = G(it)
54 xo1 = xvrai(i)
55 iter = 1
56
57 do
58 xvrai(i) = xvrai(i) - (Xf1 - abs_y) / Gvrai(i)
59 IF (ABS(xvrai(i) - xo1) <= my_eps .or. iter == 300) exit
60 xo1 = xvrai(i)
61 Xf1 = a0 + xvrai(i) * (a1 + xvrai(i) * (a2 + xvrai(i) * a3))
62 Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
63 end DO
64
65 if (ABS(xvrai(i) - xo1) > my_eps) then
66 ! iter == 300
67 print *, 'Pas de solution.'
68 print *, i, abs_y
69 STOP 1
70 end if
71
72 if (y < 0d0) xvrai(i) = - xvrai(i)
73 end DO
74
75 DO i = 1, iim -1
76 IF (xvrai(i + 1) < xvrai(i)) THEN
77 print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'
78 STOP 1
79 END IF
80 END DO
81
82 xlon = xvrai + clon
83 xprim = twopi_d / (iim * Gvrai)
84
85 end subroutine invert_zoom_x
86
87 end module invert_zoom_x_m

  ViewVC Help
Powered by ViewVC 1.1.21