- Timestamp:
- Apr 10, 2026, 5:31:47 PM (33 hours ago)
- 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) 1 module util_linear_interpolation_mod 2 3 implicit none 4 5 contains 6 7 subroutine interp_line(x1,y1,len1,x2,y2,len2) 2 8 implicit none 3 9 !----------------------------------------------------------------------- … … 10 16 ! the nearby x1(:) point 11 17 ! 12 c-----------------------------------------------------------------------18 !----------------------------------------------------------------------- 13 19 ! arguments 14 20 ! --------- 15 21 ! inputs: 16 real x1(len1) ! ordered list of abscissas17 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 done20 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(:) 21 27 ! outputs: 22 real y2(len2) ! interpolated values28 real,intent(out) :: y2(len2) ! interpolated values 23 29 !----------------------------------------------------------------------- 24 30 25 31 ! local variables: 26 integer i,j32 integer :: i,j 27 33 28 34 29 35 do i=1,len2 30 36 ! 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)))) then37 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 33 39 ! set y2(i) to y1(1) or y1(len1) 34 40 if (abs(x2(i)-x1(1)).lt.abs(x2(i)-x1(len1))) then … … 43 49 ! find the nearest neigbours and do a linear interpolation 44 50 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)))) then47 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) 49 55 endif 50 56 enddo 51 57 endif 52 58 53 enddo 59 enddo ! of do i=1,len2 54 60 55 end 61 end subroutine interp_line 62 63 end module util_linear_interpolation_mod
Note: See TracChangeset
for help on using the changeset viewer.
