source: trunk/LMDZ.GENERIC/libf/phystd/pindex.F90 @ 3523

Last change on this file since 3523 was 3233, checked in by gmilcareck, 9 months ago

Add the possibility to use a fixed vertical molar fraction profile for the
collision-induced absorption like the correlated-k.

File size: 2.2 KB
RevLine 
[3233]1module pindex_mod
2
3implicit none
4
5contains
6
7subroutine pindex(pti,vti,pref,nvarlev,nlev,nv)
8
9INTEGER, INTENT(IN) :: nvarlev
10INTEGER, INTENT(IN) :: nlev
11REAL, INTENT(IN) :: pti(nvarlev)
12REAL, INTENT(IN) :: vti(nvarlev)
13REAL, INTENT(IN) :: pref(nlev)
14
15REAL, INTENT(OUT) :: nv(nlev)
16
17REAL :: a,b
18INTEGER :: n,l
19
20!==================================================================
21!     
22!     Purpose
23!     -------
24!     Interpolate only data to the given P values.
25!
26!     Inputs
27!     ------
28!     pti                - The pressure to interpolate to
29!     vti                - The variable to interpolate to
30!     pref               - The tabulated pressure grid
31!     nvarlev            - Number of levels of pti
32!     nlev               - Number of levels of pref
33!   
34!     Outputs
35!     -------
36!     nv                 - Interpolated variable
37!
38!     Local variables
39!     -------
40!     n,l                - integer
41!     m                  - coefficient for interpolation
42!     k                  - coefficient for interpolation
43!                       
44!==================================================================
45
46DO n=1,nlev
47  IF(pti(1) .LT. pti(nvarlev)) THEN
48    IF(pref(n).LE.pti(1)) THEN
49      nv(n)=vti(1)
50    ELSEIF(pref(n).GE.pti(nvarlev)) THEN
51      nv(n)=vti(nvarlev)
52    ELSE
53      DO l=2,nvarlev
54        IF(pti(l-1).LE.pref(n) .AND. pref(n).LE.pti(l)) THEN
55          IF(abs(vti(l-1)-vti(l)).LT.1e-6) THEN
56            nv(n) = (vti(l-1) + vti(l))/2.
57          ELSE
58            a = (pti(l)-pti(l-1))/(vti(l)-vti(l-1))
59            b = (vti(l)*pti(l-1)-vti(l-1)*pti(l))/(vti(l)-vti(l-1))
60            nv(n) = (pref(n)-b)/a
61          ENDIF
62        ENDIF
63      ENDDO
64    END IF
65  ELSE
66    IF(pref(n).LE.pti(nvarlev)) THEN
67      nv(n)=vti(nvarlev)
68    ELSEIF(pref(n).GE.pti(1)) THEN
69      nv(n)=vti(1)
70    ELSE
71      DO l=2,nvarlev
72        IF(pti(l-1).GE.pref(n) .AND. pref(n).GE.pti(l)) THEN
73          IF(abs(vti(l-1)-vti(l)).LT.1e-6) THEN
74            nv(n) = (vti(l-1) + vti(l))/2.
75          ELSE
76            a = (pti(l)-pti(l-1))/(vti(l)-vti(l-1))
77            b = (vti(l)*pti(l-1)-vti(l-1)*pti(l))/(vti(l)-vti(l-1))
78            nv(n) = (pref(n)-b)/a
79          ENDIF
80        ENDIF
81      ENDDO
82    END IF
83  END IF
84ENDDO
85
86end subroutine pindex
87end module pindex_mod
Note: See TracBrowser for help on using the repository browser.