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 |
|
|
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 |
|
|
33 |
real xo1, Xf1 |
DOUBLE PRECISION 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 |
|
|
36 |
!------------------------------------------------------------------ |
!------------------------------------------------------------------ |
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, xtild(it), xtild(it + 1), xacc = 1d-6) |
52 |
Xf1 = Xf(it) |
Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3) |
|
Gvrai(i) = G(it) |
|
|
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 |
|
|
|
|
53 |
if (y < 0d0) xvrai(i) = - xvrai(i) |
if (y < 0d0) xvrai(i) = - xvrai(i) |
54 |
end DO |
end DO |
55 |
|
|
65 |
|
|
66 |
end subroutine invert_zoom_x |
end subroutine invert_zoom_x |
67 |
|
|
68 |
|
!********************************************************************** |
69 |
|
|
70 |
|
SUBROUTINE funcd(x, fval, fderiv) |
71 |
|
|
72 |
|
use coefpoly_m, only: a0, a1, a2, a3 |
73 |
|
|
74 |
|
DOUBLE PRECISION, INTENT(IN):: x |
75 |
|
DOUBLE PRECISION, INTENT(OUT):: fval, fderiv |
76 |
|
|
77 |
|
fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y |
78 |
|
fderiv = a1 + x * (2d0 * a2 + x * 3d0 * a3) |
79 |
|
|
80 |
|
END SUBROUTINE funcd |
81 |
|
|
82 |
end module invert_zoom_x_m |
end module invert_zoom_x_m |