1 |
module fxhyp_loop_ik_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
INTEGER, PARAMETER:: nmax = 30000 |
6 |
|
7 |
contains |
8 |
|
9 |
subroutine fxhyp_loop_ik(ik, decalx, xf, xtild, Xprimt, xzoom, xlon, & |
10 |
xprimm, xuv) |
11 |
|
12 |
use coefpoly_m, only: coefpoly |
13 |
USE dimens_m, ONLY: iim |
14 |
use nr_util, only: pi_d, twopi_d, twopi |
15 |
use serre, only: grossismx |
16 |
|
17 |
INTEGER, intent(in):: ik |
18 |
DOUBLE PRECISION, intent(in):: decalx |
19 |
DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), Xprimt(0:) ! (0:2 * nmax) |
20 |
DOUBLE PRECISION, intent(in):: xzoom |
21 |
real, intent(out):: xlon(:), xprimm(:) ! (iim + 1) |
22 |
|
23 |
DOUBLE PRECISION, intent(in):: xuv |
24 |
! 0. si calcul aux points scalaires |
25 |
! 0.5 si calcul aux points U |
26 |
|
27 |
! Local: |
28 |
|
29 |
DOUBLE PRECISION xo1, Xfi, xi, a0, a1, a2, a3, Xf1, Xprimin |
30 |
integer ii1, ii2, i, it, iter, idif, ii |
31 |
DOUBLE PRECISION, parameter:: my_eps = 1e-6 |
32 |
DOUBLE PRECISION xxprim(iim + 1), xvrai(iim + 1) |
33 |
INTEGER:: is2 = 0 |
34 |
|
35 |
!------------------------------------------------------------------ |
36 |
|
37 |
IF (ik == 1 .and. grossismx == 1.) THEN |
38 |
ii1 = 2 |
39 |
ii2 = iim + 1 |
40 |
else |
41 |
ii1=1 |
42 |
ii2=iim |
43 |
END IF |
44 |
|
45 |
DO i = ii1, ii2 |
46 |
Xfi = - pi_d + (i + xuv - decalx) * twopi_d / iim |
47 |
|
48 |
it = 2 * nmax |
49 |
do while (xfi < xf(it) .and. it >= 1) |
50 |
it = it - 1 |
51 |
end do |
52 |
|
53 |
! Calcul de Xf(xi) |
54 |
|
55 |
xi = xtild(it) |
56 |
|
57 |
IF (it == 2 * nmax) THEN |
58 |
it = 2 * nmax -1 |
59 |
END IF |
60 |
|
61 |
CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), & |
62 |
xtild(it), xtild(it + 1), a0, a1, a2, a3) |
63 |
Xf1 = Xf(it) |
64 |
Xprimin = a1 + xi * (2d0 * a2 + xi * 3d0 * a3) |
65 |
xo1 = xi |
66 |
iter = 1 |
67 |
|
68 |
do |
69 |
xi = xi - (Xf1 - Xfi) / Xprimin |
70 |
IF (ABS(xi - xo1) <= my_eps .or. iter == 300) exit |
71 |
xo1 = xi |
72 |
Xf1 = a0 + xi * (a1 + xi * (a2 + xi * a3)) |
73 |
Xprimin = a1 + xi * (2d0 * a2 + xi * 3d0 * a3) |
74 |
end DO |
75 |
|
76 |
if (ABS(xi - xo1) > my_eps) then |
77 |
! iter == 300 |
78 |
print *, 'Pas de solution.' |
79 |
print *, i, xfi |
80 |
STOP 1 |
81 |
end if |
82 |
|
83 |
xxprim(i) = twopi_d / (REAL(iim) * Xprimin) |
84 |
xvrai(i) = xi + xzoom |
85 |
end DO |
86 |
|
87 |
IF (ik == 1 .and. grossismx == 1.) THEN |
88 |
xvrai(1) = xvrai(iim + 1)-twopi_d |
89 |
xxprim(1) = xxprim(iim + 1) |
90 |
END IF |
91 |
|
92 |
DO i = 1, iim |
93 |
xlon(i) = xvrai(i) |
94 |
xprimm(i) = xxprim(i) |
95 |
END DO |
96 |
|
97 |
DO i = 1, iim -1 |
98 |
IF (xvrai(i + 1) < xvrai(i)) THEN |
99 |
print *, 'rlonu(', i + 1, ') < rlonu(', i, ')' |
100 |
STOP 1 |
101 |
END IF |
102 |
END DO |
103 |
|
104 |
IF (ik == 1 .and. (MINval(xvrai(:iim)) < - pi_d - 0.1d0 & |
105 |
.or. MAXval(xvrai(:iim)) > pi_d + 0.1d0)) THEN |
106 |
IF (xzoom <= 0.) THEN |
107 |
i = 1 |
108 |
|
109 |
do while (xvrai(i) < - pi_d .and. i < iim) |
110 |
i = i + 1 |
111 |
end do |
112 |
|
113 |
if (xvrai(i) < - pi_d) then |
114 |
print *, 'Xvrai plus petit que - pi !' |
115 |
STOP 1 |
116 |
end if |
117 |
|
118 |
is2 = i |
119 |
ELSE |
120 |
i = iim |
121 |
|
122 |
do while (xvrai(i) > pi_d .and. i > 1) |
123 |
i = i - 1 |
124 |
end do |
125 |
|
126 |
if (xvrai(i) > pi_d) then |
127 |
print *, 'Xvrai plus grand que pi !' |
128 |
STOP 1 |
129 |
end if |
130 |
|
131 |
is2 = i |
132 |
END IF |
133 |
END IF |
134 |
|
135 |
if (is2 /= 0) then |
136 |
IF (xzoom <= 0.) THEN |
137 |
IF (is2 /= 1) THEN |
138 |
DO ii = is2, iim |
139 |
xlon(ii-is2 + 1) = xvrai(ii) |
140 |
xprimm(ii-is2 + 1) = xxprim(ii) |
141 |
END DO |
142 |
DO ii = 1, is2 -1 |
143 |
xlon(ii + iim-is2 + 1) = xvrai(ii) + twopi_d |
144 |
xprimm(ii + iim-is2 + 1) = xxprim(ii) |
145 |
END DO |
146 |
END IF |
147 |
else |
148 |
idif = iim -is2 |
149 |
|
150 |
DO ii = 1, is2 |
151 |
xlon(ii + idif) = xvrai(ii) |
152 |
xprimm(ii + idif) = xxprim(ii) |
153 |
END DO |
154 |
|
155 |
DO ii = 1, idif |
156 |
xlon(ii) = xvrai(ii + is2) - twopi_d |
157 |
xprimm(ii) = xxprim(ii + is2) |
158 |
END DO |
159 |
end IF |
160 |
end if |
161 |
|
162 |
xlon(iim + 1) = xlon(1) + twopi |
163 |
xprimm(iim + 1) = xprimm(1) |
164 |
|
165 |
end subroutine fxhyp_loop_ik |
166 |
|
167 |
end module fxhyp_loop_ik_m |