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

Last change on this file since 2810 was 2582, checked in by emillour, 4 years ago

Generic GCM:
Fixed bug in tpindex (for low temperatures, between first and second
reference temperatures, temperature was wrongly set to tref(1)).
The input temperature was also allowed to be modified by the routine,
which is probably not a good idea and no longer the case.
Took this opportunity to turn tpindex into a module.
GM+EM

File size: 6.0 KB
RevLine 
[2582]1      module tpindex_mod
2     
3      implicit none
4     
5      contains
6     
[135]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
[2582]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
[135]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
[2582]53!     NVAR               -  index (bottom left )
54!                          of bounding box
55!     LCOEF(4)           - Interpolation coefficients
56!     WRATIO             - weight for the interpolation of abundance
[135]57!
58!     Authors
59!     -------
60!     Adapted from the NASA Ames code by R. Wordsworth (2009)
61!     
62!==================================================================
63
[2582]64      use radinc_h, only: L_PINT, L_NTREF, L_REFVAR
[135]65
66      implicit none
67
[2582]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
[135]72
[2582]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
[135]81
[2582]82!     1. For the interpolation along temperature
[135]83!     Get the upper and lower temperature grid indicies that bound the
[2582]84!     requested temperature. If the sought temperature is below
85!     the T-grid minimum, assume that lower value. If the sought teperature
86!     is above the T-grid maximum, extrapolate.
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
[135]91
[2582]92      TWL = TW
93
94      IF(TWL.LT.TREF(1)) THEN
95        ! Below lowest tabulated temperature
[135]96        MT = 1
[2582]97        TWL=TREF(1)*1.00
98        write(*,*) 'tpindex: Caution! Temperature of upper levels lower
[1616]99     $ than ref temperature for k-coef: k-coeff fixed for upper levels'
[2582]100        write(*,*) "         TW=",TWL
101        write(*,*) "         TREF(1)=",TREF(1)
102      ELSEIF(TWL.GE.TREF(L_NTREF)) THEN
103       ! Above highest tabulated temperature
104       MT = L_NTREF-1
105      ELSE
106        ! Regular case, find encompassing tabulated temperatures
[135]107        do n=1,L_NTREF-1
[2582]108          IF (TWL.GE.TREF(N) .and. TWL.LT.TREF(N+1)) THEN
[135]109            MT = n
[2582]110            EXIT
111          END IF
112        END DO
[135]113     
114      END IF
115
[2582]116!     weight for the temperature interpolation
117      U = (TWL-TREF(MT))/(TREF(MT+1)-TREF(MT))
[1600]118
[2582]119!     2. For interpolation along pressure
[135]120!     Get the upper and lower pressure grid indicies that bound the
[2582]121!     requested pressure. If the sought pressure is below
122!     the P-grid minimum, assume that lower value. If the sought pressure
123!     is above the P-grid maximum, extrapolate.
[135]124
[2582]125      PWL = log10(PW)
[135]126
[2582]127      IF(PWL.LT.PREF(1)) THEN
128        ! Below lowest tabulated pressure
129        MP = 1
130        PWL=Pref(1)*1.00
131        write(*,*) 'tpindex: Caution! Pressure of upper levels lower
132     $than ref pressure for k-coef: k-coeff fixed for upper levels'
133        write(*,*) "         PWL=",PWL
134        write(*,*) "         PREF(1)=",PREF(1)
135      ELSEIF(PWL.GE.PREF(L_PINT)) THEN
136        ! Above highest tabulated pressure
137        MP = L_PINT-1
138      ELSE
139        ! Regular case, find encompassing tabulated pressures
140        do n=1,L_PINT-1
141          IF (PWL.GE.PREF(N) .and. PWL.LT.PREF(N+1)) THEN
142            MP = n
143            EXIT
144          END IF
145        END DO
[1600]146     
[2582]147      END IF
[135]148
[2582]149!     weight for the pressure interpolation
[135]150      T = (PWL-PREF(MP))/(PREF(MP+1)-PREF(MP))
151
[2582]152!  Build the interpolation coefficients (bilinear in temperature-log(pressure))
[135]153      LCOEF(1) = (1.0-T)*(1.0-U)
154      LCOEF(2) = T*(1.0-U)
155      LCOEF(3) = T*U
156      LCOEF(4) = (1.0-T)*U
157
[2582]158!  3. For interpolation along abundances
[135]159!  Get the indicies for abundance of the varying species. There are 10 sets of
160!  k-coefficients with differing amounts of variable vs. constant gas.
161
[2582]162      IF(QVAR.LE.WREFVAR(1)) THEN
163        ! Below lowest tabulated abundance
[135]164        NVAR   = 1
[1600]165        WRATIO = 0.0D0      ! put all the weight on the first point
[2582]166      ELSEIF(QVAR.GE.WREFVAR(L_REFVAR)) THEN
167        ! Above highest tabulated abundance
[1600]168        NVAR   = L_REFVAR-1 ! TB16 in order to not oversize NVAr when doing
169                                  !NVAR+1
170        WRATIO = 1.00D0     ! put all the weight on the last point
[135]171      ELSE
[2582]172        ! Regular case, find encompassing tabulated abundances
173        ! and compute corresponding weight
[135]174        DO N=2,L_REFVAR
[2582]175          IF(QVAR.GE.WREFVAR(N-1) .and. QVAR.LT.WREFVAR(N)) THEN
[135]176            NVAR   = N-1
177            WRATIO = (QVAR - WREFVAR(N-1))/(WREFVAR(N) - WREFVAR(N-1))
[2582]178            EXIT
[135]179          END IF
180        END DO
181      END IF
182
[2582]183      end subroutine tpindex
184     
185      end module tpindex_mod
186     
Note: See TracBrowser for help on using the repository browser.