--- trunk/Sources/dyn3d/invert_zoom_x.f 2015/06/17 16:40:24 148 +++ trunk/dyn3d/invert_zoom_x.f 2018/03/20 09:35:59 265 @@ -3,18 +3,21 @@ implicit none INTEGER, PARAMETER:: nmax = 30000 + DOUBLE PRECISION abs_y + + private abs_y, funcd contains - subroutine invert_zoom_x(xf, xtild, G, xlon, xprim, xuv) + subroutine invert_zoom_x(beta, xf, xtild, G, xlon, xprim, xuv) - use coefpoly_m, only: coefpoly - USE dimens_m, ONLY: iim - use dynetat0_m, only: clon + use coefpoly_m, only: coefpoly, a1, a2, a3 + USE dimensions, ONLY: iim + use dynetat0_m, only: clon, grossismx 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) + DOUBLE PRECISION, intent(in):: beta, Xf(0:), xtild(0:), G(0:) ! (0:nmax) real, intent(out):: xlon(:), xprim(:) ! (iim) @@ -24,49 +27,43 @@ ! 0.5 si calcul aux points U ! Local: - DOUBLE PRECISION Y, abs_y, a0, a1, a2, a3 - integer i, it, iter - real, parameter:: my_eps = 1e-6 - - real xo1, Xf1 - real xvrai(iim), Gvrai(iim) - ! intermediary variables because xlon and xprim are simple precision + DOUBLE PRECISION Y + DOUBLE PRECISION h ! step of the uniform grid + integer i, it + + DOUBLE PRECISION xvrai(iim), Gvrai(iim) + ! intermediary variables because xlon and xprim are single precision !------------------------------------------------------------------ + print *, "Call sequence information: invert_zoom_x" it = 0 ! initial guess + h = twopi_d / iim DO i = 1, iim - Y = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim + Y = - pi_d + (i + xuv - 0.75d0) * h ! - pi <= y < pi abs_y = abs(y) - call hunt(xf, abs_y, it, my_lbound = 0) - ! {0 <= it <= nmax - 1} - - ! 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) - 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)) + ! Distinguish boundaries in order to avoid roundoff error. + ! funcd should be exactly equal to 0 at xtild(it) or xtild(it + + ! 1) and could be very small with the wrong sign so rtsafe + ! would fail. + if (abs_y == 0d0) then + xvrai(i) = 0d0 + gvrai(i) = grossismx + else if (abs_y == pi_d) then + xvrai(i) = pi_d + gvrai(i) = 2d0 * beta - grossismx + else + call hunt(xf, abs_y, it, my_lbound = 0) + ! {0 <= it <= nmax - 1} + + ! Calcul de xvrai(i) et Gvrai(i) + CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), & + 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) - end DO - - if (ABS(xvrai(i) - xo1) > my_eps) then - ! iter == 300 - print *, 'Pas de solution.' - print *, i, abs_y - STOP 1 end if if (y < 0d0) xvrai(i) = - xvrai(i) @@ -80,8 +77,22 @@ END DO xlon = xvrai + clon - xprim = twopi_d / (iim * Gvrai) + xprim = h / 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