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

Annotation of /trunk/filtrez/inifilr.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 313 - (hide annotations)
Mon Dec 10 15:54:30 2018 UTC (5 years, 6 months ago) by guez
Original Path: trunk/filtrez/inifilr.f
File size: 3933 byte(s)
Remove module temps. Move variable itau_dyn from module temps to
module dynetat0_m, where it is defined.

Split module dynetat0_m into dynetat0_m and dynetat0_chosen_m. The
motivation is to create smaller modules. Procedures principal_cshift
and invert_zoomx had to stay in dynetat0_m because of circular
dependency. Now we will be able to move them away. Module variables
which are chosen by the user, not computed, in program ce0l go to
dynetat0_chosen_m: day_ref, annee_ref, clon, clat, grossismx,
grossismy, dzoomx, dzoomy, taux, tauy.

Move variable "pa" from module disvert_m to module
dynetat0_chosen_m. Define "pa" in dynetat0_chosen rather than etat0.

Define day_ref and annee_ref in procedure read_serre rather than
etat0.

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 313 ! Libraries:
43 guez 143 use jumble, only: new_unit
44 guez 169 use nr_util, only: pi, ifirstloc, assert
45 guez 3
46 guez 313 USE dimensions, ONLY: iim, jjm
47     USE dynetat0_m, ONLY: rlatu, rlatv, xprimu
48     USE dynetat0_chosen_m, ONLY: grossismx
49     use inifgn_m, only: inifgn
50     use inifilr_hemisph_m, only: inifilr_hemisph
51    
52 guez 54 ! Local:
53 guez 164
54 guez 132 REAL dlatu(jjm)
55 guez 166 REAL rlamda(2:iim) ! > 0, in descending order
56 guez 168 real eignvl(iim) ! eigenvalues (<= 0) sorted in descending order
57 guez 166 INTEGER j, unit
58 guez 151 REAL colat0 ! > 0
59 guez 169 integer j1 ! index of negative latitude closest to the equator
60 guez 3
61 guez 154 real eignfnu(iim, iim), eignfnv(iim, iim)
62 guez 169 ! eigenvectors of the discrete second derivative with respect to
63     ! longitude, at rlon[uv] longitudes
64 guez 154
65 guez 54 !-----------------------------------------------------------
66 guez 3
67 guez 54 print *, "Call sequence information: inifilr"
68 guez 3
69 guez 154 CALL inifgn(eignvl, eignfnu, eignfnv)
70 guez 3
71 guez 54 ! Calcul de colat0
72 guez 143 forall (j = 1:jjm) dlatu(j) = rlatu(j) - rlatu(j + 1)
73     colat0 = min(0.5, minval(dlatu) / minval(xprimu(:iim)))
74 guez 54 PRINT *, 'colat0 = ', colat0
75 guez 3
76 guez 169 rlamda = iim / pi / colat0 * grossismx / sqrt(- eignvl(2: iim))
77 guez 168 print *, "1 / rlamda(iim) = ", 1. / rlamda(iim)
78 guez 169 ! This is demonstrated in the notes but just to be sure:
79 guez 175 call assert(rlamda(iim) * colat0 >= 1. - 2. * epsilon(0.), &
80 guez 169 "inifilr rlamda(iim) * colat0")
81    
82 guez 166 call new_unit(unit)
83     open(unit, file = "modfrst.csv", status = "replace", action = "write")
84     write(unit, fmt = *) '"rlat (degrees)" modfrst' ! title line
85 guez 3
86 guez 167 j1 = ifirstloc(rlatu <= 0.)
87 guez 3
88 guez 169 call inifilr_hemisph(rlatu(j1 - 1:2:- 1), rlamda, unit, eignfnv, jfiltnu, &
89     matriceun, matrinvn)
90 guez 167 jfiltnu = j1 - jfiltnu
91 guez 166 matriceun = matriceun(:, :, jfiltnu - 1:1:- 1)
92     matrinvn = matrinvn(:, :, jfiltnu - 1:1:- 1)
93 guez 3
94 guez 169 call inifilr_hemisph(- rlatu(j1:jjm), rlamda, unit, eignfnv, jfiltsu, &
95     matriceus, matrinvs)
96 guez 167 jfiltsu = j1 - 1 + jfiltsu
97 guez 3
98 guez 167 j1 = ifirstloc(rlatv <= 0.)
99 guez 3
100 guez 169 call inifilr_hemisph(rlatv(j1 - 1:1:- 1), rlamda, unit, eignfnu, jfiltnv, &
101     matricevn)
102 guez 167 jfiltnv = j1 - jfiltnv
103 guez 166 matricevn = matricevn(:, :, jfiltnv:1:- 1)
104 guez 3
105 guez 169 call inifilr_hemisph(- rlatv(j1:jjm), rlamda, unit, eignfnu, jfiltsv, &
106     matricevs)
107 guez 167 jfiltsv = j1 - 1 + jfiltsv
108 guez 3
109 guez 166 close(unit)
110 guez 151 PRINT *, 'jfiltnu =', jfiltnu
111     PRINT *, 'jfiltsu =', jfiltsu
112     PRINT *, 'jfiltnv =', jfiltnv
113     PRINT *, 'jfiltsv =', jfiltsv
114 guez 3
115 guez 54 END SUBROUTINE inifilr
116 guez 32
117 guez 54 end module inifilr_m

  ViewVC Help
Powered by ViewVC 1.1.21