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

Contents of /trunk/filtrez/inifilr.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 168 - (show annotations)
Wed Sep 9 10:41:47 2015 UTC (8 years, 9 months ago) by guez
Original Path: trunk/Sources/filtrez/inifilr.f
File size: 3797 byte(s)
In order to be able to choose finer resolutions, set large memory
model in compiler options and use dynamic libraries.

Variables rlatd, rlond, cuphy and cvphy of module comgeomphy were
never used. (In LMDZ, they are used only for Orchid.)

There is a bug in PGI Fortran 13.10 that does not accept the
combination of forall, pack and spread in regr_pr_av and
regr_pr_int. In order to circumvent this bug, created the function
gr_dyn_phy.

In program test_inifilr, use a single latitude coordinate for north
and south.

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 (<= 0) sorted in descending order
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 print *, "1 / rlamda(iim) = ", 1. / rlamda(iim)
81 call new_unit(unit)
82 open(unit, file = "modfrst.csv", status = "replace", action = "write")
83 write(unit, fmt = *) '"rlat (degrees)" modfrst' ! title line
84
85 j1 = ifirstloc(rlatu <= 0.)
86
87 call inifilr_hemisph(rlatu(j1 - 1:2:- 1), colat0, rlamda, unit, eignfnv, &
88 jfiltnu, matriceun, matrinvn)
89 jfiltnu = j1 - jfiltnu
90 matriceun = matriceun(:, :, jfiltnu - 1:1:- 1)
91 matrinvn = matrinvn(:, :, jfiltnu - 1:1:- 1)
92
93 call inifilr_hemisph(- rlatu(j1:jjm), colat0, rlamda, unit, eignfnv, &
94 jfiltsu, matriceus, matrinvs)
95 jfiltsu = j1 - 1 + jfiltsu
96
97 j1 = ifirstloc(rlatv <= 0.)
98
99 call inifilr_hemisph(rlatv(j1 - 1:1:- 1), colat0, rlamda, unit, eignfnu, &
100 jfiltnv, matricevn)
101 jfiltnv = j1 - jfiltnv
102 matricevn = matricevn(:, :, jfiltnv:1:- 1)
103
104 call inifilr_hemisph(- rlatv(j1:jjm), colat0, rlamda, unit, eignfnu, &
105 jfiltsv, matricevs)
106 jfiltsv = j1 - 1 + jfiltsv
107
108 close(unit)
109 PRINT *, 'jfiltnu =', jfiltnu
110 PRINT *, 'jfiltsu =', jfiltsu
111 PRINT *, 'jfiltnv =', jfiltnv
112 PRINT *, 'jfiltsv =', jfiltsv
113
114 END SUBROUTINE inifilr
115
116 end module inifilr_m

  ViewVC Help
Powered by ViewVC 1.1.21