/[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 139 - (hide annotations)
Tue May 26 17:46:03 2015 UTC (9 years ago) by guez
File size: 2239 byte(s)
dynetat0 read rlonu, rlatu, rlonv, rlatv, cu_2d, cv_2d, aire_2d from
"start.nc" and then these variables were overwritten by
inigeom. Corrected this. Now, inigeom does not compute rlonu, rlatu,
rlonv and rlatv. Moreover, cu_2d, cv_2d, aire_2d are not written to
"restart.nc". Since xprimu, xprimv, xprimm025, xprimp025, rlatu1,
rlatu2, yprimu1, yprimu2 are computed at the same time as rlonu,
rlatu, rlonv, rlatv, and since it would not be convenient to separate
those computations, we decide to write xprimu, xprimv, xprimm025,
xprimp025, rlatu1, rlatu2, yprimu1, yprimu2 into "restart.nc", read
them from "start.nc" and not compute them in inigeom. So, in summary,
"start.nc" contains all the coordinates and their derivatives, and
inigeom only computes the 2D-variables.

Technical details:

Moved variables rlatu, rlonv, rlonu, rlatv, xprimu, xprimv from module
comgeom to module dynetat0_m. Upgraded local variables rlatu1,
yprimu1, rlatu2, yprimu2, xprimm025, xprimp025 of procedure inigeom to
variables of module dynetat0_m.

Removed unused local variable yprimu of procedure inigeom and
corresponding argument yyprimu of fyhyp.

Moved variables clat, clon, grossismx, grossismy, dzoomx, dzoomy,
taux, tauy from module serre to module dynetat0_m (since they are read
from "start.nc"). The default values are now defined in read_serre
instead of in the declarations. Changed name of module serre to
read_serre_m, no more module variable here.

The calls to fxhyp and fyhyp are moved from inigeom to etat0.

Side effects in programs other than gcm: etat0 and read_serre write
variables of module dynetat0; the programs test_fxyp and
test_inter_barxy need more source files.

Removed unused arguments len and nd of cv3_tracer. Removed unused
argument PPSOL of LWU.

Bug fix in test_inter_barxy: forgotten call to read_serre.

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 121
16     DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), Xprimt(0:) ! (0:2 * nmax)
17 guez 124 real, intent(out):: xlon(:), xprimm(:) ! (iim)
18 guez 121
19     DOUBLE PRECISION, intent(in):: xuv
20     ! 0. si calcul aux points scalaires
21     ! 0.5 si calcul aux points U
22    
23     ! Local:
24 guez 126 DOUBLE PRECISION xo1, Xfi, a0, a1, a2, a3, Xf1, Xprimin
25     integer i, it, iter
26     DOUBLE PRECISION, parameter:: my_eps = 1d-6
27 guez 121
28 guez 126 DOUBLE PRECISION xxprim(iim), xvrai(iim)
29     ! intermediary variables because xlon and xprimm are simple precision
30    
31 guez 121 !------------------------------------------------------------------
32    
33 guez 126 DO i = 1, iim
34     Xfi = - pi_d + (i + xuv - 0.75d0) * twopi_d / iim
35 guez 121
36     it = 2 * nmax
37     do while (xfi < xf(it) .and. it >= 1)
38     it = it - 1
39     end do
40    
41 guez 126 ! Calcul de Xf(xvrai(i))
42 guez 121
43 guez 126 xvrai(i) = xtild(it)
44 guez 121
45     IF (it == 2 * nmax) THEN
46     it = 2 * nmax -1
47     END IF
48    
49     CALL coefpoly(Xf(it), Xf(it + 1), Xprimt(it), Xprimt(it + 1), &
50     xtild(it), xtild(it + 1), a0, a1, a2, a3)
51     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