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

Last change on this file since 910 was 903, checked in by aslmd, 12 years ago

LMDZ.GENERIC. Corrected a small bug in bilinearbig. Added support for ciclad in makegcm_ifort. Raised the value for capcal in case soil turned off.

  • 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   !! AS: important to optimize here because the array is quite large
25   !! ... and actually calculations only need to be done once
26   !! IF ind=-9999 we have not calculated yet
27   if ( ind == -9999) then
28      !1st check we're within the wavenumber range
29      if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then
30         ind=-1
31      else
32        i=1
33        x2=x_arr(i)
34        do while ( x2 .le. x )
35          x1=x2
36          i=i+1
37          x2=x_arr(i)
38          ind=i-1
39        end do
40      endif
41   endif
42
43   !! Either we already saw we are out of wavenumber range
44   !! ... and we just have to set f=0 and exit
45   if ( ind == -1) then
46      f=0.0D+0
47      return
48   !! Or we already determined ind -- so we just proceed
49   else
50      x1=x_arr(ind)
51      x2=x_arr(ind+1)
52   endif
53
54!     ... and for y within the temperature range
55      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
56         write(*,*) 'Warning from bilinearH2H2:'
57         write(*,*) 'Outside continuum temperature range!'
58         if(y.lt.y_arr(1))then
59            y=y_arr(1)+0.01
60         endif
61         if(y.gt.y_arr(nY))then
62            y=y_arr(nY)-0.01
63         endif
64      else
65        j=1
66        y2=y_arr(j)
67        do while ( y2 .le. y )
68          y1=y2
69          j=j+1
70          y2=y_arr(j)
71          b=j-1
72        end do
73      endif
74
75      f11=f2d_arr(ind,b)
76      f21=f2d_arr(ind+1,b)
77      f12=f2d_arr(ind,b+1)
78      f22=f2d_arr(ind+1,b+1)
79
80      call bilinear(f,f11,f21,f12,f22,x,x1,x2,y,y1,y2)
81
82      return
83    end subroutine bilinearbig
Note: See TracBrowser for help on using the repository browser.