ARTEMIS
sample.H
Go to the documentation of this file.
1 /* Copyright 2022 Edoardo Zoni, Remi Lehe, David Grote, Axel Huebl
2  *
3  * This file is part of ABLASTR.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef ABLASTR_COARSEN_SAMPLE_H_
8 #define ABLASTR_COARSEN_SAMPLE_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 {
50  amrex::Real Interp (
51  amrex::Array4<amrex::Real const> const& arr_src,
52  amrex::GpuArray<int,3> const& sf,
53  amrex::GpuArray<int,3> const& sc,
55  const int i,
56  const int j,
57  const int k,
58  const int comp
59  )
60  {
61  using namespace amrex::literals;
62 
63  // Indices of destination array (coarse)
64  const int ic[3] = { i, j, k };
65 
66  // Number of points and starting indices of source array (fine)
67  int np[3], idx_min[3];
68 
69  // Compute number of points
70  for ( int l = 0; l < 3; ++l ) {
71  if ( cr[l] == 1 ) np[l] = 1+amrex::Math::abs(sf[l]-sc[l]); // no coarsening
72  else np[l] = 2-sf[l];
73  }
74 
75  // Compute starting indices of source array (fine)
76  for ( int l = 0; l < 3; ++l ) {
77  if ( cr[l] == 1 ) idx_min[l] = ic[l]-sc[l]*(1-sf[l]); // no coarsening
78  else idx_min[l] = ic[l]*cr[l]+static_cast<int>(cr[l]/2)*(1-sc[l])-(1-sf[l]);
79  }
80 
81  // Auxiliary integer variables
82  const int numx = np[0];
83  const int numy = np[1];
84  const int numz = np[2];
85  const int imin = idx_min[0];
86  const int jmin = idx_min[1];
87  const int kmin = idx_min[2];
88  int ii, jj, kk;
89  amrex::Real const wx = 1.0_rt / static_cast<amrex::Real>(numx);
90  amrex::Real const wy = 1.0_rt / static_cast<amrex::Real>(numy);
91  amrex::Real const wz = 1.0_rt / static_cast<amrex::Real>(numz);
92 
93  // Interpolate over points computed above
94  amrex::Real c = 0.0_rt;
95  for (int kref = 0; kref < numz; ++kref) {
96  for (int jref = 0; jref < numy; ++jref) {
97  for (int iref = 0; iref < numx; ++iref) {
98  ii = imin+iref;
99  jj = jmin+jref;
100  kk = kmin+kref;
101  c += wx*wy*wz*arr_src(ii,jj,kk,comp);
102  }
103  }
104  }
105  return c;
106  }
107 
127  void Loop ( amrex::MultiFab& mf_dst,
128  const amrex::MultiFab& mf_src,
129  const int dcomp,
130  const int scomp,
131  const int ncomp,
132  const amrex::IntVect ngrow,
133  const amrex::IntVect crse_ratio=amrex::IntVect(1) );
134 
154  void Coarsen ( amrex::MultiFab& mf_dst,
155  const amrex::MultiFab& mf_src,
156  const int dcomp,
157  const int scomp,
158  const int ncomp,
159  const int ngrow,
160  const amrex::IntVect crse_ratio=amrex::IntVect(1) );
161  void Coarsen ( amrex::MultiFab& mf_dst,
162  const amrex::MultiFab& mf_src,
163  const int dcomp,
164  const int scomp,
165  const int ncomp,
166  const amrex::IntVect ngrowvect,
167  const amrex::IntVect crse_ratio=amrex::IntVect(1) );
168 
169 } // namespace ablastr::coarsen::sample
170 
171 #endif // ABLASTR_COARSEN_SAMPLE_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
Definition: sample.cpp:23
void Loop(amrex::MultiFab &mf_dst, const amrex::MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const amrex::IntVect ngrowvect, const amrex::IntVect crse_ratio)
Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contai...
Definition: sample.cpp:25
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, 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: sample.H:50
void Coarsen(amrex::MultiFab &mf_dst, const amrex::MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const int ngrow, const amrex::IntVect crse_ratio)
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in th...
Definition: sample.cpp:107
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