source: LMDZ5/branches/AI-cosp/libf/dyn3d_common/invert_zoom_x_m.F90 @ 5440

Last change on this file since 5440 was 2228, checked in by lguez, 10 years ago

Correcting a problem from revision 2218. The type double precision
with option "-fdefault-real-8" of gfortran is promoted to 16-byte
precision and there is no specific procedure in arth with this
precision. Could not add a specific procedure in arth with double
precision because, with ifort, the option "-real-size 64" does not
promote the double precision, so that would make two identical
specific procedures in arth.

In module nrtype, replaced double precision by a parameterized real
kind so that the effective precision does not depend on a compiler
option.

In coefpoly, fxhyp, fyhyp and invert_zoom_x, use the parameterized
real kind defined in nrtype, instead of double precision.

Also, in module nrtype, removed unused derived types sprs2_sp and
sprs2_dp.

File size: 2.2 KB
Line 
1module invert_zoom_x_m
2
3  implicit none
4
5  INTEGER, PARAMETER:: nmax = 30000
6
7contains
8
9  subroutine invert_zoom_x(xf, xtild, Xprimt, xlon, xprimm, xuv)
10
11    use coefpoly_m, only: coefpoly
12    use nrtype, only: pi, pi_d, twopi_d, k8
13
14    include "dimensions.h"
15    ! for iim
16
17    include "serre.h"
18    ! for clon
19
20    REAL(K8), intent(in):: Xf(0:), xtild(0:), Xprimt(0:) ! (0:2 * nmax)
21    real, intent(out):: xlon(:), xprimm(:) ! (iim)
22
23    REAL(K8), intent(in):: xuv
24    ! 0. si calcul aux points scalaires
25    ! 0.5 si calcul aux points U
26
27    ! Local:
28    REAL(K8) xo1, Xfi, a0, a1, a2, a3, Xf1, Xprimin
29    integer i, it, iter
30    REAL(K8), parameter:: my_eps = 1e-6_k8
31
32    REAL(K8) xxprim(iim), xvrai(iim)
33    ! intermediary variables because xlon and xprimm are simple precision
34
35    !------------------------------------------------------------------
36
37    DO i = 1, iim
38       Xfi = - pi_d + (i + xuv - 0.75_k8) * twopi_d / iim
39
40       it = 2 * nmax
41       do while (xfi < xf(it) .and. it >= 1)
42          it = it - 1
43       end do
44
45       ! Calcul de Xf(xvrai(i))
46
47       xvrai(i) = xtild(it)
48
49       IF (it == 2 * nmax) THEN
50          it = 2 * nmax -1
51       END IF
52
53       CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), &
54            xtild(it), xtild(it + 1), a0, a1, a2, a3)
55       Xf1 = Xf(it)
56       Xprimin = a1 + xvrai(i) * (2._k8 * a2 + xvrai(i) * 3._k8 * a3)
57       xo1 = xvrai(i)
58       iter = 1
59
60       do
61          xvrai(i) = xvrai(i) - (Xf1 - Xfi) / Xprimin
62          IF (ABS(xvrai(i) - xo1) <= my_eps .or. iter == 300) exit
63          xo1 = xvrai(i)
64          Xf1 = a0 + xvrai(i) * (a1 + xvrai(i) * (a2 + xvrai(i) * a3))
65          Xprimin = a1 + xvrai(i) * (2._k8 * a2 + xvrai(i) * 3._k8 * a3)
66       end DO
67
68       if (ABS(xvrai(i) - xo1) > my_eps) then
69          ! iter == 300
70          print *, 'Pas de solution.'
71          print *, i, xfi
72          STOP 1
73       end if
74
75       xxprim(i) = twopi_d / (iim * Xprimin)
76    end DO
77
78    DO i = 1, iim -1
79       IF (xvrai(i + 1) < xvrai(i)) THEN
80          print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'
81          STOP 1
82       END IF
83    END DO
84
85    xlon = xvrai + clon / 180. * pi
86    xprimm = xxprim
87
88  end subroutine invert_zoom_x
89
90end module invert_zoom_x_m
Note: See TracBrowser for help on using the repository browser.