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

Annotation of /trunk/filtrez/inifilr.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 169 - (hide annotations)
Mon Sep 14 17:13:16 2015 UTC (8 years, 9 months ago) by guez
Original Path: trunk/Sources/filtrez/inifilr.f
File size: 3879 byte(s)
In inifilr_hemisph, colat0 is necessarily >= 1. / rlamda(iim) (see
notes) so we simplify the definition of jfilt. No need to keep modfrst
values at other latitudes than the current one, and we can have one
loop on latitudes instead of two.

Just encapsulated transp into a module.

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

  ViewVC Help
Powered by ViewVC 1.1.21