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