source: trunk/LMDZ.GENERIC/libf/phystd/bilinearbig.F90 @ 879

Last change on this file since 879 was 878, checked in by jleconte, 12 years ago

12/02/2013 == JL

  • Follows previous commit by Aymeric about bilinear interpolations:
    • Extended to all existing continua
    • generalized bilinearbig to work for various size inputs
  • because N2 and H2O continua databases are smaller, improvement around 15% for

an earth case.

  • Property svn:executable set to *
File size: 2.0 KB
Line 
1      subroutine bilinearbig(nX,nY,x_arr,y_arr,f2d_arr,x_in,y_in,f,ind)
2
3!     Necessary for interpolation of continuum data
4!     optimized by A. Spiga 01/2013
5
6      implicit none
7
8      integer nX,nY,i,j,ind,b
9
10      real*8 x_in,y_in,x1,x2,y1,y2
11      real*8 f,f11,f12,f21,f22,fA,fB
12      real*8 x_arr(nX)
13      real*8 y_arr(nY)
14      real*8 f2d_arr(nX,nY)
15      real*8,save :: x,y
16
17      integer strlen
18      character*100 label
19      label='subroutine bilinear'
20
21      x=x_in
22      y=y_in
23
24
25   !! AS: important to optimize here because the array is quite large
26   !! ... and actually calculations only need to be done once
27   !! --> Case 1 : we have not calculated yet
28   if ( ind == -9999) then
29      !1st check we're within the wavenumber range
30      if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then
31         f=0.0D+0
32         ind=-1
33      else
34        i=1
35        x2=x_arr(i)
36        do while ( x2 .le. x )
37          x1=x2
38          i=i+1
39          x2=x_arr(i)
40          ind=i-1
41        end do
42      endif
43   !! --> case 2 : we already saw we are out of wavenumber range
44   else if ( ind == -1) then
45      f=0.0D+0
46      return
47   !! --> case 3 : we already determined a -- so we just proceed
48   else
49      x1=x_arr(ind)
50      x2=x_arr(ind+1)
51   endif
52
53!     ... and for y within the temperature range
54      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
55         write(*,*) 'Warning from bilinearH2H2:'
56         write(*,*) 'Outside continuum temperature range!'
57         if(y.lt.y_arr(1))then
58            y=y_arr(1)+0.01
59         endif
60         if(y.gt.y_arr(nY))then
61            y=y_arr(nY)-0.01
62         endif
63      else
64        j=1
65        y2=y_arr(j)
66        do while ( y2 .le. y )
67          y1=y2
68          j=j+1
69          y2=y_arr(j)
70          b=j-1
71        end do
72      endif
73
74      f11=f2d_arr(ind,b)
75      f21=f2d_arr(ind+1,b)
76      f12=f2d_arr(ind,b+1)
77      f22=f2d_arr(ind+1,b+1)
78
79      call bilinear(f,f11,f21,f12,f22,x,x1,x2,y,y1,y2)
80
81      return
82    end subroutine bilinearbig
Note: See TracBrowser for help on using the repository browser.