LMDZ
invert_zoom_x_m.F90
Go to the documentation of this file.
2 
3  implicit none
4 
5  INTEGER, PARAMETER:: nmax = 30000
6 
7 contains
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 
90 end module invert_zoom_x_m
subroutine coefpoly(xf1, xf2, xprim1, xprim2, xtild1, xtild2, a0, a1, a2, a3)
Definition: coefpoly_m.F90:8
!$Header!c!c!c include serre h!c REAL clon
Definition: serre.h:8
integer, parameter k8
Definition: nrtype.F90:5
real(k8), parameter twopi_d
Definition: nrtype.F90:19
!$Id mode_top_bound COMMON comconstr && pi
Definition: comconst.h:7
subroutine invert_zoom_x(xf, xtild, Xprimt, xlon, xprimm, xuv)
Definition: nrtype.F90:1
integer, parameter nmax
c c zjulian c cym CALL iim cym klev iim
Definition: ini_bilKP_ave.h:24
real(k8), parameter pi_d
Definition: nrtype.F90:15