source: trunk/LMDZ.MARS/libf/phymars/interp_line.F @ 3902

Last change on this file since 3902 was 3902, checked in by emillour, 3 months ago

Mars PCM:
More code tidying: turn geticecover, interp_line, orbite, simpleclouds, tabfi
and tcondco2 into modules.
EM

File size: 2.0 KB
Line 
1      module interp_line_mod
2     
3      implicit none
4     
5      contains
6     
7      subroutine interp_line(x1,y1,len1,x2,y2,len2)
8      implicit none
9!-----------------------------------------------------------------------
10!
11!  Purpose: Do a series of linear interpolations
12!  Data sets are organized as vectors (see below)
13!  If x2(:), and abscissa at which interpolation is requiered, lies
14!  outside of the interval covered by x1(:), instead of doing an
15!  extrapolation, y2() is set to the value y1() corresponding to
16!  the nearby x1(:) point
17
18c-----------------------------------------------------------------------
19!  arguments
20!  ---------
21!  inputs:
22      real,intent(in) :: x1(len1) ! ordered list of abscissas
23      real,intent(in) :: y1(len1) ! values at x1(:)
24      integer,intent(in) :: len1  ! length of x1(:) and y1(:)
25      real,intent(in) :: x2(len2) !ordered list of abscissas at which interpolation is done
26      integer,intent(in) :: len2  ! length of x2(:) and y2(:)
27!  outputs:
28      real,intent(out) :: y2(len2) ! interpolated values
29!-----------------------------------------------------------------------
30
31! local variables:
32      integer i,j
33     
34
35      do i=1,len2
36        ! check if x2(i) lies outside of the interval covered by x1()
37        if(((x2(i).le.x1(1)).and.(x2(i).le.x1(len1))).or.
38     &     ((x2(i).ge.x1(1)).and.(x2(i).ge.x1(len1)))) then
39          ! set y2(i) to y1(1) or y1(len1)
40          if (abs(x2(i)-x1(1)).lt.abs(x2(i)-x1(len1))) then
41            ! x2(i) lies closest to x1(1)
42            y2(i)=y1(1)
43          else
44            ! x2(i) lies closest to x1(len1)
45            y2(i)=y1(len1)
46          endif
47
48        else
49        ! find the nearest neigbours and do a linear interpolation
50         do j=1,len1-1
51          if(((x2(i).ge.x1(j)).and.(x2(i).le.x1(j+1))).or.
52     &       ((x2(i).le.x1(j)).and.(x2(i).ge.x1(j+1)))) then
53            y2(i)=((x2(i)-x1(j))/(x1(j+1)-x1(j)))*y1(j+1)+
54     &            ((x2(i)-x1(j+1))/(x1(j)-x1(j+1)))*y1(j)
55          endif
56         enddo
57        endif
58
59      enddo
60
61      end subroutine interp_line
62
63      end module interp_line_mod
Note: See TracBrowser for help on using the repository browser.