source: LMDZ5/branches/LMDZ5_SPLA/libf/phylmd/plevel.F90 @ 5425

Last change on this file since 5425 was 1992, checked in by lguez, 11 years ago

Converted to free source form files in libf/phylmd which were still in
fixed source form. The conversion was done using the polish mode of
the NAG Fortran Compiler.

In addition to converting to free source form, the processing of the
files also:

-- indented the code (including comments);

-- set Fortran keywords to uppercase, and set all other identifiers
to lower case;

-- added qualifiers to end statements (for example "end subroutine
conflx", instead of "end");

-- changed the terminating statements of all DO loops so that each
loop ends with an ENDDO statement (instead of a labeled continue).

-- replaced #include by include.

  • 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
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 3.3 KB
Line 
1
2! $Header$
3
4! ================================================================
5! ================================================================
6SUBROUTINE plevel(ilon, ilev, lnew, pgcm, pres, qgcm, qpres)
7  ! ================================================================
8  ! ================================================================
9  USE netcdf
10  USE dimphy
11  IMPLICIT NONE
12
13  ! ym#include "dimensions.h"
14  ! y#include "dimphy.h"
15
16  ! ================================================================
17
18  ! Interpoler des champs 3-D u, v et g du modele a un niveau de
19  ! pression donnee (pres)
20
21  ! INPUT:  ilon ----- nombre de points
22  ! ilev ----- nombre de couches
23  ! lnew ----- true si on doit reinitialiser les poids
24  ! pgcm ----- pressions modeles
25  ! pres ----- pression vers laquelle on interpolle
26  ! Qgcm ----- champ GCM
27  ! Qpres ---- champ interpolle au niveau pres
28
29  ! ================================================================
30
31  ! arguments :
32  ! -----------
33
34  INTEGER ilon, ilev
35  LOGICAL lnew
36
37  REAL pgcm(ilon, ilev)
38  REAL qgcm(ilon, ilev)
39  REAL pres
40  REAL qpres(ilon)
41
42  ! local :
43  ! -------
44
45  ! ym      INTEGER lt(klon), lb(klon)
46  ! ym      REAL ptop, pbot, aist(klon), aisb(klon)
47
48  ! ym      save lt,lb,ptop,pbot,aist,aisb
49  INTEGER, ALLOCATABLE, SAVE, DIMENSION (:) :: lt, lb
50  REAL, ALLOCATABLE, SAVE, DIMENSION (:) :: aist, aisb
51  !$OMP THREADPRIVATE(lt,lb,aist,aisb)
52  REAL, SAVE :: ptop, pbot
53  !$OMP THREADPRIVATE(ptop, pbot)
54  LOGICAL, SAVE :: first = .TRUE.
55  !$OMP THREADPRIVATE(first)
56  INTEGER i, k
57
58  REAL missing_val
59
60  missing_val = nf90_fill_real
61
62  IF (first) THEN
63    ALLOCATE (lt(klon), lb(klon), aist(klon), aisb(klon))
64    first = .FALSE.
65  END IF
66
67  ! =====================================================================
68  IF (lnew) THEN
69    ! on r�nitialise les r�ndicages et les poids
70    ! =====================================================================
71
72
73    ! Chercher les 2 couches les plus proches du niveau a obtenir
74
75    ! Eventuellement, faire l'extrapolation a partir des deux couches
76    ! les plus basses ou les deux couches les plus hautes:
77    DO i = 1, klon
78      IF (abs(pres-pgcm(i,ilev))<abs(pres-pgcm(i,1))) THEN
79        lt(i) = ilev ! 2
80        lb(i) = ilev - 1 ! 1
81      ELSE
82        lt(i) = 2
83        lb(i) = 1
84      END IF
85    END DO
86    DO k = 1, ilev - 1
87      DO i = 1, klon
88        pbot = pgcm(i, k)
89        ptop = pgcm(i, k+1)
90        IF (ptop<=pres .AND. pbot>=pres) THEN
91          lt(i) = k + 1
92          lb(i) = k
93        END IF
94      END DO
95    END DO
96
97    ! Interpolation lineaire:
98
99    DO i = 1, klon
100      ! interpolation en logarithme de pression:
101
102      ! ...   Modif . P. Le Van    ( 20/01/98) ....
103      ! Modif Fr��ic Hourdin (3/01/02)
104
105      aist(i) = log(pgcm(i,lb(i))/pres)/log(pgcm(i,lb(i))/pgcm(i,lt(i)))
106      aisb(i) = log(pres/pgcm(i,lt(i)))/log(pgcm(i,lb(i))/pgcm(i,lt(i)))
107    END DO
108
109
110  END IF ! lnew
111
112  ! ======================================================================
113  ! inteprollation
114  ! ======================================================================
115
116  DO i = 1, klon
117    qpres(i) = qgcm(i, lb(i))*aisb(i) + qgcm(i, lt(i))*aist(i)
118  END DO
119
120  ! Je mets les vents a zero quand je rencontre une montagne
121  DO i = 1, klon
122    IF (pgcm(i,1)<pres) THEN
123      qpres(i) = missing_val
124    END IF
125  END DO
126
127
128  RETURN
129END SUBROUTINE plevel
Note: See TracBrowser for help on using the repository browser.