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

Annotation of /trunk/filtrez/inifilr.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 166 - (hide annotations)
Wed Jul 29 14:32:55 2015 UTC (8 years, 10 months ago) by guez
Original Path: trunk/Sources/filtrez/inifilr.f
File size: 3968 byte(s)
Split ppm3d.f into files containing a single procedure.

Factorized computations of filtering matrices into a procedure
inifilr_hemisph. Had then to change the matrices from allocatable to
pointer and from customized lower bound to lower bound 1. The change
in lower bounds does not matter because the matrices are only used as
a whole as actual arguments.

Also, in infilr, instead of finding jfilt[ns][uv] from approximately
jjm /2, start at index j1 that corresponds to the equator. This is not
the same if there is a zoom in latitude.

Also, the test (rlamda(modfrst[ns][uv](j)) * cos(rlat[uv](j)) < 1) in
the loops on filtered latitudes is not useful now that we start from
j1: it is necessarily true. See notes.

Just encapsulated lwvn into a module and removed unused argument ktraer.

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

  ViewVC Help
Powered by ViewVC 1.1.21