/[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 168 - (hide annotations)
Wed Sep 9 10:41:47 2015 UTC (8 years, 8 months ago) by guez
File size: 735 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 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     ! filtering matrix
16    
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