source: trunk/LMDZ.GENERIC/libf/phystd/tpindex.F @ 3436

Last change on this file since 3436 was 2822, checked in by aslmd, 2 years ago

Fixed a real conversion

File size: 5.7 KB
Line 
1      module tpindex_mod
2     
3      implicit none
4     
5      contains
6     
7      subroutine tpindex(pw,tw,qvar,pref,tref,wrefvar,LCOEF,MT,MP,
8     &     NVAR,wratio)
9
10!==================================================================
11!     
12!     Purpose
13!     -------
14!     Interpolate K-coefficients to the given P,T and Qvar values.
15!
16!     Notes
17!     -----
18!     The interpolation is the usual one in two dimensions given
19!     in "Numerical Recipes", where the "X" are P, the "Y" are
20!     T, and the F(X,Y) are the CO2 K-coefficients.
21!
22!     The interpolating box is:
23!
24!           (PL,TU)                        (PR,TU)
25!
26!                          (TW,PW)
27!
28!           
29!           (PL,TL)                        (PR,TL)
30!
31!      PL  - Pressure left
32!      PR  - Pressure right
33!      TL  - Temperature lower
34!      TU  - Temperature upper
35!      PW  - Pressure wanted
36!      TW  - Temperature wanted
37!
38!     Inputs
39!     ------
40!     PW                 - The pressure to interpolate to
41!     TW                 - The temperature to interpolate to
42!     Qvar               - The abundance to interpolate to
43!     Pref               - The tabulated (log) pressure grid
44!     Tref               - The tabulated temperature grid
45!     WREFVAR            - The tabulated abundance grid
46!   
47!     Outputs
48!     -------
49!     MT                 - Temperature index (bottom left temperature)
50!                          of bounding box
51!     MP                 - Pressure index (bottom left pressure)
52!                          of bounding box
53!     NVAR               -  index (bottom left )
54!                          of bounding box
55!     LCOEF(4)           - Interpolation coefficients
56!     WRATIO             - weight for the interpolation of abundance
57!
58!     Authors
59!     -------
60!     Adapted from the NASA Ames code by R. Wordsworth (2009)
61!     
62!==================================================================
63
64      use radinc_h, only: L_PINT, L_NTREF, L_REFVAR
65
66      implicit none
67
68      REAL,INTENT(IN) :: Tref(L_NTREF)
69      REAL,INTENT(IN) :: pref(L_PINT)
70      REAL,INTENT(IN) :: wrefvar(L_REFVAR)
71      REAL,INTENT(IN) :: PW, TW, Qvar
72
73      INTEGER,INTENT(OUT) :: MT, MP, NVAR
74      REAL,INTENT(OUT) :: LCOEF(4), wratio
75     
76      INTEGER :: N
77      REAL :: PWL ! local value for log10(PW)
78      REAL :: TWL ! local value for TW
79      REAL :: T ! linear interpolation weight along log(pressure)
80      REAL :: U ! linear interpolation weight along temperature
81
82!     1. For the interpolation along temperature
83!     Get the upper and lower temperature grid indicies that bound the
84!     requested temperature. If the sought temperature is below
85!     the T-grid minimum, assume that lower value. If the sought temperature
86!     is above the T-grid maximum, assume that maximum value.
87!     TW : temperature to interpolate to
88!     TREF : reference temperature array
89!     MT : index of TREF for (lower) bounding temperature
90!     U : weight for interpolation in temperature
91
92      TWL = TW
93
94      IF(TWL.LT.TREF(1)) THEN
95        ! Below lowest tabulated temperature
96        MT = 1
97        TWL=TREF(1)*1.00
98      ELSEIF(TWL.GE.TREF(L_NTREF)) THEN
99       ! Above highest tabulated temperature
100        MT=L_NTREF-1
101        TWL = TREF(L_NTREF)*1.00
102      ELSE
103        ! Regular case, find encompassing tabulated temperatures
104        do n=1,L_NTREF-1
105          IF (TWL.GE.TREF(N) .and. TWL.LT.TREF(N+1)) THEN
106            MT = n
107            EXIT
108          END IF
109        END DO
110     
111      END IF
112
113!     weight for the temperature interpolation
114      U = (TWL-TREF(MT))/(TREF(MT+1)-TREF(MT))
115
116!     2. For interpolation along pressure
117!     Get the upper and lower pressure grid indicies that bound the
118!     requested pressure. If the sought pressure is below
119!     the P-grid minimum, assume that lower value. If the sought pressure
120!     is above the P-grid maximum, assume that maximum value.
121
122      PWL = log10(PW)
123
124      IF(PWL.LT.PREF(1)) THEN
125        ! Below lowest tabulated pressure
126        MP = 1
127        PWL=Pref(1)*1.00
128      ELSEIF(PWL.GE.PREF(L_PINT)) THEN
129        ! Above highest tabulated pressure
130        MP = L_PINT-1
131        PWL = PREF(L_PINT)*1.00
132      ELSE
133        ! Regular case, find encompassing tabulated pressures
134        do n=1,L_PINT-1
135          IF (PWL.GE.PREF(N) .and. PWL.LT.PREF(N+1)) THEN
136            MP = n
137            EXIT
138          END IF
139        END DO
140     
141      END IF
142
143!     weight for the pressure interpolation
144      T = (PWL-PREF(MP))/(PREF(MP+1)-PREF(MP))
145
146!  Build the interpolation coefficients (bilinear in temperature-log(pressure))
147      LCOEF(1) = (1.0-T)*(1.0-U)
148      LCOEF(2) = T*(1.0-U)
149      LCOEF(3) = T*U
150      LCOEF(4) = (1.0-T)*U
151
152!  3. For interpolation along abundances
153!  Get the indicies for abundance of the varying species. There are 10 sets of
154!  k-coefficients with differing amounts of variable vs. constant gas.
155
156      IF(QVAR.LE.WREFVAR(1)) THEN
157        ! Below lowest tabulated abundance
158        NVAR   = 1
159        WRATIO = 0.0D0      ! put all the weight on the first point
160      ELSEIF(QVAR.GE.WREFVAR(L_REFVAR)) THEN
161        ! Above highest tabulated abundance
162        NVAR   = L_REFVAR-1 ! TB16 in order to not oversize NVAr when doing
163                                  !NVAR+1
164        WRATIO = 1.00D0     ! put all the weight on the last point
165      ELSE
166        ! Regular case, find encompassing tabulated abundances
167        ! and compute corresponding weight
168        DO N=2,L_REFVAR
169          IF(QVAR.GE.WREFVAR(N-1) .and. QVAR.LT.WREFVAR(N)) THEN
170            NVAR   = N-1
171            WRATIO = (QVAR - WREFVAR(N-1))/(WREFVAR(N) - WREFVAR(N-1))
172            EXIT
173          END IF
174        END DO
175      END IF
176
177      end subroutine tpindex
178     
179      end module tpindex_mod
180     
Note: See TracBrowser for help on using the repository browser.