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

Annotation of /trunk/Sources/filtrez/filtreg_hemisph.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 169 - (hide annotations)
Mon Sep 14 17:13:16 2015 UTC (8 years, 8 months ago) by guez
File size: 763 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 136 module filtreg_hemisph_m
2    
3     implicit none
4    
5     contains
6    
7 guez 141 subroutine filtreg_hemisph(champ, sdd, matri)
8 guez 136
9     USE dimens_m, ONLY: iim
10    
11     REAL, intent(inout):: champ(:, :, :) ! (iim + 1, :, :)
12 guez 168 REAL, intent(in):: sdd(:) ! (iim) xprim[uv]^{\pm 1/2}
13 guez 136
14 guez 168 real, intent(in), dimension(:, :, :):: matri ! (iim, iim, :)
15 guez 169 ! filtering matrix, last dimension is latitude
16 guez 168
17 guez 136 ! Local:
18     integer l, j
19    
20     !-----------------------------------------------------------------
21    
22 guez 142 forall (j = 1:size(champ, 2), l = 1:size(champ, 3))
23     champ(:iim, j, l) = champ(:iim, j, l) &
24     + matmul(matri(:, :, j), champ(:iim, j, l) * sdd) / sdd
25     END forall
26 guez 136
27 guez 142 champ(iim + 1, :, :) = champ(1, :, :)
28    
29 guez 136 end subroutine filtreg_hemisph
30    
31     end module filtreg_hemisph_m

  ViewVC Help
Powered by ViewVC 1.1.21