source: trunk/LMDZ.TITAN/libf/phytitan/tpindex.F @ 3533

Last change on this file since 3533 was 3090, checked in by slebonnois, 15 months ago

BdeBatz? : Cleans microphysics and makes few corrections for physics

File size: 4.4 KB
RevLine 
[1648]1      subroutine tpindex(pw,tw,pref,tref,LCOEF,MT,MP)
[135]2
3!==================================================================
4!     
5!     Purpose
6!     -------
[1648]7!     Interpolate K-coefficients to the given P,T.
[135]8!
9!     Notes
10!     -----
11!     The interpolation is the usual one in two dimensions given
12!     in "Numerical Recipes", where the "X" are P, the "Y" are
13!     T, and the F(X,Y) are the CO2 K-coefficients.
14!
15!     The interpolating box is:
16!
17!           (PL,TU)                        (PR,TU)
18!
19!                          (TW,PW)
20!
21!           
22!           (PL,TL)                        (PR,TL)
23!
24!      PL  - Pressure left
25!      PR  - Pressure right
26!      TL  - Temperature lower
27!      TU  - Temperature upper
28!      PW  - Pressure wanted
29!      TW  - Temperature wanted
30!
31!     Inputs
32!     ------
33!     PW                 - The pressure to interpolate to
34!     TW                 - The temperature to interpolate to
35!     Pref(NP)           - The pressure grid array
36!     Tref(NT)           - The temperature grid array
37!   
38!     Outputs
39!     -------
40!     TI                 - Interpolation term (pressure)
41!     UI                 - Interpolation term (temperature)
42!     MT                 - Temperature index (bottom left temperature)
43!                          of bounding box
44!     MP                 - Pressure index (bottom left pressure)
45!                          of bounding box
46!
47!     Authors
48!     -------
49!     Adapted from the NASA Ames code by R. Wordsworth (2009)
[3090]50!     Last update : B. de Batz de Trenquelléon (2023)
[135]51!     
52!==================================================================
53
54      use radinc_h
55
56      implicit none
57
[3083]58      ! INPUTS :
59      !=========
60      REAL,INTENT(IN) :: Tref(L_NTREF)
61      REAL,INTENT(IN) :: pref(L_PINT)
62      REAL,INTENT(IN) :: PW, TW
[135]63
[3083]64      ! OUTPUTS :
65      !==========
66      INTEGER,INTENT(OUT) :: MT, MP
67      REAL,INTENT(OUT) :: LCOEF(4)
[135]68
[3083]69      ! LOCAL :
70      !========
71      INTEGER :: N
72      REAL :: PWL ! local value for log10(PW)
73      REAL :: TWL ! local value for TW
74      REAL :: T ! linear interpolation weight along log(pressure)
75      REAL :: U ! linear interpolation weight along temperature
[135]76
[3083]77
78!     1. For the interpolation along temperature
79!        Get the upper and lower temperature grid indicies that bound the
80!        requested temperature. If the sought temperature is below
81!        the T-grid minimum, assume that lower value. If the sought temperature
82!        is above the T-grid maximum, assume that maximum value.
83!        TW : temperature to interpolate to
84!        TREF : reference temperature array
85!        MT : index of TREF for (lower) bounding temperature
86!        U : weight for interpolation in temperature
87
88      TWL = TW
89      IF(TWL.LT.TREF(1)) THEN
90
91        ! Below lowest tabulated temperature
[135]92        MT = 1
[3083]93        TWL=TREF(1)*1.00
94      ELSEIF(TWL.GE.TREF(L_NTREF)) THEN
95       ! Above highest tabulated temperature
96        MT=L_NTREF-1
97        TWL = TREF(L_NTREF)*1.00
98      ELSE
99        ! Regular case, find encompassing tabulated temperatures
[135]100        do n=1,L_NTREF-1
[3083]101          IF (TWL.GE.TREF(N) .and. TWL.LT.TREF(N+1)) THEN
[135]102            MT = n
[3083]103            EXIT
104          END IF
105        END DO
[135]106     
107      END IF
[3083]108!     weight for the temperature interpolation
109      U = (TWL-TREF(MT))/(TREF(MT+1)-TREF(MT))
[135]110
[1600]111
[3083]112!     2. For interpolation along pressure
113!        Get the upper and lower pressure grid indicies that bound the
114!        requested pressure. If the sought pressure is below
115!        the P-grid minimum, assume that lower value. If the sought pressure
116!        is above the P-grid maximum, assume that maximum value.
[135]117
[3083]118      PWL = log10(PW)
119      IF(PWL.LT.PREF(1)) THEN
120        ! Below lowest tabulated pressure
121        MP = 1
[1616]122        PWL=Pref(1)*1.00
[3083]123      ELSEIF(PWL.GE.PREF(L_PINT)) THEN
124        ! Above highest tabulated pressure
125        MP = L_PINT-1
126        PWL = PREF(L_PINT)*1.00
127      ELSE
128        ! Regular case, find encompassing tabulated pressures
129        do n=1,L_PINT-1
130          IF (PWL.GE.PREF(N) .and. PWL.LT.PREF(N+1)) THEN
131            MP = n
132            EXIT
133          END IF
134        END DO
[1600]135      ENDIF
[135]136
[3083]137!     weight for the pressure interpolation
[135]138      T = (PWL-PREF(MP))/(PREF(MP+1)-PREF(MP))
139
[3083]140!  Build the interpolation coefficients (bilinear in temperature-log(pressure))
[135]141      LCOEF(1) = (1.0-T)*(1.0-U)
142      LCOEF(2) = T*(1.0-U)
143      LCOEF(3) = T*U
144      LCOEF(4) = (1.0-T)*U
145
[3083]146      end subroutine tpindex
Note: See TracBrowser for help on using the repository browser.