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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 154 - (show annotations)
Tue Jul 7 17:49:23 2015 UTC (8 years, 10 months ago) by guez
File size: 2148 byte(s)
Removed argument dtphys of physiq. Use it directly from comconst in
physiq instead.

Donwgraded variables eignfnu, eignfnv of module inifgn_m to dummy
arguments of SUBROUTINE inifgn. They were not used elsewhere than in
the calling procedure inifilr. Renamed argument dv of inifgn to eignval_v.

Made alboc and alboc_cd independent of the size of arguments. Now we
can call them only at indices knindex in interfsurf_hq, where we need
them. Fixed a bug in alboc_cd: rmu0 was modified, and the
corresponding actual argument in interfsurf_hq is an intent(in)
argument of interfsurf_hq.

Variables of size knon instead of klon in interfsur_lim and interfsurf_hq.

Removed argument alb_new of interfsurf_hq because it was the same than
alblw. Simplified test on cycle_diurne, following LMDZ.

Moved tests on nbapp_rad from physiq to read_clesphys2. No need for
separate counter itaprad, we can use itap. Define lmt_pas and radpas
from integer input parameters instead of real-type computed values.

1 module invert_zoom_x_m
2
3 implicit none
4
5 INTEGER, PARAMETER:: nmax = 30000
6 DOUBLE PRECISION abs_y
7
8 private abs_y, funcd
9
10 contains
11
12 subroutine invert_zoom_x(xf, xtild, G, xlon, xprim, xuv)
13
14 use coefpoly_m, only: coefpoly, a1, a2, a3
15 USE dimens_m, ONLY: iim
16 use dynetat0_m, only: clon
17 use nr_util, only: pi_d, twopi_d
18 use numer_rec_95, only: hunt, rtsafe
19
20 DOUBLE PRECISION, intent(in):: Xf(0:), xtild(0:), G(0:) ! (0:nmax)
21
22 real, intent(out):: xlon(:), xprim(:) ! (iim)
23
24 DOUBLE PRECISION, intent(in):: xuv
25 ! between - 0.25 and 0.5
26 ! 0. si calcul aux points scalaires
27 ! 0.5 si calcul aux points U
28
29 ! Local:
30 DOUBLE PRECISION Y
31 DOUBLE PRECISION h ! step of the uniform grid
32 integer i, it
33
34 DOUBLE PRECISION xvrai(iim), Gvrai(iim)
35 ! intermediary variables because xlon and xprim are simple precision
36
37 !------------------------------------------------------------------
38
39 it = 0 ! initial guess
40 h = twopi_d / iim
41
42 DO i = 1, iim
43 Y = - pi_d + (i + xuv - 0.75d0) * h
44 ! - pi <= y < pi
45 abs_y = abs(y)
46
47 call hunt(xf, abs_y, it, my_lbound = 0)
48 ! {0 <= it <= nmax - 1}
49
50 ! Calcul de xvrai(i) et Gvrai(i)
51 CALL coefpoly(Xf(it), Xf(it + 1), G(it), G(it + 1), xtild(it), &
52 xtild(it + 1))
53 xvrai(i) = rtsafe(funcd, xtild(it), xtild(it + 1), xacc = 1d-6)
54 Gvrai(i) = a1 + xvrai(i) * (2d0 * a2 + xvrai(i) * 3d0 * a3)
55 if (y < 0d0) xvrai(i) = - xvrai(i)
56 end DO
57
58 DO i = 1, iim -1
59 IF (xvrai(i + 1) < xvrai(i)) THEN
60 print *, 'xvrai(', i + 1, ') < xvrai(', i, ')'
61 STOP 1
62 END IF
63 END DO
64
65 xlon = xvrai + clon
66 xprim = h / Gvrai
67
68 end subroutine invert_zoom_x
69
70 !**********************************************************************
71
72 SUBROUTINE funcd(x, fval, fderiv)
73
74 use coefpoly_m, only: a0, a1, a2, a3
75
76 DOUBLE PRECISION, INTENT(IN):: x
77 DOUBLE PRECISION, INTENT(OUT):: fval, fderiv
78
79 fval = a0 + x * (a1 + x * (a2 + x * a3)) - abs_y
80 fderiv = a1 + x * (2d0 * a2 + x * 3d0 * a3)
81
82 END SUBROUTINE funcd
83
84 end module invert_zoom_x_m

  ViewVC Help
Powered by ViewVC 1.1.21