Ignore:
Timestamp:
Apr 10, 2026, 5:31:47 PM (33 hours ago)
Author:
emillour
Message:

Generic PCM:
Clean up code: put inter_line.F in a module util_linear_interpolation.F90 and
also turn soil.F and soil_settings.F into modules.
EM

File:
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phygeneric/util_linear_interpolation.F90

    r4178 r4179  
    1       subroutine interp_line(x1,y1,len1,x2,y2,len2)
     1module util_linear_interpolation_mod
     2
     3implicit none
     4
     5contains
     6
     7    subroutine interp_line(x1,y1,len1,x2,y2,len2)
    28      implicit none
    39!-----------------------------------------------------------------------
     
    1016!  the nearby x1(:) point
    1117
    12 c-----------------------------------------------------------------------
     18!-----------------------------------------------------------------------
    1319!  arguments
    1420!  ---------
    1521!  inputs:
    16       real x1(len1) ! ordered list of abscissas
    17       real y1(len1) ! values at x1(:)
    18       integer len1  ! length of x1(:) and y1(:)
    19       real x2(len2) !ordered list of abscissas at which interpolation is done
    20       integer len2  ! length of x2(:) and y2(:)
     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(:)
    2127!  outputs:
    22       real y2(len2) ! interpolated values
     28      real,intent(out) :: y2(len2) ! interpolated values
    2329!-----------------------------------------------------------------------
    2430
    2531! local variables:
    26       integer i,j
     32      integer :: i,j
    2733     
    2834
    2935      do i=1,len2
    3036        ! check if x2(i) lies outside of the interval covered by x1()
    31         if(((x2(i).le.x1(1)).and.(x2(i).le.x1(len1))).or.
    32      &     ((x2(i).ge.x1(1)).and.(x2(i).ge.x1(len1)))) then
     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
    3339          ! set y2(i) to y1(1) or y1(len1)
    3440          if (abs(x2(i)-x1(1)).lt.abs(x2(i)-x1(len1))) then
     
    4349        ! find the nearest neigbours and do a linear interpolation
    4450         do j=1,len1-1
    45           if(((x2(i).ge.x1(j)).and.(x2(i).le.x1(j+1))).or.
    46      &       ((x2(i).le.x1(j)).and.(x2(i).ge.x1(j+1)))) then
    47             y2(i)=((x2(i)-x1(j))/(x1(j+1)-x1(j)))*y1(j+1)+
    48      &            ((x2(i)-x1(j+1))/(x1(j)-x1(j+1)))*y1(j)
     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)
    4955          endif
    5056         enddo
    5157        endif
    5258
    53       enddo
     59      enddo ! of do i=1,len2
    5460
    55       end
     61    end subroutine interp_line
     62
     63end module util_linear_interpolation_mod
Note: See TracChangeset for help on using the changeset viewer.