/[lmdze]/trunk/Sources/dyn3d/invert_zoom_x.f
ViewVC logotype

Annotation of /trunk/Sources/dyn3d/invert_zoom_x.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 145 - (hide annotations)
Tue Jun 16 15:23:29 2015 UTC (8 years, 11 months ago) by guez
File size: 2360 byte(s)
Renamed bibio to misc.

In procedure fxhyp, use the fact that xf is an odd function of xtild.

In procedure invert_zoom_x, replace linear search in xf by
bisection. Also, use result from previous loop iteration as initial
guess. Variable "it" cannot be equal to 2 * nmax after search.

Unused arguments: hm of cv3_feed; ph, qnk, tv,tvp of cv3_mixing; ppsol
of lw; rconst, temp of vdif_kcay; rconst, plev, temp, ustar, l_mix of
yamada.

1 guez 131 module invert_zoom_x_m
2 guez 121
3     implicit none
4    
5     INTEGER, PARAMETER:: nmax = 30000
6    
7     contains
8    
9 guez 131 subroutine invert_zoom_x(xf, xtild, Xprimt, xlon, xprimm, xuv)
10 guez 121
11     use coefpoly_m, only: coefpoly
12     USE dimens_m, ONLY: iim
13 guez 139 use dynetat0_m, only: clon
14 guez 124 use nr_util, only: pi_d, twopi_d
15 guez 145 use numer_rec_95, only: hunt
16 guez 121
17     DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), Xprimt(0:) ! (0:2 * nmax)
18 guez 124 real, intent(out):: xlon(:), xprimm(:) ! (iim)
19 guez 121
20     DOUBLE PRECISION, intent(in):: xuv
21 guez 145 ! between - 0.25 and 0.5
22 guez 121 ! 0. si calcul aux points scalaires
23 guez 145 ! 0.5 si calcul aux points U
24 guez 121
25     ! Local:
26 guez 126 DOUBLE PRECISION xo1, Xfi, a0, a1, a2, a3, Xf1, Xprimin
27     integer i, it, iter
28     DOUBLE PRECISION, parameter:: my_eps = 1d-6
29 guez 121
30 guez 126 DOUBLE PRECISION xxprim(iim), xvrai(iim)
31     ! intermediary variables because xlon and xprimm are simple precision
32    
33 guez 121 !------------------------------------------------------------------
34    
35 guez 145 it = 0 ! initial guess
36    
37 guez 126 DO i = 1, iim
38     Xfi = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim
39 guez 145 ! - pi <= xfi < pi
40 guez 121
41 guez 145 call hunt(xf, xfi, it, my_lbound = 0)
42     it = max(0, it)
43     ! In principle, xfi >= xf(0), max with 0 just in case of
44     ! roundoff error. {0 <= it <= 2 * nmax - 1}
45 guez 121
46 guez 126 ! Calcul de Xf(xvrai(i))
47 guez 121
48     CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), &
49     xtild(it), xtild(it + 1), a0, a1, a2, a3)
50 guez 145 xvrai(i) = xtild(it)
51 guez 121 Xf1 = Xf(it)
52 guez 126 Xprimin = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
53     xo1 = xvrai(i)
54 guez 121 iter = 1
55    
56     do
57 guez 126 xvrai(i) = xvrai(i) - (Xf1 - Xfi) / Xprimin
58     IF (ABS(xvrai(i) - xo1) <= my_eps .or. iter == 300) exit
59     xo1 = xvrai(i)
60     Xf1 = a0 + xvrai(i) * (a1 + xvrai(i) * (a2 + xvrai(i) * a3))
61     Xprimin = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
62 guez 121 end DO
63    
64 guez 126 if (ABS(xvrai(i) - xo1) > my_eps) then
65 guez 121 ! iter == 300
66     print *, 'Pas de solution.'
67     print *, i, xfi
68     STOP 1
69     end if
70    
71 guez 126 xxprim(i) = twopi_d / (iim * Xprimin)
72 guez 121 end DO
73    
74     DO i = 1, iim -1
75     IF (xvrai(i + 1) < xvrai(i)) THEN
76 guez 124 print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'
77 guez 121 STOP 1
78     END IF
79     END DO
80    
81 guez 127 xlon = xvrai + clon
82 guez 126 xprimm = xxprim
83 guez 121
84 guez 131 end subroutine invert_zoom_x
85 guez 121
86 guez 131 end module invert_zoom_x_m

  ViewVC Help
Powered by ViewVC 1.1.21