Afivo  0.3
m_a2_interp.f90
1 
5 module m_a2_interp
6  use m_a2_types
7 
8  implicit none
9  private
10 
11  public :: a2_interp0
12  public :: a2_interp1
13  public :: a2_interp0_to_grid
14 
15 contains
16 
18  function a2_interp0(tree, rr, ivs, n_var) result(vals)
19  use m_a2_utils, only: a2_get_loc
20  type(a2_t), intent(in) :: tree
21  real(dp), intent(in) :: rr(2)
22  integer, intent(in) :: n_var
23  integer, intent(in) :: ivs(n_var)
24  real(dp) :: vals(n_var)
25  type(a2_loc_t) :: loc
26 
27  loc = a2_get_loc(tree, rr)
28 
29  if (loc%id == -1) then
30  print *, "a2_interp0 error, no box at ", rr
31  stop
32  end if
33 
34  vals = tree%boxes(loc%id)%cc(loc%ix(1), loc%ix(2), ivs)
35  end function a2_interp0
36 
38  function a2_interp1(tree, r, ivs, n_var, id_guess) result(vals)
39  use m_a2_utils, only: a2_get_id_at
40  type(a2_t), intent(in) :: tree
41  real(dp), intent(in) :: r(2)
42  integer, intent(in) :: n_var
43  integer, intent(in) :: ivs(n_var)
44  integer, intent(in), optional :: id_guess
45  real(dp) :: vals(n_var)
46  integer :: i, iv, id, ix(2)
47  real(dp) :: r_loc(2), dvec(2), ovec(2), w(2, 2)
48 
49  id = a2_get_id_at(tree, r, id_guess)
50 
51  if (id <= af_no_box) then
52  print *, "a2_interp1: point outside domain", r
53  stop
54  end if
55 
56  ! Compute ix such that r lies between cell centers at ix and ix + 1
57  ix = nint((r - tree%boxes(id)%r_min) / tree%boxes(id)%dr)
58  r_loc = a2_r_cc(tree%boxes(id), ix)
59  dvec = r - r_loc
60 
61  ! Normalize dvec to a value [0, 1]
62  dvec = dvec / tree%boxes(id)%dr
63  ovec = 1 - dvec
64 
65  ! Compute weights of linear interpolation
66  w(1, 1) = ovec(1) * ovec(2)
67  w(2, 1) = dvec(1) * ovec(2)
68  w(1, 2) = ovec(1) * dvec(2)
69  w(2, 2) = dvec(1) * dvec(2)
70 
71  do i = 1, size(ivs)
72  iv = ivs(i)
73  vals(i) = sum(w * tree%boxes(id)%cc(ix(1):ix(1)+1, &
74  ix(2):ix(2)+1, iv))
75  end do
76  end function a2_interp1
77 
79  subroutine a2_interp0_to_grid(tree, rr, iv, amount, to_density)
80  use m_a2_utils, only: a2_get_loc
81  type(a2_t), intent(inout) :: tree
82  integer, intent(in) :: iv
83  real(dp), intent(in) :: rr(2)
84  real(dp), intent(in) :: amount
85  logical, intent(in) :: to_density
86  real(dp) :: actual_amount
87  type(a2_loc_t) :: loc
88  integer :: id, ix(2)
89 
90  loc = a2_get_loc(tree, rr)
91 
92  if (loc%id == -1) then
93  print *, "a2_interp0_to_grid error, no box at ", rr
94  stop
95  end if
96 
97  id = loc%id
98  ix = loc%ix
99 
101  if (to_density) then
102  actual_amount = amount / tree%boxes(id)%dr**2
103  else
104  actual_amount = amount
105  end if
106 
107  tree%boxes(id)%cc(ix(1), ix(2), iv) = &
108  tree%boxes(id)%cc(ix(1), ix(2), iv) + &
109  actual_amount
110  end subroutine a2_interp0_to_grid
111 
112 end module m_a2_interp
This module contains all kinds of different &#39;helper&#39; routines for Afivo. If the number of routines fo...
Definition: m_a2_utils.f90:5
Type which stores all the boxes and levels, as well as some information about the number of boxes...
Definition: m_a2_types.f90:93
This module contains routines related to interpolation, which can interpolate &#39;to&#39; the grid and &#39;from...
Definition: m_a2_interp.f90:5
This module contains the basic types and constants that are used in the 2-dimensional version of Afiv...
Definition: m_a2_types.f90:5
Type specifying the location of a cell.
Definition: m_a2_types.f90:124