WarpX
CoarsenMR.H
Go to the documentation of this file.
1 #ifndef WARPX_COARSEN_MR_H_
2 #define WARPX_COARSEN_MR_H_
3 
4 #include <AMReX_Array.H>
5 #include <AMReX_Array4.H>
6 #include <AMReX_Extension.H>
7 #include <AMReX_GpuQualifiers.H>
8 #include <AMReX_REAL.H>
9 
10 #include <AMReX_BaseFwd.H>
11 
12 #include <cstdlib>
13 
14 namespace CoarsenMR{
15 
16  using namespace amrex;
17 
37  Real Interp ( Array4<Real const> const& arr_src,
38  GpuArray<int,3> const& sf,
39  GpuArray<int,3> const& sc,
40  GpuArray<int,3> const& cr,
41  const int i,
42  const int j,
43  const int k,
44  const int comp )
45  {
46  // Indices of destination array (coarse)
47  const int ic[3] = { i, j, k };
48 
49  // Number of points and starting indices of source array (fine)
50  int np[3], idx_min[3];
51 
52  // Compute number of points
53  for ( int l = 0; l < 3; ++l ) {
54  if ( cr[l] == 1 ) np[l] = 1; // no coarsening
55  else np[l] = cr[l]*(1-sf[l])*(1-sc[l]) // cell-centered
56  +(2*(cr[l]-1)+1)*sf[l]*sc[l]; // nodal
57  }
58 
59  // Compute starting indices of source array (fine)
60  for ( int l = 0; l < 3; ++l ) {
61  if ( cr[l] == 1 ) idx_min[l] = ic[l]; // no coarsening
62  else idx_min[l] = ic[l]*cr[l]*(1-sf[l])*(1-sc[l]) // cell-centered
63  +(ic[l]*cr[l]-cr[l]+1)*sf[l]*sc[l]; // nodal
64  }
65 
66  // Auxiliary integer variables
67  const int numx = np[0];
68  const int numy = np[1];
69  const int numz = np[2];
70  const int imin = idx_min[0];
71  const int jmin = idx_min[1];
72  const int kmin = idx_min[2];
73  const int sfx = sf[0];
74  const int sfy = sf[1];
75  const int sfz = sf[2];
76  const int scx = sc[0];
77  const int scy = sc[1];
78  const int scz = sc[2];
79  const int crx = cr[0];
80  const int cry = cr[1];
81  const int crz = cr[2];
82  int ii, jj, kk;
83  Real wx, wy, wz;
84 
85  // Add neutral elements (=0) beyond guard cells in source array (fine)
86  auto const arr_src_safe = [arr_src]
87  AMREX_GPU_DEVICE (int const ix, int const iy, int const iz, int const n) noexcept
88  {
89  return arr_src.contains( ix, iy, iz ) ? arr_src(ix,iy,iz,n) : 0.0_rt;
90  };
91 
92  // Interpolate over points computed above. Weights are computed in order
93  // to guarantee total charge conservation for both cell-centered data
94  // (equal weights) and nodal data (weights depend on distance between
95  // points on fine and coarse grids). Terms multiplied by (1-sf)*(1-sc)
96  // are ON for cell-centered data and OFF for nodal data, while terms
97  // multiplied by sf*sc are ON for nodal data and OFF for cell-centered data.
98  // Python script Source/Utils/check_interp_points_and_weights.py can be
99  // used to check interpolation points and weights in 1D.
100  Real c = 0.0_rt;
101  for (int kref = 0; kref < numz; ++kref) {
102  for (int jref = 0; jref < numy; ++jref) {
103  for (int iref = 0; iref < numx; ++iref) {
104  ii = imin+iref;
105  jj = jmin+jref;
106  kk = kmin+kref;
107  wx = (1.0_rt/static_cast<Real>(numx))*(1-sfx)*(1-scx) // if cell-centered
108  +((amrex::Math::abs(crx-amrex::Math::abs(ii-i*crx)))/static_cast<Real>(crx*crx))*sfx*scx; // if nodal
109  wy = (1.0_rt/static_cast<Real>(numy))*(1-sfy)*(1-scy) // if cell-centered
110  +((amrex::Math::abs(cry-amrex::Math::abs(jj-j*cry)))/static_cast<Real>(cry*cry))*sfy*scy; // if nodal
111  wz = (1.0_rt/static_cast<Real>(numz))*(1-sfz)*(1-scz) // if cell-centered
112  +((amrex::Math::abs(crz-amrex::Math::abs(kk-k*crz)))/static_cast<Real>(crz*crz))*sfz*scz; // if nodal
113  c += wx*wy*wz*arr_src_safe(ii,jj,kk,comp);
114  }
115  }
116  }
117  return c;
118  }
119 
133  void Loop ( MultiFab& mf_dst,
134  const MultiFab& mf_src,
135  const int ncomp,
136  const IntVect ngrow,
137  const IntVect crse_ratio );
138 
149  void Coarsen ( MultiFab& mf_dst,
150  const MultiFab& mf_src,
151  const IntVect crse_ratio );
152 }
153 
154 #endif // WARPX_COARSEN_MR_H_
void Loop(MultiFab &mf_dst, const MultiFab &mf_src, const int ncomp, const IntVect ngrow, const IntVect crse_ratio)
Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contai...
Definition: CoarsenMR.cpp:19
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
imin
Definition: check_interp_points_and_weights.py:123
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
ii
Definition: check_interp_points_and_weights.py:147
int n
Definition: run_libensemble_on_warpx.py:67
void Coarsen(MultiFab &mf_dst, const MultiFab &mf_src, const IntVect crse_ratio)
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in th...
Definition: CoarsenMR.cpp:90
i
Definition: check_interp_points_and_weights.py:173
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool contains(int i, int j, int k) const noexcept
jj
Definition: check_interp_points_and_weights.py:159
cr
TODO Coarsening for IO: interpolation points and weights def coarsening_points_and_weights_for_IO( i...
Definition: check_interp_points_and_weights.py:100
Definition: CoarsenMR.H:14
AMREX_GPU_DEVICE AMREX_FORCE_INLINE Real Interp(Array4< Real const > const &arr_src, GpuArray< int, 3 > const &sf, GpuArray< int, 3 > const &sc, GpuArray< int, 3 > const &cr, const int i, const int j, const int k, const int comp)
Interpolates the floating point data contained in the source Array4 arr_src, extracted from a fine Mu...
Definition: CoarsenMR.H:37