1 |
! |
module laplacien_m |
|
! $Header: /home/cvsroot/LMDZ4/libf/dyn3d/laplacien.F,v 1.1.1.1 2004/05/19 12:53:06 lmdzadmin Exp $ |
|
|
! |
|
|
SUBROUTINE laplacien ( klevel, teta, divgra ) |
|
|
c |
|
|
c P. Le Van |
|
|
c |
|
|
c ************************************************************ |
|
|
c .... calcul de (div( grad )) de teta ..... |
|
|
c ************************************************************ |
|
|
c klevel et teta sont des arguments d'entree pour le s-prog |
|
|
c divgra est un argument de sortie pour le s-prog |
|
|
c |
|
|
use grad_m, only: grad |
|
|
use dimens_m |
|
|
use paramet_m |
|
|
use comgeom |
|
|
use filtreg_m, only: filtreg |
|
|
use divergf_m, only: divergf |
|
|
|
|
|
IMPLICIT NONE |
|
|
c |
|
|
|
|
|
c |
|
|
c ......... variables en arguments .............. |
|
|
c |
|
|
INTEGER, intent(in):: klevel |
|
|
REAL teta( ip1jmp1,klevel ), divgra( ip1jmp1,klevel ) |
|
|
c |
|
|
c ............ variables locales .............. |
|
|
c |
|
|
REAL ghy(ip1jm,llm), ghx(ip1jmp1,llm) |
|
|
c ....................................................... |
|
|
|
|
|
|
|
|
c |
|
|
CALL SCOPY ( ip1jmp1 * klevel, teta, 1, divgra, 1 ) |
|
|
|
|
|
CALL filtreg( divgra, jjp1, klevel, 2, 1, .TRUE., 1 ) |
|
|
CALL grad ( klevel,divgra, ghx , ghy ) |
|
|
CALL divergf ( klevel, ghx , ghy , divgra ) |
|
2 |
|
|
3 |
RETURN |
IMPLICIT NONE |
4 |
END |
|
5 |
|
contains |
6 |
|
|
7 |
|
SUBROUTINE laplacien(klevel, teta) |
8 |
|
|
9 |
|
! From LMDZ4/libf/dyn3d/laplacien.F, version 1.1.1.1 2004/05/19 12:53:06 |
10 |
|
! P. Le Van |
11 |
|
! Calcul de div(grad) de teta. |
12 |
|
|
13 |
|
use grad_m, only: grad |
14 |
|
use filtreg_m, only: filtreg |
15 |
|
use divergf_m, only: divergf |
16 |
|
USE dimens_m, ONLY: llm |
17 |
|
USE paramet_m, ONLY: ip1jm, ip1jmp1, jjp1 |
18 |
|
|
19 |
|
INTEGER, intent(in):: klevel |
20 |
|
REAL, intent(inout):: teta(ip1jmp1, klevel) |
21 |
|
|
22 |
|
! Variables locales: |
23 |
|
REAL ghy(ip1jm, llm), ghx(ip1jmp1, llm) |
24 |
|
|
25 |
|
!----------------------------------------------------------------- |
26 |
|
|
27 |
|
CALL filtreg(teta, jjp1, klevel, 2, 1, .TRUE.) |
28 |
|
CALL grad(klevel, teta, ghx, ghy) |
29 |
|
CALL divergf(klevel, ghx, ghy, teta) |
30 |
|
|
31 |
|
END SUBROUTINE laplacien |
32 |
|
|
33 |
|
end module laplacien_m |