source: LMDZ6/trunk/libf/dyn3d_common/invert_zoom_x_m.f90 @ 5348

Last change on this file since 5348 was 5271, checked in by abarral, 5 weeks ago

Move dimensions.h into a module
Nb: doesn't compile yet

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