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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 136 - (show annotations)
Thu Apr 30 18:35:49 2015 UTC (9 years ago) by guez
File size: 2247 byte(s)
Clarified the logic in filtreg by creating a procedure
filtreg_hemisph. It was terrible with a loop on hemispheres and tests
on hemisphere inside the loop, plus maddening indirections on latitude
bounds, plus repeated code. Went from 126 lines to much clearer 74 +
32 = 106 lines.

In module inifilr_m, finally made the arrays matrice[uv][ns],
matrinv[ns] dynamic (following LMDZ). Changed the lower bound of
matriceun and matrinvn in the 3rd dimension: 2 instead of 1, the index
1 was not defined (nor used).

In module inifilr_m, changed the bounds of matriceus and matrinvs in
the 3rd dimension: jfiltsu:jjm instead of 1:jjm - jfiltsu + 1. Changed
the bounds of matricevs in the 3rd dimension: jfiltsv:jjm instead of
1:jjm - jfiltsv + 1. It is a little simpler and clearer this way in
procedure inifilr.


1 module filtreg_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE filtreg(champ, direct, intensive)
8
9 ! From filtrez/filtreg.F, version 1.1.1.1, 2004/05/19 12:53:09
10 ! Author: P. Le Van
11 ! Objet : filtre matriciel longitudinal, avec les matrices pr\'ecalcul\'ees
12 ! pour l'op\'erateur filtre.
13
14 USE coefils, ONLY: sddu, sddv, unsddu, unsddv
15 USE dimens_m, ONLY: iim, jjm
16 use filtreg_hemisph_m, only: filtreg_hemisph
17 use inifilr_m, only: jfiltnu, jfiltnv, jfiltsu, jfiltsv, matriceun, &
18 matriceus, matricevn, matricevs, matrinvn, matrinvs
19 use nr_util, only: assert
20
21 REAL, intent(inout):: champ(:, :, :) ! (iim + 1, nlat, :)
22 ! en entr\'ee : champ \`a filtrer, en sortie : champ filtr\'e
23
24 logical, intent(in):: direct ! filtre direct ou inverse
25
26 logical, intent(in):: intensive
27 ! champ intensif ou extensif (pond\'er\'e par les aires)
28
29 ! Local:
30 INTEGER nlat ! nombre de latitudes \`a filtrer
31 REAL sdd1(iim), sdd2(iim)
32
33 !-----------------------------------------------------------
34
35 call assert(size(champ, 1) == iim + 1, "filtreg iim + 1")
36 nlat = size(champ, 2)
37 call assert(nlat == jjm .or. nlat == jjm + 1, "filtreg nlat")
38
39 if (nlat == jjm + 1) then
40 IF (intensive) THEN
41 sdd1 = sddv
42 sdd2 = unsddv
43 ELSE
44 sdd1 = unsddv
45 sdd2 = sddv
46 END IF
47 if (direct) then
48 call filtreg_hemisph(champ(:, 2:jfiltnu, :), sdd1, sdd2, matriceun)
49 call filtreg_hemisph(champ(:, jfiltsu:jjm, :), sdd1, sdd2, matriceus)
50 else
51 call filtreg_hemisph(champ(:, 2:jfiltnu, :), sdd1, sdd2, - matrinvn)
52 call filtreg_hemisph(champ(:, jfiltsu:jjm, :), sdd1, sdd2, - matrinvs)
53 end if
54 else
55 IF (intensive) THEN
56 sdd1 = sddu
57 sdd2 = unsddu
58 ELSE
59 sdd1 = unsddu
60 sdd2 = sddu
61 END IF
62 if (direct) then
63 call filtreg_hemisph(champ(:, :jfiltnv, :), sdd1, sdd2, matricevn)
64 call filtreg_hemisph(champ(:, jfiltsv:jjm, :), sdd1, sdd2, matricevs)
65 else
66 PRINT *, 'filtreg: inverse filter on scalar grid only'
67 STOP 1
68 END IF
69 end if
70
71 END SUBROUTINE filtreg
72
73 end module filtreg_m

  ViewVC Help
Powered by ViewVC 1.1.21