/[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 121 - (show annotations)
Wed Jan 28 16:10:02 2015 UTC (9 years, 3 months ago) by guez
Original Path: trunk/dyn3d/fxhyp_loop_ik.f
File size: 4388 byte(s)
In procedure fxhyp, extracted the body of the loop on ik into a new
procedure:  fxhyp_loop_ik.

dzoomx and dzoomy must now be fractions of the entire range, they
cannot be ranges in degrees or rad.

In fxhyp, force Xf(2 * nmax) = pi_d instead of possibly doing it in
fxhyp_loop_ik.

In fxhyp_loop_ik, when testing whether xvrai is between -pi and pi,
changed the boundaries from -pi - 0.1 to - pi_d - 1d-5 and from pi +
0.1 to pi_d + 1d-5. This reveals a misconception of the
code. Therefore, this version does not work.

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, xi2
30 integer ii1, ii2, i, it, iter, idif, ii
31 DOUBLE PRECISION, parameter:: my_eps = 1e-3
32 DOUBLE PRECISION xxprim(iim + 1), xvrai(iim + 1)
33 INTEGER:: is2 = 0
34
35 !------------------------------------------------------------------
36
37 xo1 = 0.
38
39 IF (ik == 1 .and. grossismx == 1.) THEN
40 ii1 = 2
41 ii2 = iim + 1
42 else
43 ii1=1
44 ii2=iim
45 END IF
46
47 DO i = ii1, ii2
48 Xfi = - pi_d + (i + xuv - decalx) * twopi_d / iim
49
50 it = 2 * nmax
51 do while (xfi < xf(it) .and. it >= 1)
52 it = it - 1
53 end do
54
55 ! Calcul de Xf(xi)
56
57 xi = xtild(it)
58
59 IF (it == 2 * nmax) THEN
60 it = 2 * nmax -1
61 END IF
62
63 ! Appel de la routine qui calcule les coefficients a0, a1,
64 ! a2, a3 d'un polynome de degre 3 qui passe par les points
65 ! (Xf(it), xtild(it)) et (Xf(it + 1), xtild(it + 1))
66
67 CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), &
68 xtild(it), xtild(it + 1), a0, a1, a2, a3)
69
70 Xf1 = Xf(it)
71 Xprimin = a1 + 2. * a2 * xi + 3. * a3 * xi**2
72
73 iter = 1
74
75 do
76 xi = xi - (Xf1 - Xfi) / Xprimin
77 IF (ABS(xi - xo1) <= my_eps .or. iter == 300) exit
78 xo1 = xi
79 xi2 = xi * xi
80 Xf1 = a0 + a1 * xi + a2 * xi2 + a3 * xi2 * xi
81 Xprimin = a1 + 2. * a2 * xi + 3. * a3 * xi2
82 end DO
83
84 if (ABS(xi - xo1) > my_eps) then
85 ! iter == 300
86 print *, 'Pas de solution.'
87 print *, i, xfi
88 STOP 1
89 end if
90
91 xxprim(i) = twopi_d / (REAL(iim) * Xprimin)
92 xvrai(i) = xi + xzoom
93 end DO
94
95 IF (ik == 1 .and. grossismx == 1.) THEN
96 xvrai(1) = xvrai(iim + 1)-twopi_d
97 xxprim(1) = xxprim(iim + 1)
98 END IF
99
100 DO i = 1, iim
101 xlon(i) = xvrai(i)
102 xprimm(i) = xxprim(i)
103 END DO
104
105 DO i = 1, iim -1
106 IF (xvrai(i + 1) < xvrai(i)) THEN
107 print *, 'rlonu(', i + 1, ') < rlonu(', i, ')'
108 STOP 1
109 END IF
110 END DO
111
112 IF (.not. (MINval(xvrai(:iim)) >= - pi_d - 1d-5 &
113 .and. MAXval(xvrai(:iim)) <= pi_d + 1d-5)) THEN
114 IF (xzoom <= 0.) THEN
115 IF (ik == 1) THEN
116 i = 1
117
118 do while (xvrai(i) < - pi_d .and. i < iim)
119 i = i + 1
120 end do
121
122 if (xvrai(i) < - pi_d) then
123 print *, 'Xvrai plus petit que - pi !'
124 STOP 1
125 end if
126
127 is2 = i
128 END IF
129
130 IF (is2 /= 1) THEN
131 DO ii = is2, iim
132 xlon(ii-is2 + 1) = xvrai(ii)
133 xprimm(ii-is2 + 1) = xxprim(ii)
134 END DO
135 DO ii = 1, is2 -1
136 xlon(ii + iim-is2 + 1) = xvrai(ii) + twopi_d
137 xprimm(ii + iim-is2 + 1) = xxprim(ii)
138 END DO
139 END IF
140 ELSE
141 IF (ik == 1) THEN
142 i = iim
143
144 do while (xvrai(i) > pi_d .and. i > 1)
145 i = i - 1
146 end do
147
148 if (xvrai(i) > pi_d) then
149 print *, 'Xvrai plus grand que pi !'
150 STOP 1
151 end if
152
153 is2 = i
154 END IF
155
156 idif = iim -is2
157
158 DO ii = 1, is2
159 xlon(ii + idif) = xvrai(ii)
160 xprimm(ii + idif) = xxprim(ii)
161 END DO
162
163 DO ii = 1, idif
164 xlon(ii) = xvrai(ii + is2) - twopi_d
165 xprimm(ii) = xxprim(ii + is2)
166 END DO
167 END IF
168 END IF
169
170 xlon(iim + 1) = xlon(1) + twopi
171 xprimm(iim + 1) = xprimm(1)
172
173 end subroutine fxhyp_loop_ik
174
175 end module fxhyp_loop_ik_m

  ViewVC Help
Powered by ViewVC 1.1.21