/[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 149 - (show annotations)
Thu Jun 18 12:23:44 2015 UTC (9 years ago) by guez
File size: 2090 byte(s)
In invert_zoom_x, call rtsafe instead of the equivalent coding that
was there. funcd needs to access a[0-4] and abs_y so we upgrade a[0-4]
from arguments of coefpoly to variables of module coefpoly_m and abs_y
from local variable of invert_zoom_x to private variable of module
invert_zoom_x_m.

Removed unused arguments t10m and q10m of hbtm.

1 module invert_zoom_x_m
2
3 implicit none
4
5 INTEGER, PARAMETER:: nmax = 30000
6 DOUBLE PRECISION abs_y
7
8 private abs_y, funcd
9
10 contains
11
12 subroutine invert_zoom_x(xf, xtild, G, xlon, xprim, xuv)
13
14 use coefpoly_m, only: coefpoly, a1, a2, a3
15 USE dimens_m, ONLY: iim
16 use dynetat0_m, only: clon
17 use nr_util, only: pi_d, twopi_d
18 use numer_rec_95, only: hunt, rtsafe
19
20 DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), G(0:) ! (0:nmax)
21
22 real, intent(out):: xlon(:), xprim(:) ! (iim)
23
24 DOUBLE PRECISION, intent(in):: xuv
25 ! between - 0.25 and 0.5
26 ! 0. si calcul aux points scalaires
27 ! 0.5 si calcul aux points U
28
29 ! Local:
30 DOUBLE PRECISION Y
31 integer i, it
32
33 real xvrai(iim), Gvrai(iim)
34 ! intermediary variables because xlon and xprim are simple precision
35
36 !------------------------------------------------------------------
37
38 it = 0 ! initial guess
39
40 DO i = 1, iim
41 Y = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim
42 ! - pi <= y < pi
43 abs_y = abs(y)
44
45 call hunt(xf, abs_y, it, my_lbound = 0)
46 ! {0 <= it <= nmax - 1}
47
48 ! Calcul de xvrai(i) et Gvrai(i)
49 CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &
50 xtild(it + 1))
51 xvrai(i) = rtsafe(funcd, real(xtild(it)), real(xtild(it + 1)), &
52 xacc = 1e-6)
53 Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
54 if (y < 0d0) xvrai(i) = - xvrai(i)
55 end DO
56
57 DO i = 1, iim -1
58 IF (xvrai(i + 1) < xvrai(i)) THEN
59 print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'
60 STOP 1
61 END IF
62 END DO
63
64 xlon = xvrai + clon
65 xprim = twopi_d / (iim * Gvrai)
66
67 end subroutine invert_zoom_x
68
69 !**********************************************************************
70
71 SUBROUTINE funcd(x, fval, fderiv)
72
73 use coefpoly_m, only: a0, a1, a2, a3
74
75 REAL, INTENT(IN):: x
76 REAL, INTENT(OUT):: fval, fderiv
77
78 fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y
79 fderiv = a1 + x * (2. * a2 + x * 3. * a3)
80
81 END SUBROUTINE funcd
82
83 end module invert_zoom_x_m

  ViewVC Help
Powered by ViewVC 1.1.21