ARTEMIS
average.H
Go to the documentation of this file.
1 /* Copyright 2022 Edoardo Zoni, Remi Lehe, Prabhat Kumar, Axel Huebl
2  *
3  * This file is part of ABLASTR.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef ABLASTR_COARSEN_AVERAGE_H_
8 #define ABLASTR_COARSEN_AVERAGE_H_
9 
10 
11 #include <AMReX_Array.H>
12 #include <AMReX_Array4.H>
13 #include <AMReX_BLassert.H>
14 #include <AMReX_Extension.H>
15 #include <AMReX_GpuQualifiers.H>
16 #include <AMReX_IntVect.H>
17 #include <AMReX_Math.H>
18 #include <AMReX_REAL.H>
19 
20 #include <AMReX_BaseFwd.H>
21 
22 #include <cstdlib>
23 
24 
30 {
52  amrex::Real
54  amrex::Array4<amrex::Real const> const &arr_src,
55  amrex::GpuArray<int, 3> const &sf,
56  amrex::GpuArray<int, 3> const &sc,
58  int const i,
59  int const j,
60  int const k,
61  int const comp
62  )
63  {
64  using namespace amrex::literals;
65 
66  AMREX_ASSERT_WITH_MESSAGE(sf[0] == sc[0], "Interp: Staggering for component 0 does not match!");
67  AMREX_ASSERT_WITH_MESSAGE(sf[1] == sc[1], "Interp: Staggering for component 1 does not match!");
68  AMREX_ASSERT_WITH_MESSAGE(sf[2] == sc[2], "Interp: Staggering for component 2 does not match!");
69 
70  // Indices of destination array (coarse)
71  int const ic[3] = {i, j, k};
72 
73  // Number of points and starting indices of source array (fine)
74  int np[3], idx_min[3];
75 
76  // Compute number of points
77  for (int l = 0; l < 3; ++l) {
78  if (cr[l] == 1)
79  np[l] = 1; // no coarsening
80  else
81  np[l] = cr[l] * (1 - sf[l]) * (1 - sc[l]) // cell-centered
82  + (2 * (cr[l] - 1) + 1) * sf[l] * sc[l]; // nodal
83  }
84 
85  // Compute starting indices of source array (fine)
86  for (int l = 0; l < 3; ++l) {
87  if (cr[l] == 1)
88  idx_min[l] = ic[l]; // no coarsening
89  else
90  idx_min[l] = ic[l] * cr[l] * (1 - sf[l]) * (1 - sc[l]) // cell-centered
91  + (ic[l] * cr[l] - cr[l] + 1) * sf[l] * sc[l]; // nodal
92  }
93 
94  // Auxiliary integer variables
95  int const numx = np[0];
96  int const numy = np[1];
97  int const numz = np[2];
98  int const imin = idx_min[0];
99  int const jmin = idx_min[1];
100  int const kmin = idx_min[2];
101  int const sfx = sf[0];
102  int const sfy = sf[1];
103  int const sfz = sf[2];
104  int const scx = sc[0];
105  int const scy = sc[1];
106  int const scz = sc[2];
107  int const crx = cr[0];
108  int const cry = cr[1];
109  int const crz = cr[2];
110  int ii, jj, kk;
111  amrex::Real wx, wy, wz;
112 
113  // Add neutral elements (=0) beyond guard cells in source array (fine)
114  auto const arr_src_safe = [arr_src]
115  AMREX_GPU_DEVICE(int const ix, int const iy, int const iz, int const n) noexcept {
116  return arr_src.contains(ix, iy, iz) ? arr_src(ix, iy, iz, n) : 0.0_rt;
117  };
118 
119  // Interpolate over points computed above. Weights are computed in order
120  // to guarantee total charge conservation for both cell-centered data
121  // (equal weights) and nodal data (weights depend on distance between
122  // points on fine and coarse grids). Terms multiplied by (1-sf)*(1-sc)
123  // are ON for cell-centered data and OFF for nodal data, while terms
124  // multiplied by sf*sc are ON for nodal data and OFF for cell-centered data.
125  // Python script Source/Utils/check_interp_points_and_weights.py can be
126  // used to check interpolation points and weights in 1D.
127  amrex::Real c = 0.0_rt;
128  for (int kref = 0; kref < numz; ++kref) {
129  for (int jref = 0; jref < numy; ++jref) {
130  for (int iref = 0; iref < numx; ++iref) {
131  ii = imin + iref;
132  jj = jmin + jref;
133  kk = kmin + kref;
134  wx = (1.0_rt / static_cast<amrex::Real>(numx)) * (1 - sfx) * (1 - scx) // if cell-centered
135  + ((amrex::Math::abs(crx - amrex::Math::abs(ii - i * crx))) /
136  static_cast<amrex::Real>(crx * crx)) * sfx * scx; // if nodal
137  wy = (1.0_rt / static_cast<amrex::Real>(numy)) * (1 - sfy) * (1 - scy) // if cell-centered
138  + ((amrex::Math::abs(cry - amrex::Math::abs(jj - j * cry))) /
139  static_cast<amrex::Real>(cry * cry)) * sfy * scy; // if nodal
140  wz = (1.0_rt / static_cast<amrex::Real>(numz)) * (1 - sfz) * (1 - scz) // if cell-centered
141  + ((amrex::Math::abs(crz - amrex::Math::abs(kk - k * crz))) /
142  static_cast<amrex::Real>(crz * crz)) * sfz * scz; // if nodal
143  c += wx * wy * wz * arr_src_safe(ii, jj, kk, comp);
144  }
145  }
146  }
147  return c;
148  }
149 
163  void
164  Loop (
165  amrex::MultiFab & mf_dst,
166  amrex::MultiFab const & mf_src,
167  int const ncomp,
168  amrex::IntVect const ngrow,
169  amrex::IntVect const crse_ratio
170  );
171 
182  void
183  Coarsen (
184  amrex::MultiFab & mf_dst,
185  amrex::MultiFab const & mf_src,
186  amrex::IntVect const crse_ratio
187  );
188 
189 } // namespace ablastr::coarsen::average
190 
191 #endif // ABLASTR_COARSEN_AVERAGE_H_
#define AMREX_ASSERT_WITH_MESSAGE(EX, MSG)
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
Definition: average.cpp:23
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Interp(amrex::Array4< amrex::Real const > const &arr_src, amrex::GpuArray< int, 3 > const &sf, amrex::GpuArray< int, 3 > const &sc, amrex::GpuArray< int, 3 > const &cr, int const i, int const j, int const k, int const comp)
Interpolates the floating point data contained in the source Array4 arr_src, extracted from a fine Mu...
Definition: average.H:53
void Coarsen(amrex::MultiFab &mf_dst, amrex::MultiFab const &mf_src, amrex::IntVect const crse_ratio)
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in th...
Definition: average.cpp:96
void Loop(amrex::MultiFab &mf_dst, amrex::MultiFab const &mf_src, int const ncomp, amrex::IntVect const ngrow, amrex::IntVect const crse_ratio)
Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contai...
Definition: average.cpp:25
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
i
Definition: check_interp_points_and_weights.py:174
ii
Definition: check_interp_points_and_weights.py:148
jj
Definition: check_interp_points_and_weights.py:160
imin
Definition: check_interp_points_and_weights.py:124
cr
TODO Coarsening for IO: interpolation points and weights def coarsening_points_and_weights_for_IO( i,...
Definition: check_interp_points_and_weights.py:101
int n
Definition: run_libensemble_on_warpx.py:67
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool contains(int i, int j, int k) const noexcept