source: branches/iLoveclim/SOURCES/Old-sources/moy_mxmy_shift.f90 @ 77

Last change on this file since 77 was 77, checked in by dumas, 8 years ago

Merge branche iLOVECLIM sur rev 76

File size: 2.2 KB
Line 
1!> \file moy_mxmy_shift.f90
2!!  Fait la moyenne ponderee d'un tableau X2D sur les mailles staggered
3!<
4
5!> SUBROUTINE: moy_mxmy()
6!! \author ...
7!! \date ...
8!! @note Cette routine permet de faire  la moyenne ponderee d'un tableau X2D sur les mailles staggered
9!! \param n1   [in]    dimension des tableaux
10!! \param n2   [in]    dimension des tableaux
11!! \param X2D  [in]    tableau sur les noeuds majeurs
12!! \param X_mx [out]  tableau sur les noeuds mineurs x
13!! \param X_my [out]  tableau sur les noeuds mineurs y
14!! \return X_mx, X_my
15!<
16subroutine moy_mxmy(n1,n2,X2D,X_mx,X_my)
17! fait la moyenne ponderee d'un tableau X2D sur les mailles staggered
18
19implicit none
20integer, intent(in) :: n1,n2   !< dimension des tableaux
21real, dimension(n1,n2),intent(in) :: X2D   !< tableau sur les noeuds majeurs
22real, dimension(n1,n2),intent(out) :: X_mx  !< tableau sur les noeuds mineurs x
23real, dimension(n1,n2),intent(out) :: X_my  !< tableau sur les noeuds mineurs y
24
25! tableaux de travail
26
27real, dimension(n1,n2) :: X_sud
28real, dimension(n1,n2) :: X_nord
29real, dimension(n1,n2) :: X_west
30real, dimension(n1,n2) :: X_est
31real, dimension(n1,n2) :: X_SW
32real, dimension(n1,n2) :: X_SE
33real, dimension(n1,n2) :: X_NW
34real, dimension(n1,n2) :: X_NE
35
36
37
38integer :: i,j
39integer :: i_moins1,j_moins1,i_plus1,j_plus1
40
41! l'orientation correspond a la direction vers laquelle on se deplace
42!---------------------------------------------------------------------
43X_west = eoshift(X2D(:,:),shift=-1,boundary=0.,dim=1)  ! X_west(i,j) = X(i-1,j)
44X_est  = eoshift(X2D(:,:),shift=+1,boundary=0.,dim=1)  ! X_est(i,j)  = X(i+1,j)
45X_sud  = eoshift(X2D(:,:),shift=-1,boundary=0.,dim=2)  ! X_sud(i,j)  = X(i,j-1)
46X_nord = eoshift(X2D(:,:),shift=+1,boundary=0.,dim=2)  ! X_nord(i,j) = X(i,j+1)
47
48X_se   = eoshift( X_sud(:,:),shift=+1,boundary=0.,dim=1)  ! X_se(i,j) = X(i+1,j-1)
49X_sw   = eoshift( X_sud(:,:),shift=-1,boundary=0.,dim=1)  ! X_sw(i,j) = X(i-1,j-1)
50X_ne   = eoshift(X_nord(:,:),shift=+1,boundary=0.,dim=1)  ! X_ne(i,j) = X(i+1,j+1)
51X_nw   = eoshift(X_nord(:,:),shift=-1,boundary=0.,dim=1)  ! X_nw(i,j) = X(i-1,j+1)
52
53X_mx = 0.25*(X2D+X_west)+0.125*((X_nw+X_nord)+(X_sw+X_sud))
54X_my = 0.25*(X2D+X_sud)+0.125*((X_est+X_se)+(X_sw+X_west))
55
56
57return
58end subroutine moy_mxmy
59
60
Note: See TracBrowser for help on using the repository browser.