--- trunk/Sources/dyn3d/invert_zoom_x.f 2015/06/16 17:27:33 146 +++ trunk/Sources/dyn3d/invert_zoom_x.f 2015/06/18 13:49:26 150 @@ -3,20 +3,23 @@ implicit none INTEGER, PARAMETER:: nmax = 30000 + DOUBLE PRECISION abs_y + + private abs_y, funcd contains - subroutine invert_zoom_x(xf, xtild, G, xlon, xprimm, xuv) + subroutine invert_zoom_x(xf, xtild, G, xlon, xprim, xuv) - use coefpoly_m, only: coefpoly + use coefpoly_m, only: coefpoly, a1, a2, a3 USE dimens_m, ONLY: iim use dynetat0_m, only: clon use nr_util, only: pi_d, twopi_d - use numer_rec_95, only: hunt + use numer_rec_95, only: hunt, rtsafe DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), G(0:) ! (0:nmax) - real, intent(out):: xlon(:), xprimm(:) ! (iim) + real, intent(out):: xlon(:), xprim(:) ! (iim) DOUBLE PRECISION, intent(in):: xuv ! between - 0.25 and 0.5 @@ -24,52 +27,30 @@ ! 0.5 si calcul aux points U ! Local: - DOUBLE PRECISION xo1, Xfi, abs_xfi, a0, a1, a2, a3, Xf1, Xprimin - integer i, it, iter - DOUBLE PRECISION, parameter:: my_eps = 1d-6 + DOUBLE PRECISION Y + integer i, it - DOUBLE PRECISION xxprim(iim), xvrai(iim) - ! intermediary variables because xlon and xprimm are simple precision + DOUBLE PRECISION xvrai(iim), Gvrai(iim) + ! intermediary variables because xlon and xprim are simple precision !------------------------------------------------------------------ it = 0 ! initial guess DO i = 1, iim - Xfi = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim - ! - pi <= xfi < pi - abs_xfi = abs(xfi) + Y = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim + ! - pi <= y < pi + abs_y = abs(y) - call hunt(xf, abs_xfi, it, my_lbound = 0) + call hunt(xf, abs_y, it, my_lbound = 0) ! {0 <= it <= nmax - 1} - ! Calcul de xvrai(i) et xxprim(i) - + ! Calcul de xvrai(i) et Gvrai(i) CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), & - xtild(it + 1), a0, a1, a2, a3) - xvrai(i) = xtild(it) - Xf1 = Xf(it) - Xprimin = G(it) - xo1 = xvrai(i) - iter = 1 - - do - xvrai(i) = xvrai(i) - (Xf1 - abs_xfi) / Xprimin - 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)) - Xprimin = 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_xfi - STOP 1 - end if - - if (xfi < 0) xvrai(i) = - xvrai(i) - xxprim(i) = twopi_d / (iim * Xprimin) + xtild(it + 1)) + xvrai(i) = rtsafe(funcd, xtild(it), xtild(it + 1), xacc = 1d-6) + Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3) + if (y < 0d0) xvrai(i) = - xvrai(i) end DO DO i = 1, iim -1 @@ -80,8 +61,22 @@ END DO xlon = xvrai + clon - xprimm = xxprim + xprim = twopi_d / (iim * Gvrai) end subroutine invert_zoom_x + !********************************************************************** + + SUBROUTINE funcd(x, fval, fderiv) + + use coefpoly_m, only: a0, a1, a2, a3 + + DOUBLE PRECISION, INTENT(IN):: x + DOUBLE PRECISION, INTENT(OUT):: fval, fderiv + + fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y + fderiv = a1 + x * (2d0 * a2 + x * 3d0 * a3) + + END SUBROUTINE funcd + end module invert_zoom_x_m