/[lmdze]/trunk/Sources/filtrez/inifilr.f
ViewVC logotype

Contents of /trunk/Sources/filtrez/inifilr.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 167 - (show annotations)
Mon Aug 24 16:30:33 2015 UTC (8 years, 10 months ago) by guez
File size: 3786 byte(s)
Added program test_inifilr.

Encapsulated ppm3d into a module and added implicit none. Removed
unused argument dum.

Encountered a problem in procedure invert_zoom_x. With grossismx=2.9,
DZOOMX=0.3, taux=5, for xuv = -0.25, for i = 1, rtsafe fails because
fval is about 1e-16 instead of 0 at xval = pi. So distinguished the
cases abs_y = 0 or pi. Needed then to add argument beta to
invert_zoom_x.

Moved the output of eignvalues of differentiation matrix from inifilr
to inifgn, where they are computed.

Simpler definition of j1 in inifilr.

1 module inifilr_m
2
3 IMPLICIT NONE
4
5 ! North:
6
7 INTEGER jfiltnu, jfiltnv
8 ! index of the last scalar line filtered in northern hemisphere
9
10 real, pointer:: matriceun(:, :, :) ! (iim, iim, jfiltnu - 1)
11 ! matrice filtre pour les champs situes sur la grille scalaire
12
13 real, pointer:: matrinvn(:, :, :) ! (iim, iim, jfiltnu - 1)
14 ! matrice filtre pour les champs situes sur la grille scalaire, pour
15 ! le filtre inverse
16
17 real, pointer:: matricevn(:, :, :) ! (iim, iim, jfiltnv)
18 ! matrice filtre pour les champs situes sur la grille de V ou de Z
19
20 ! South:
21
22 integer jfiltsu, jfiltsv
23 ! index of the first line filtered in southern hemisphere
24
25 real, pointer:: matriceus(:, :, :) ! (iim, iim, jjm - jfiltsu + 1)
26 ! matrice filtre pour les champs situes sur la grille scalaire
27
28 real, pointer:: matrinvs(:, :, :) ! (iim, iim, jjm - jfiltsu + 1)
29 ! matrice filtre pour les champs situes sur la grille scalaire, pour
30 ! le filtre inverse
31
32 real, pointer:: matricevs(:, :, :) ! (iim, iim, jjm - jfiltsv + 1)
33 ! matrice filtre pour les champs situes sur la grille de V ou de Z
34
35 contains
36
37 SUBROUTINE inifilr
38
39 ! From filtrez/inifilr.F, version 1.1.1.1, 2004/05/19 12:53:09
40 ! H. Upadhyaya, O. Sharma
41
42 ! This procedure computes the filtering coefficients for scalar
43 ! lines and meridional wind v lines. The modes are filtered from
44 ! modfrst to iim. We filter all those latitude lines where coefil
45 ! < 1. No filtering at poles. colat0 is to be used when alpha
46 ! (stretching coefficient) is set equal to zero for the regular
47 ! grid case.
48
49 USE dimens_m, ONLY : iim, jjm
50 USE dynetat0_m, ONLY : rlatu, rlatv, xprimu, grossismx
51 use inifgn_m, only: inifgn
52 use inifilr_hemisph_m, only: inifilr_hemisph
53 use jumble, only: new_unit
54 use nr_util, only: pi, ifirstloc
55
56 ! Local:
57
58 REAL dlatu(jjm)
59 REAL rlamda(2:iim) ! > 0, in descending order
60 real eignvl(iim) ! eigenvalues sorted in descending order (<= 0)
61 INTEGER j, unit
62 REAL colat0 ! > 0
63 integer j1 ! index of smallest positive latitude
64
65 real eignfnu(iim, iim), eignfnv(iim, iim)
66 ! eigenvectors of the discrete second derivative with respect to longitude
67
68 !-----------------------------------------------------------
69
70 print *, "Call sequence information: inifilr"
71
72 CALL inifgn(eignvl, eignfnu, eignfnv)
73
74 ! Calcul de colat0
75 forall (j = 1:jjm) dlatu(j) = rlatu(j) - rlatu(j + 1)
76 colat0 = min(0.5, minval(dlatu) / minval(xprimu(:iim)))
77 PRINT *, 'colat0 = ', colat0
78
79 rlamda = iim / (pi * colat0 / grossismx) / sqrt(- eignvl(2: iim))
80 call new_unit(unit)
81 open(unit, file = "modfrst.csv", status = "replace", action = "write")
82 write(unit, fmt = *) '"rlat (degrees)" modfrst' ! title line
83
84 ! D\'etermination de jfilt[ns][uv] :
85
86 j1 = ifirstloc(rlatu <= 0.)
87
88 call inifilr_hemisph(rlatu(j1 - 1:2:- 1), colat0, rlamda, unit, eignfnv, &
89 jfiltnu, matriceun, matrinvn)
90 jfiltnu = j1 - jfiltnu
91 matriceun = matriceun(:, :, jfiltnu - 1:1:- 1)
92 matrinvn = matrinvn(:, :, jfiltnu - 1:1:- 1)
93
94 call inifilr_hemisph(- rlatu(j1:jjm), colat0, rlamda, unit, eignfnv, &
95 jfiltsu, matriceus, matrinvs)
96 jfiltsu = j1 - 1 + jfiltsu
97
98 j1 = ifirstloc(rlatv <= 0.)
99
100 call inifilr_hemisph(rlatv(j1 - 1:1:- 1), colat0, rlamda, unit, eignfnu, &
101 jfiltnv, matricevn)
102 jfiltnv = j1 - jfiltnv
103 matricevn = matricevn(:, :, jfiltnv:1:- 1)
104
105 call inifilr_hemisph(- rlatv(j1:jjm), colat0, rlamda, unit, eignfnu, &
106 jfiltsv, matricevs)
107 jfiltsv = j1 - 1 + jfiltsv
108
109 close(unit)
110 PRINT *, 'jfiltnu =', jfiltnu
111 PRINT *, 'jfiltsu =', jfiltsu
112 PRINT *, 'jfiltnv =', jfiltnv
113 PRINT *, 'jfiltsv =', jfiltsv
114
115 END SUBROUTINE inifilr
116
117 end module inifilr_m

  ViewVC Help
Powered by ViewVC 1.1.21