source: LMDZ5/trunk/libf/dyn3d/pres2lev_mod.F90 @ 1939

Last change on this file since 1939 was 1932, checked in by lguez, 11 years ago

Declaration of variables appearing in array bounds must be specified
before the declaration of the array.

Variables used in FCTTRE.h are declared in YOMCST.h and YOETHF.h so
YOMCST.h and YOETHF.h must be included before FCTTRE.h.

(There were compilation errors with gfortran and debugging options.)

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 2.0 KB
Line 
1! $Id: pres2lev.F 1179 2009-06-11 14:18:47Z jghattas $
2!
3MODULE pres2lev_mod
4
5CONTAINS
6
7!******************************************************
8SUBROUTINE pres2lev(varo,varn,lmo,lmn,po,pn,ni,nj,ok_invertp)
9!
10! interpolation lineaire pour passer
11! a une nouvelle discretisation verticale pour
12! les variables de GCM
13! Francois Forget (01/1995)
14! MOdif remy roca 12/97 pour passer de pres2sig
15! Modif F.Codron 07/08 po en 3D
16!**********************************************************
17
18  IMPLICIT NONE
19
20!   Declarations:
21! ==============
22!
23!  ARGUMENTS
24!  """""""""
25  LOGICAL, INTENT(IN) :: ok_invertp
26  INTEGER, INTENT(IN) :: lmo ! dimensions ancienne couches
27  INTEGER, INTENT(IN) :: lmn ! dimensions nouvelle couches
28 
29  INTEGER, INTENT(IN) :: ni,nj ! nombre de point horizontal
30  REAL, INTENT(IN) :: po(ni*nj,lmo) ! niveau de pression ancienne grille
31  REAL, INTENT(IN) :: pn(ni*nj,lmn) ! niveau de pression nouvelle grille
32
33  REAL, INTENT(IN)  :: varo(ni*nj,lmo) ! var dans l'ancienne grille
34  REAL, INTENT(OUT) :: varn(ni*nj,lmn) ! var dans la nouvelle grille
35
36  REAL :: zvaro(ni*nj,lmo),zpo(ni*nj,lmo)
37
38! Autres variables
39! """"""""""""""""
40  INTEGER ::  ln ,lo, k
41  REAL    :: coef
42
43
44! Inversion de l'ordre des niveaux verticaux
45  IF (ok_invertp) THEN
46    DO lo=1,lmo
47      DO k=1,ni*nj
48        zpo(k,lo)=po(k,lmo+1-lo)
49        zvaro(k,lo)=varo(k,lmo+1-lo)
50      ENDDO
51    ENDDO
52  ELSE
53    DO lo=1,lmo
54      DO k=1,ni*nj
55        zpo(k,lo)=po(k,lo)
56        zvaro(k,lo)=varo(k,lo)
57      ENDDO
58    ENDDO
59  ENDIF
60
61  DO ln=1,lmn
62    DO lo=1,lmo-1
63      DO k=1,ni*nj
64        IF (pn(k,ln) >= zpo(k,1) ) THEN
65          varn(k,ln) = zvaro(k,1)
66        ELSE IF (pn(k,ln) <= zpo(k,lmo)) THEN
67          varn(k,ln) = zvaro(k,lmo)
68        ELSE IF ( pn(k,ln) <= zpo(k,lo) .AND. pn(k,ln) > zpo(k,lo+1) ) THEN
69          coef = (pn(k,ln)-zpo(k,lo)) / (zpo(k,lo+1)-zpo(k,lo))
70          varn(k,ln) = zvaro(k,lo) + coef*(zvaro(k,lo+1)-zvaro(k,lo))
71        ENDIF
72         
73      ENDDO 
74    ENDDO
75  ENDDO               
76
77END SUBROUTINE pres2lev   
78
79END MODULE pres2lev_mod
Note: See TracBrowser for help on using the repository browser.