Afivo  0.3
m_af_interp.f90
1 #include "cpp_macros.h"
2 
5 module m_af_interp
6  use m_af_types
7 
8  implicit none
9  private
10 
11  public :: af_interp0
12  public :: af_interp1
13  public :: af_interp0_to_grid
14  public :: af_interp1_fc
15 
16 contains
17 
19  function af_interp0(tree, r, ivs, success, id_guess) result(vals)
20  use m_af_utils, only: af_get_loc
21  type(af_t), intent(in) :: tree
22  real(dp), intent(in) :: r(ndim)
23  integer, intent(in) :: ivs(:)
24  logical, intent(out) :: success
25  integer, intent(inout), optional :: id_guess
26  real(dp) :: vals(size(ivs))
27  type(af_loc_t) :: loc
28 
29  loc = af_get_loc(tree, r, guess=id_guess)
30  ! Update guess
31  if (present(id_guess)) id_guess = loc%id
32 
33  if (loc%id == -1) then
34  success = .false.
35  vals = 0.0_dp
36  else
37  success = .true.
38  vals = tree%boxes(loc%id)%cc(dindex(loc%ix), ivs)
39  end if
40  end function af_interp0
41 
43  function af_interp1(tree, r, ivs, success, id_guess) result(vals)
44  use m_af_utils, only: af_get_id_at
45  type(af_t), intent(in) :: tree
46  real(dp), intent(in) :: r(ndim)
47  integer, intent(in) :: ivs(:)
48  logical, intent(out) :: success
49  integer, intent(inout), optional :: id_guess
50  real(dp) :: vals(size(ivs))
51  integer :: i, iv, id, ix(ndim)
52  real(dp) :: r_loc(ndim), dvec(ndim), ovec(ndim), w(dtimes(2))
53 
54  id = af_get_id_at(tree, r, guess=id_guess)
55  ! Update guess
56  if (present(id_guess)) id_guess = id
57 
58  if (id <= af_no_box) then
59  success = .false.
60  vals = 0.0_dp
61  else
62  success = .true.
63  ! Compute ix such that r lies between cell centers at ix and ix + 1
64  ix = nint((r - tree%boxes(id)%r_min) / tree%boxes(id)%dr)
65  r_loc = af_r_cc(tree%boxes(id), ix)
66  dvec = r - r_loc
67 
68  ! Normalize dvec to a value [0, 1]
69  dvec = dvec / tree%boxes(id)%dr
70  ovec = 1 - dvec
71 
72  ! Compute weights of linear interpolation
73 #if NDIM == 1
74  w(1) = ovec(1)
75  w(2) = dvec(1)
76 #elif NDIM == 2
77  w(1, 1) = ovec(1) * ovec(2)
78  w(2, 1) = dvec(1) * ovec(2)
79  w(1, 2) = ovec(1) * dvec(2)
80  w(2, 2) = dvec(1) * dvec(2)
81 #elif NDIM == 3
82  w(1, 1, 1) = ovec(1) * ovec(2) * ovec(3)
83  w(2, 1, 1) = dvec(1) * ovec(2) * ovec(3)
84  w(1, 2, 1) = ovec(1) * dvec(2) * ovec(3)
85  w(2, 2, 1) = dvec(1) * dvec(2) * ovec(3)
86  w(1, 1, 2) = ovec(1) * ovec(2) * dvec(3)
87  w(2, 1, 2) = dvec(1) * ovec(2) * dvec(3)
88  w(1, 2, 2) = ovec(1) * dvec(2) * dvec(3)
89  w(2, 2, 2) = dvec(1) * dvec(2) * dvec(3)
90 #endif
91 
92  do i = 1, size(ivs)
93  iv = ivs(i)
94 #if NDIM == 1
95  vals(i) = sum(w * tree%boxes(id)%cc(ix(1):ix(1)+1, iv))
96 #elif NDIM == 2
97  vals(i) = sum(w * tree%boxes(id)%cc(ix(1):ix(1)+1, &
98  ix(2):ix(2)+1, iv))
99 #elif NDIM == 3
100  vals(i) = sum(w * tree%boxes(id)%cc(ix(1):ix(1)+1, &
101  ix(2):ix(2)+1, ix(3):ix(3)+1, iv))
102 #endif
103  end do
104  end if
105  end function af_interp1
106 
108  subroutine af_interp0_to_grid(tree, r, iv, amount, to_density)
109  use m_af_utils, only: af_get_loc
110  type(af_t), intent(inout) :: tree
111  integer, intent(in) :: iv
112  real(dp), intent(in) :: r(ndim)
113  real(dp), intent(in) :: amount
114  logical, intent(in) :: to_density
115  real(dp) :: actual_amount
116  type(af_loc_t) :: loc
117  integer :: id, ix(ndim)
118 
119  loc = af_get_loc(tree, r)
120 
121  if (loc%id == -1) then
122  print *, "af_interp0_to_grid error, no box at ", r
123  stop
124  end if
125 
126  id = loc%id
127  ix = loc%ix
128 
130  if (to_density) then
131  actual_amount = amount / product(tree%boxes(id)%dr)
132  else
133  actual_amount = amount
134  end if
135 
136  tree%boxes(id)%cc(dindex(ix), iv) = &
137  tree%boxes(id)%cc(dindex(ix), iv) + &
138  actual_amount
139  end subroutine af_interp0_to_grid
140 
142  function af_interp1_fc(tree, r, ifc, success, id_guess) result(vals)
143  use m_af_utils, only: af_get_id_at
144  type(af_t), intent(in) :: tree
145  real(dp), intent(in) :: r(ndim)
146  integer, intent(in) :: ifc
147  logical, intent(out) :: success
148  integer, intent(inout), optional :: id_guess
149  real(dp) :: vals(ndim), inv_dr(ndim)
150  integer :: id, ix(ndim)
151  real(dp) :: ix_frac(ndim), r_rel(ndim)
152 
153  id = af_get_id_at(tree, r, guess=id_guess)
154 
155  ! Update guess
156  if (present(id_guess)) id_guess = id
157 
158  if (id <= af_no_box) then
159  success = .false.
160  vals = 0.0_dp
161  else
162  success = .true.
163  inv_dr = 1/tree%boxes(id)%dr
164 
165  r_rel = r - tree%boxes(id)%r_min
166  ix_frac = r_rel * inv_dr + 1
167  ix = floor(ix_frac)
168  where (ix < 1) ix = 1
169  where (ix > tree%n_cell) ix = tree%n_cell
170  ix_frac = ix_frac - ix
171 
172 #if NDIM == 2
173  vals(1) = (1 - ix_frac(1)) * tree%boxes(id)%fc(ix(1), ix(2), 1, ifc) + &
174  ix_frac(1) * tree%boxes(id)%fc(ix(1)+1, ix(2), 1, ifc)
175  vals(2) = (1 - ix_frac(2)) * tree%boxes(id)%fc(ix(1), ix(2), 2, ifc) + &
176  ix_frac(2) * tree%boxes(id)%fc(ix(1), ix(2)+1, 2, ifc)
177 #elif NDIM == 3
178  vals(1) = (1 - ix_frac(1)) * tree%boxes(id)%fc(ix(1), ix(2), ix(3), 1, ifc) + &
179  ix_frac(1) * tree%boxes(id)%fc(ix(1)+1, ix(2), ix(3), 1, ifc)
180  vals(2) = (1 - ix_frac(2)) * tree%boxes(id)%fc(ix(1), ix(2), ix(3), 2, ifc) + &
181  ix_frac(2) * tree%boxes(id)%fc(ix(1), ix(2)+1, ix(3), 2, ifc)
182  vals(3) = (1 - ix_frac(3)) * tree%boxes(id)%fc(ix(1), ix(2), ix(3), 3, ifc) + &
183  ix_frac(3) * tree%boxes(id)%fc(ix(1), ix(2), ix(3)+1, 3, ifc)
184 #endif
185  end if
186  end function af_interp1_fc
187 
188 end module m_af_interp
This module contains routines related to interpolation, which can interpolate 'to' the grid and 'from...
Definition: m_af_interp.f90:5
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Definition: m_af_types.f90:4
This module contains all kinds of different 'helper' routines for Afivo. If the number of routines fo...
Definition: m_af_utils.f90:5
Type specifying the location of a cell.
Definition: m_af_types.f90:396
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
Definition: m_af_types.f90:326