MONC
Functions/Subroutines
interpolation_mod Module Reference

Functions/Subroutines

subroutine piecewise_linear_1d (zvals, vals, zgrid, field)
 Does a simple 1d piecewise linear interpolation. More...
 
subroutine interpolate_point_linear_1d (zvals, vals, z, f, extrapolate)
 Does a simple 1d linear interpolation to a point. More...
 

Function/Subroutine Documentation

◆ interpolate_point_linear_1d()

subroutine interpolation_mod::interpolate_point_linear_1d ( real(kind=default_precision), dimension(:), intent(in)  zvals,
real(kind=default_precision), dimension(:), intent(in)  vals,
real(kind=default_precision), intent(in)  z,
real(kind=default_precision), intent(out)  f,
character(*), intent(in), optional  extrapolate 
)

Does a simple 1d linear interpolation to a point.

Parameters
zvalsinput z nodes
valsinput nodal values
zlocation to interpolate onto
foutput interpolated value

Definition at line 72 of file interpolation.F90.

72 
73  real(kind=DEFAULT_PRECISION), intent(in) :: zvals(:), vals(:)
74  real(kind=DEFAULT_PRECISION), intent(in) :: z
75  real(kind=DEFAULT_PRECISION), intent(out) :: f
76  character(*), intent(in), optional :: extrapolate
77 
78  integer :: nn ! loop counter
79  integer :: nnodes ! number of input values
80 
81  integer, parameter :: maxchars=20
82  character(MAXCHARS) :: ext_type
83 
84  nnodes=size(zvals)
85  ! suggested fix from JME
86  do
87  if (nnodes == 1) exit
88  if (zvals(nnodes-1) < zvals(nnodes)) exit
89  nnodes = nnodes - 1
90  enddo
91  !
92  ext_type='linear'
93  if (present(extrapolate))ext_type=trim(extrapolate)
94 
95  if (z < zvals(1))then
96  nn=1
97  select case (trim(ext_type))
98  case ('linear')
99  f = vals(nn) + (vals(nn+1) - vals(nn))/(zvals(nn+1) - zvals(nn)) &
100  *(z - zvals(nn))
101  case ('constant')
102  f = vals(nn)
103  end select
104  return
105  end if
106 
107  if (present(extrapolate))ext_type=trim(extrapolate)
108 
109  if (z >= zvals(nnodes))then
110  nn=nnodes
111  select case (trim(ext_type))
112  case ('linear')
113  f = vals(nn) + (vals(nn-1) - vals(nn))/(zvals(nn-1) - zvals(nn)) &
114  *(z - zvals(nn))
115  case ('constant')
116  f = vals(nn)
117  end select
118  return
119  end if
120 
121  do nn=1,nnodes-1
122  if (zvals(nn) <= z .and. z < zvals(nn+1))then
123  f = vals(nn) + (vals(nn+1) - vals(nn))/(zvals(nn+1) - zvals(nn)) &
124  *(z - zvals(nn))
125  exit
126  end if
127  end do
128 
129 
Here is the caller graph for this function:

◆ piecewise_linear_1d()

subroutine interpolation_mod::piecewise_linear_1d ( real(kind=default_precision), dimension(:), intent(in)  zvals,
real(kind=default_precision), dimension(:), intent(in)  vals,
real(kind=default_precision), dimension(:), intent(inout)  zgrid,
real(kind=default_precision), dimension(:), intent(inout)  field 
)

Does a simple 1d piecewise linear interpolation.

Parameters
zvalsinput z nodes
valsinput nodal values
zgridgrid to interpolate onto
fieldoutput interpolated values

Definition at line 15 of file interpolation.F90.

15 
16  real(kind=DEFAULT_PRECISION), intent(in) :: zvals(:), vals(:)
17  real(kind=DEFAULT_PRECISION), intent(inout) :: zgrid(:)
18  real(kind=DEFAULT_PRECISION), intent(inout) :: field(:)
19 
20  real(kind=DEFAULT_PRECISION) :: initgd_lem(size(zvals)+1)
21  real(kind=DEFAULT_PRECISION) :: zngd_lem(size(zvals)+1)
22  real(kind=DEFAULT_PRECISION) :: field_lem(size(field))
23 
24  real(kind=DEFAULT_PRECISION) :: verylow, veryhigh
25 
26  integer :: nn,k ! loop counters
27  integer :: nnodes ! number of input values
28 
29  nnodes=size(zvals)
30 
31 !!$ do k=1, size(field)
32 !!$ do nn=2,nnodes
33 !!$ if (zgrid(k) < zvals(nn))then
34 !!$ field(k) = vals(nn-1) + (vals(nn) - vals(nn-1))/(zvals(nn) - zvals(nn-1)) &
35 !!$ *(zgrid(k) - zvals(nn-1))
36 !!$ exit
37 !!$ end if
38 !!$ end do
39 !!$ end do
40 
41  ! now repeat the code of the LEM as a check
42  verylow = -1.0e5
43 
44  initgd_lem(1) = vals(1)
45  zngd_lem(1) = verylow
46  DO nn=2,nnodes+1
47  initgd_lem(nn) = vals(nn-1)
48  zngd_lem(nn) = zvals(nn-1)
49  ENDDO
50 
51  do nn=1,nnodes
52  DO k=1,size(field)
53  IF( zngd_lem(nn).LE.zgrid(k) .AND. zgrid(k).LT.zngd_lem(nn+1)) then
54  field(k) = initgd_lem(nn) + ( initgd_lem(nn+1) - initgd_lem(nn))*(zgrid(k)- zngd_lem(nn))/ &
55  (zngd_lem(nn+1) - zngd_lem(nn))
56  end IF
57  end DO
58  end do
59 
60  !do k=1, size(field)
61  ! print *, 'field_monc, field_lem =', field(k), field_lem(k), k
62  !enddo
63 
Here is the caller graph for this function: