1#ifndef AMREX_BASEFAB_H_
2#define AMREX_BASEFAB_H_
3#include <AMReX_Config.H>
44extern std::atomic<Long> atomic_total_bytes_allocated_in_fabs;
45extern std::atomic<Long> atomic_total_bytes_allocated_in_fabs_hwm;
46extern std::atomic<Long> atomic_total_cells_allocated_in_fabs;
47extern std::atomic<Long> atomic_total_cells_allocated_in_fabs_hwm;
48extern Long private_total_bytes_allocated_in_fabs;
49extern Long private_total_bytes_allocated_in_fabs_hwm;
50extern Long private_total_cells_allocated_in_fabs;
51extern Long private_total_cells_allocated_in_fabs_hwm;
53#pragma omp threadprivate(private_total_bytes_allocated_in_fabs)
54#pragma omp threadprivate(private_total_bytes_allocated_in_fabs_hwm)
55#pragma omp threadprivate(private_total_cells_allocated_in_fabs)
56#pragma omp threadprivate(private_total_cells_allocated_in_fabs_hwm)
73 explicit SrcComp (
int ai) noexcept : i(ai) {}
98requires (std::is_arithmetic_v<T>)
104requires (std::is_trivially_default_constructible_v<T> && !std::is_arithmetic_v<T>)
108 for (
Long i = 0; i < n; ++i) {
114requires (!std::is_trivially_default_constructible_v<T>)
125requires (std::is_trivially_destructible_v<T>)
131requires (!std::is_trivially_destructible_v<T>)
210 bool shared = false,
Arena* ar =
nullptr);
261 requires (std::is_trivially_destructible_v<U>)
274 [[nodiscard]] std::
size_t nBytes () const noexcept {
return this->
truesize*
sizeof(T); }
281 [[nodiscard]] std::size_t
nBytes (
const Box& bx,
int ncomps)
const noexcept
282 {
return bx.numPts() *
sizeof(T) * ncomps; }
285 [[nodiscard]]
int nComp () const noexcept {
return this->
nvar; }
288 [[nodiscard]]
const int*
nCompPtr() const noexcept {
289 return &(this->
nvar);
299 [[nodiscard]]
const Box&
box () const noexcept {
return this->
domain; }
360 [[nodiscard]] T*
dataPtr (
int n = 0) noexcept {
369 [[nodiscard]]
const T*
dataPtr (
int n = 0) const noexcept {
465 void getVal (T* data,
const IntVect& pos,
int N,
int numcomp)
const noexcept;
469 template <RunOn run_on AMREX_DEFAULT_RUNON>
470 requires (std::is_same_v<T,float> || std::is_same_v<T,double>)
479 template <RunOn run_on AMREX_DEFAULT_RUNON>
480 void setVal (T
const&
x,
const Box& bx,
int dcomp,
int ncomp)
noexcept;
482 template <RunOn run_on AMREX_DEFAULT_RUNON>
517 const
Box& destbox,
int destcomp,
int numcomp) noexcept;
527 int numcomp = 1) noexcept;
540 int numcomp,
void* dst) const noexcept;
545 int numcomp, const
void* src) noexcept;
550 int numcomp, const
void* src) noexcept;
587 [[nodiscard]]
Real norm (
int p,
int scomp = 0,
int numcomp = 1) const;
591 [[nodiscard]]
Real norm (const
Box& subbox,
int p,
int scomp = 0,
int numcomp = 1) const;
597 void abs (
int comp,
int numcomp=1) noexcept;
602 void abs (const
Box& subbox,
int comp = 0,
int numcomp=1) noexcept;
607 [[nodiscard]] T
min (
int comp = 0) const noexcept;
612 [[nodiscard]] T
min (const
Box& subbox,
int comp = 0) const noexcept;
617 [[nodiscard]] T
max (
int comp = 0) const noexcept;
622 [[nodiscard]] T
max (const
Box& subbox,
int comp = 0) const noexcept;
627 [[nodiscard]] std::pair<T,T>
minmax (
int comp = 0) const noexcept;
632 [[nodiscard]] std::pair<T,T>
minmax (const
Box& subbox,
int comp = 0) const noexcept;
637 [[nodiscard]] T
maxabs (
int comp = 0) const noexcept;
642 [[nodiscard]] T
maxabs (const
Box& subbox,
int comp = 0) const noexcept;
710 [[nodiscard]] T
sum (
int comp,
int numcomp = 1) const noexcept;
713 [[nodiscard]] T
sum (const
Box& subbox,
int comp,
int numcomp = 1) const noexcept;
757 int srccomp,
int destcomp,
int numcomp=1) noexcept;
777 int numcomp=1) noexcept;
784 int srccomp,
int destcomp,
int numcomp=1) noexcept;
793 int srccomp,
int destcomp,
int numcomp) noexcept;
798 int srccomp,
int destcomp,
int numcomp=1) noexcept;
806 int srccomp,
int destcomp,
int numcomp=1) noexcept;
811 const
BaseFab<T>& src1,
int comp1,
812 const
BaseFab<T>& src2,
int comp2) noexcept;
828 int numcomp=1) noexcept;
835 int srccomp,
int destcomp,
int numcomp=1) noexcept;
861 int numcomp=1) noexcept;
869 int srccomp,
int destcomp,
int numcomp=1) noexcept;
894 int numcomp=1) noexcept;
901 int srccomp,
int destcomp,
int numcomp=1) noexcept;
926 int numcomp=1) noexcept;
935 int srccomp,
int destcomp,
int numcomp=1) noexcept;
949 const
BaseFab<T>& f2, const
Box& b2,
int comp2,
951 const
Box& b,
int comp,
int numcomp = 1) noexcept;
956 const
BaseFab<T>& f2,
int comp2,
958 const
Box& b,
int comp,
int numcomp = 1) noexcept;
971 const
BaseFab<T>& f2, const
Box& b2,
int comp2,
973 int comp,
int numcomp = 1) noexcept;
978 int numcomp = 1) const noexcept;
983 int numcomp) const noexcept;
993 template <RunOn run_on AMREX_DEFAULT_RUNON>
997 template <RunOn run_on AMREX_DEFAULT_RUNON>
1000 template <RunOn run_on AMREX_DEFAULT_RUNON>
1004 template <RunOn run_on AMREX_DEFAULT_RUNON>
1007 template <RunOn run_on AMREX_DEFAULT_RUNON>
1011 template <RunOn run_on AMREX_DEFAULT_RUNON>
1015 template <RunOn run_on AMREX_DEFAULT_RUNON>
1023 template <RunOn run_on AMREX_DEFAULT_RUNON>
1027 template <RunOn run_on AMREX_DEFAULT_RUNON>
1031 template <RunOn run_on AMREX_DEFAULT_RUNON>
1034 template <RunOn run_on AMREX_DEFAULT_RUNON>
1038 template <RunOn run_on AMREX_DEFAULT_RUNON>
1045 template <RunOn run_on AMREX_DEFAULT_RUNON>
1048 template <RunOn run_on AMREX_DEFAULT_RUNON>
1052 template <RunOn run_on AMREX_DEFAULT_RUNON>
1056 template <RunOn run_on AMREX_DEFAULT_RUNON>
1059 template <RunOn run_on AMREX_DEFAULT_RUNON>
1063 template <RunOn run_on AMREX_DEFAULT_RUNON>
1070 template <RunOn run_on AMREX_DEFAULT_RUNON>
1073 template <RunOn run_on AMREX_DEFAULT_RUNON>
1077 template <RunOn run_on AMREX_DEFAULT_RUNON>
1081 template <RunOn run_on AMREX_DEFAULT_RUNON>
1084 template <RunOn run_on AMREX_DEFAULT_RUNON>
1088 template <RunOn run_on AMREX_DEFAULT_RUNON>
1095 template <RunOn run_on AMREX_DEFAULT_RUNON>
1098 template <RunOn run_on AMREX_DEFAULT_RUNON>
1102 template <RunOn run_on AMREX_DEFAULT_RUNON>
1106 template <RunOn run_on AMREX_DEFAULT_RUNON>
1109 template <RunOn run_on AMREX_DEFAULT_RUNON>
1113 template <RunOn run_on AMREX_DEFAULT_RUNON>
1120 template <RunOn run_on AMREX_DEFAULT_RUNON>
1123 template <RunOn run_on AMREX_DEFAULT_RUNON>
1127 template <RunOn run_on AMREX_DEFAULT_RUNON>
1131 template <RunOn run_on AMREX_DEFAULT_RUNON>
1154 [[nodiscard]] T
dot (const
Box& bx,
int destcomp,
int numcomp) const noexcept;
1190 return this->dptr + (this->domain.index(p)+n*this->domain.numPts());
1203 return this->dptr + (this->domain.index(p)+n*this->domain.numPts());
1211 if (this->arena()->isManaged()) {
1212#if defined(AMREX_USE_SYCL)
1217#elif defined(AMREX_USE_CUDA) && !defined(_WIN32)
1219 std::size_t s =
sizeof(T)*this->nvar*this->domain.numPts();
1220#if defined(CUDART_VERSION) && (CUDART_VERSION >= 13000)
1221 cudaMemLocation location = {};
1222 location.type = cudaMemLocationTypeHost;
1231#elif defined(AMREX_USE_HIP)
1243 if (this->arena()->isManaged()) {
1244#if defined(AMREX_USE_SYCL)
1245 std::size_t s =
sizeof(T)*this->nvar*this->domain.numPts();
1246 auto& q = Gpu::Device::streamQueue();
1247 q.submit([&] (sycl::handler& h) { h.prefetch(this->dptr, s); });
1248#elif defined(AMREX_USE_CUDA) && !defined(_WIN32)
1250 std::size_t s =
sizeof(T)*this->nvar*this->domain.numPts();
1251#if defined(CUDART_VERSION) && (CUDART_VERSION >= 13000)
1252 cudaMemLocation location = {};
1253 location.type = cudaMemLocationTypeDevice;
1263#elif defined(AMREX_USE_HIP)
1280 return this->dptr[this->domain.index(p)+n*this->domain.numPts()];
1291 return this->dptr[this->domain.index(p)];
1304 return this->dptr[this->domain.index(p)+n*this->domain.numPts()];
1315 return this->dptr[this->domain.index(p)];
1323 int numcomp)
const noexcept
1325 const int loc = this->domain.index(pos);
1326 const Long sz = this->domain.numPts();
1331 for (
int k = 0; k < numcomp; k++) {
1332 data[k] = this->dptr[loc+(n+k)*sz];
1339 const IntVect& pos)
const noexcept
1341 getVal(data,pos,0,this->nvar);
1377template <RunOn run_on>
1378requires (std::is_same_v<T,float> || std::is_same_v<T,double>)
1382 amrex::fill_snan<run_on>(this->dptr, this->truesize);
1386template <RunOn run_on>
1394template <RunOn run_on>
1402template <RunOn run_on>
1410template <RunOn run_on>
1418template <RunOn run_on>
1421 const Box& destbox,
int destcomp,
int numcomp)
noexcept
1427 AMREX_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
1428 AMREX_ASSERT(destcomp >= 0 && destcomp+numcomp <= this->nvar);
1434 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
1445template <RunOn run_on>
1453template <RunOn run_on>
1467 if (this->nvar == 0) {
return; }
1468 AMREX_ASSERT(std::numeric_limits<Long>::max()/this->nvar > this->domain.numPts());
1470 this->truesize = this->nvar*this->domain.numPts();
1471 this->ptr_owner =
true;
1472 this->dptr =
static_cast<T*
>(this->alloc(this->truesize*
sizeof(T)));
1481 if constexpr (std::is_same_v<T,float> || std::is_same_v<T,double>) {
1485 this->
template fill_snan<RunOn::Device>();
1490 this->
template fill_snan<RunOn::Host>();
1510 :
DataAllocator{ar}, domain(bx), nvar(n), shared_memory(shared)
1518 dptr(const_cast<T*>(rhs.dataPtr(scomp))),
1519 domain(rhs.domain), nvar(ncomp),
1520 truesize(ncomp*rhs.domain.numPts())
1525 this->
dptr =
nullptr;
1527 this->copy<RunOn::Device>(rhs, this->
domain, scomp, this->
domain, 0, ncomp);
1537 : dptr(p), domain(bx), nvar(ncomp), truesize(bx.numPts()*ncomp)
1543 : dptr(const_cast<T*>(p)), domain(bx), nvar(ncomp), truesize(bx.numPts()*ncomp)
1552 nvar(a.nComp()), truesize(a.size())
1560 nvar(a.nComp()), truesize(a.size())
1565 : dptr(
const_cast<T*
>(a.p)),
1568 nvar(a.nComp()), truesize(a.size())
1573 : dptr(
const_cast<T*
>(a.p)),
1576 nvar(a.nComp()), truesize(a.size())
1588 dptr(rhs.dptr), domain(rhs.domain),
1589 nvar(rhs.nvar), truesize(rhs.truesize),
1590 ptr_owner(rhs.ptr_owner), shared_memory(rhs.shared_memory)
1592 , alloc_stream(rhs.alloc_stream)
1596 rhs.ptr_owner =
false;
1605 DataAllocator::operator=(rhs);
1607 domain = rhs.domain;
1609 truesize = rhs.truesize;
1610 ptr_owner = rhs.ptr_owner;
1611 shared_memory = rhs.shared_memory;
1613 alloc_stream = rhs.alloc_stream;
1617 rhs.ptr_owner =
false;
1623template <RunOn run_on>
1638 if (ar ==
nullptr) {
1647 else if (this->dptr ==
nullptr || !this->ptr_owner)
1649 if (this->shared_memory) {
1650 amrex::Abort(
"BaseFab::resize: BaseFab in shared memory cannot increase size");
1653 this->dptr =
nullptr;
1656 else if (this->nvar*this->domain.numPts() > this->truesize
1658 || (arena()->isStreamOrderedArena() && alloc_stream !=
Gpu::gpuStream())
1662 if (this->shared_memory) {
1663 amrex::Abort(
"BaseFab::resize: BaseFab in shared memory cannot increase size");
1674requires (std::is_trivially_destructible_v<U>)
1680 o = this->ptr_owner;
1681 this->ptr_owner =
false;
1682 if (o && this->dptr) {
1683 if (this->nvar > 1) {
1692 return Elixir((o ? this->dptr :
nullptr), this->arena());
1704 if (this->ptr_owner)
1706 if (this->shared_memory)
1708 amrex::Abort(
"BaseFab::clear: BaseFab cannot be owner of shared memory");
1714 this->arena()->streamOrderedFree(this->dptr, alloc_stream);
1716 this->free(this->dptr);
1719 if (this->nvar > 1) {
1726 this->dptr =
nullptr;
1732std::unique_ptr<T,DataDeleter>
1735 std::unique_ptr<T,DataDeleter> r(
nullptr,
DataDeleter{this->arena()});
1736 if (this->dptr && this->ptr_owner) {
1737 r.reset(this->dptr);
1738 this->ptr_owner =
false;
1739 if (this->nvar > 1) {
1749template <RunOn run_on>
1754 void* dst)
const noexcept
1765 d(i,j,k,n) = s(i,j,k,n+srccomp);
1767 return sizeof(T)*d.
size();
1776template <RunOn run_on,
typename BUF>
1781 const void* src)
noexcept
1793 d(i,j,k,n+dstcomp) =
static_cast<T
>(s(i,j,k,n));
1795 return sizeof(BUF)*s.
size();
1804template <RunOn run_on,
typename BUF>
1809 const void* src)
noexcept
1821 d(i,j,k,n+dstcomp) +=
static_cast<T
>(s(i,j,k,n));
1823 return sizeof(BUF)*s.
size();
1832template <RunOn run_on>
1840template <RunOn run_on>
1844 this->abs<run_on>(this->domain,0,this->nvar);
1848template <RunOn run_on>
1852 this->abs<run_on>(this->domain,comp,numcomp);
1856template <RunOn run_on>
1863 a(i,j,k,n+comp) = std::abs(a(i,j,k,n+comp));
1868template <RunOn run_on>
1871 int scomp,
int ncomp)
const noexcept
1873 BL_ASSERT(this->domain.contains(subbox));
1874 BL_ASSERT(scomp >= 0 && scomp + ncomp <= this->nvar);
1884 reduce_op.
eval(subbox, reduce_data,
1889 for (
int n = 0; n < ncomp; ++n) {
1895 ReduceTuple hv = reduce_data.
value(reduce_op);
1896 r = amrex::get<0>(hv);
1900 amrex::LoopOnCpu(subbox, ncomp, [=,&r] (
int i,
int j,
int k,
int n)
noexcept
1903 Real t =
static_cast<Real>(std::abs(a(i,j,k,n+scomp)));
1912template <RunOn run_on>
1916 return norm<run_on>(this->domain,p,comp,numcomp);
1920template <RunOn run_on>
1924 BL_ASSERT(this->domain.contains(subbox));
1925 BL_ASSERT(comp >= 0 && comp + numcomp <= this->nvar);
1935 reduce_op.
eval(subbox, reduce_data,
1939 for (
int n = 0; n < numcomp; ++n) {
1944 ReduceTuple hv = reduce_data.
value(reduce_op);
1945 nrm = amrex::get<0>(hv);
1946 }
else if (p == 1) {
1950 reduce_op.
eval(subbox, reduce_data,
1954 for (
int n = 0; n < numcomp; ++n) {
1955 t +=
static_cast<Real>(std::abs(a(i,j,k,n+comp)));
1959 ReduceTuple hv = reduce_data.
value(reduce_op);
1960 nrm = amrex::get<0>(hv);
1961 }
else if (p == 2) {
1965 reduce_op.
eval(subbox, reduce_data,
1969 for (
int n = 0; n < numcomp; ++n) {
1970 t +=
static_cast<Real>(a(i,j,k,n+comp)*a(i,j,k,n+comp));
1974 ReduceTuple hv = reduce_data.
value(reduce_op);
1975 nrm = amrex::get<0>(hv);
1983 amrex::LoopOnCpu(subbox, numcomp, [=,&nrm] (
int i,
int j,
int k,
int n)
noexcept
1985 Real t =
static_cast<Real>(std::abs(a(i,j,k,n+comp)));
1988 }
else if (p == 1) {
1989 amrex::LoopOnCpu(subbox, numcomp, [=,&nrm] (
int i,
int j,
int k,
int n)
noexcept
1991 nrm += std::abs(a(i,j,k,n+comp));
1993 }
else if (p == 2) {
1994 amrex::LoopOnCpu(subbox, numcomp, [=,&nrm] (
int i,
int j,
int k,
int n)
noexcept
1996 nrm += a(i,j,k,n+comp)*a(i,j,k,n+comp);
2007template <RunOn run_on>
2011 return this->min<run_on>(this->domain,comp);
2015template <RunOn run_on>
2024 using ReduceTuple =
typename decltype(reduce_data)::Type;
2025 reduce_op.
eval(subbox, reduce_data,
2028 return { a(i,j,k) };
2030 ReduceTuple hv = reduce_data.
value(reduce_op);
2031 return amrex::get<0>(hv);
2035 T r = std::numeric_limits<T>::max();
2045template <RunOn run_on>
2049 return this->max<run_on>(this->domain,comp);
2053template <RunOn run_on>
2062 using ReduceTuple =
typename decltype(reduce_data)::Type;
2063 reduce_op.
eval(subbox, reduce_data,
2066 return { a(i,j,k) };
2068 ReduceTuple hv = reduce_data.
value(reduce_op);
2069 return amrex::get<0>(hv);
2073 T r = std::numeric_limits<T>::lowest();
2083template <RunOn run_on>
2087 return this->minmax<run_on>(this->domain,comp);
2091template <RunOn run_on>
2100 using ReduceTuple =
typename decltype(reduce_data)::Type;
2101 reduce_op.
eval(subbox, reduce_data,
2104 auto const x = a(i,j,k);
2107 ReduceTuple hv = reduce_data.
value(reduce_op);
2108 return std::make_pair(amrex::get<0>(hv), amrex::get<1>(hv));
2112 T rmax = std::numeric_limits<T>::lowest();
2113 T rmin = std::numeric_limits<T>::max();
2116 auto const x = a(i,j,k);
2120 return std::make_pair(rmin,rmax);
2125template <RunOn run_on>
2129 return this->maxabs<run_on>(this->domain,comp);
2133template <RunOn run_on>
2142 using ReduceTuple =
typename decltype(reduce_data)::Type;
2143 reduce_op.
eval(subbox, reduce_data,
2146 return { std::abs(a(i,j,k)) };
2148 ReduceTuple hv = reduce_data.
value(reduce_op);
2149 return amrex::get<0>(hv);
2164template <RunOn run_on>
2172 std::numeric_limits<int>::lowest(),
2173 std::numeric_limits<int>::lowest())};
2180 if ((*flag == 0) && (a(i,j,k) == value)) {
2203template <RunOn run_on>
2207 return this->minIndex<run_on>(this->domain,comp);
2211template <RunOn run_on>
2215 T min_val = this->min<run_on>(subbox, comp);
2216 return this->indexFromValue<run_on>(subbox, comp, min_val);
2220template <RunOn run_on>
2224 min_val = this->min<run_on>(subbox, comp);
2225 min_idx = this->indexFromValue<run_on>(subbox, comp, min_val);
2229template <RunOn run_on>
2233 return this->maxIndex<run_on>(this->domain,comp);
2237template <RunOn run_on>
2241 T max_val = this->max<run_on>(subbox, comp);
2242 return this->indexFromValue<run_on>(subbox, comp, max_val);
2246template <RunOn run_on>
2250 max_val = this->max<run_on>(subbox, comp);
2251 max_idx = this->indexFromValue<run_on>(subbox, comp, max_val);
2255template <RunOn run_on>
2259 mask.resize(this->domain,1);
2267 using ReduceTuple =
typename decltype(reduce_data)::Type;
2268 reduce_op.
eval(this->domain, reduce_data,
2272 if (a(i,j,k) < val) {
2281 ReduceTuple hv = reduce_data.
value(reduce_op);
2282 cnt = amrex::get<0>(hv);
2288 if (a(i,j,k) < val) {
2301template <RunOn run_on>
2305 mask.resize(this->domain,1);
2313 using ReduceTuple =
typename decltype(reduce_data)::Type;
2314 reduce_op.
eval(this->domain, reduce_data,
2318 if (a(i,j,k) <= val) {
2327 ReduceTuple hv = reduce_data.
value(reduce_op);
2328 cnt = amrex::get<0>(hv);
2334 if (a(i,j,k) <= val) {
2347template <RunOn run_on>
2351 mask.resize(this->domain,1);
2359 using ReduceTuple =
typename decltype(reduce_data)::Type;
2360 reduce_op.
eval(this->domain, reduce_data,
2364 if (a(i,j,k) == val) {
2373 ReduceTuple hv = reduce_data.
value(reduce_op);
2374 cnt = amrex::get<0>(hv);
2380 if (a(i,j,k) == val) {
2393template <RunOn run_on>
2397 mask.resize(this->domain,1);
2405 using ReduceTuple =
typename decltype(reduce_data)::Type;
2406 reduce_op.
eval(this->domain, reduce_data,
2410 if (a(i,j,k) > val) {
2419 ReduceTuple hv = reduce_data.
value(reduce_op);
2420 cnt = amrex::get<0>(hv);
2426 if (a(i,j,k) > val) {
2439template <RunOn run_on>
2443 mask.resize(this->domain,1);
2451 using ReduceTuple =
typename decltype(reduce_data)::Type;
2452 reduce_op.
eval(this->domain, reduce_data,
2456 if (a(i,j,k) >= val) {
2465 ReduceTuple hv = reduce_data.
value(reduce_op);
2466 cnt = amrex::get<0>(hv);
2472 if (a(i,j,k) >= val) {
2485template <RunOn run_on>
2489 Box ovlp(this->domain);
2491 return ovlp.
ok() ? this->atomicAdd<run_on>(
x,ovlp,ovlp,0,0,this->nvar) : *
this;
2495template <RunOn run_on>
2498 int srccomp,
int destcomp,
int numcomp)
noexcept
2505 BL_ASSERT( srccomp >= 0 && srccomp+numcomp <=
x.nComp());
2512 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
2522template <RunOn run_on>
2526 Box ovlp(this->domain);
2528 return ovlp.
ok() ? saxpy<run_on>(a,
x,ovlp,ovlp,0,0,this->nvar) : *
this;
2532template <RunOn run_on>
2535 const Box& srcbox,
const Box& destbox,
2536 int srccomp,
int destcomp,
int numcomp)
noexcept
2543 BL_ASSERT( srccomp >= 0 && srccomp+numcomp <=
x.nComp());
2550 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
2553 d(i,j,k,n+destcomp) = s(i+
offset.x,j+
offset.y,k+
offset.z,n+srccomp) + a*d(i,j,k,n+destcomp);
2560template <RunOn run_on>
2568 BL_ASSERT( comp1 >= 0 && comp1+numcomp <= src1.nComp());
2569 BL_ASSERT( comp2 >= 0 && comp2+numcomp <= src2.nComp());
2577 d(i,j,k,n+destcomp) += s1(i,j,k,n+comp1) * s2(i,j,k,n+comp2);
2584template <RunOn run_on>
2589 int comp,
int numcomp)
noexcept
2599 BL_ASSERT(comp1 >= 0 && comp1+numcomp <= f1.nComp());
2600 BL_ASSERT(comp2 >= 0 && comp2+numcomp <= f2.nComp());
2609 const Dim3 off1{.
x = slo1.x-dlo.x, .y = slo1.y-dlo.y, .z = slo1.z-dlo.z};
2610 const Dim3 off2{.
x = slo2.x-dlo.x, .y = slo2.y-dlo.y, .z = slo2.z-dlo.z};
2614 d(i,j,k,n+comp) = alpha*s1(i+off1.x,j+off1.y,k+off1.z,n+comp1)
2615 + beta*s2(i+off2.x,j+off2.y,k+off2.z,n+comp2);
2621template <RunOn run_on>
2625 int numcomp)
const noexcept
2632 BL_ASSERT(ycomp >= 0 && ycomp+numcomp <=
y.nComp());
2638 const Dim3 offset{.
x = ylo.x-xlo.x, .y = ylo.y-xlo.y, .z = ylo.z-xlo.z};
2646 using ReduceTuple =
typename decltype(reduce_data)::Type;
2647 reduce_op.
eval(xbx, reduce_data,
2651 for (
int n = 0; n < numcomp; ++n) {
2656 ReduceTuple hv = reduce_data.
value(reduce_op);
2657 r = amrex::get<0>(hv);
2671template <RunOn run_on>
2675 int numcomp)
const noexcept
2682 BL_ASSERT(ycomp >= 0 && ycomp+numcomp <=
y.nComp());
2688 const Dim3 offset{.
x = ylo.x-xlo.x, .y = ylo.y-xlo.y, .z = ylo.z-xlo.z};
2698 using ReduceTuple =
typename decltype(reduce_data)::Type;
2699 reduce_op.
eval(xbx, reduce_data,
2702 int m =
static_cast<int>(
static_cast<bool>(ma(i,j,k)));
2704 for (
int n = 0; n < numcomp; ++n) {
2709 ReduceTuple hv = reduce_data.
value(reduce_op);
2710 r = amrex::get<0>(hv);
2716 int m =
static_cast<int>(
static_cast<bool>(ma(i,j,k)));
2725template <RunOn run_on>
2733template <RunOn run_on>
2741template <RunOn run_on>
2745 return this->negate<run_on>(this->domain,
DestComp{comp},
NumComps{numcomp});
2749template <RunOn run_on>
2757template <RunOn run_on>
2761 return this->invert<run_on>(r, this->domain,
DestComp{comp},
NumComps{numcomp});
2765template <RunOn run_on>
2773template <RunOn run_on>
2777 return this->plus<run_on>(r, this->domain,
DestComp{comp},
NumComps{numcomp});
2781template <RunOn run_on>
2789template <RunOn run_on>
2797template <RunOn run_on>
2801 Box ovlp(this->domain);
2803 return ovlp.
ok() ? this->atomicAdd<run_on>(src,ovlp,ovlp,srccomp,destcomp,numcomp) : *
this;
2807template <RunOn run_on>
2810 int numcomp)
noexcept
2816template <RunOn run_on>
2819 int numcomp)
noexcept
2821 Box ovlp(this->domain);
2824 return ovlp.
ok() ? this->atomicAdd<run_on>(src,ovlp,ovlp,srccomp,destcomp,numcomp) : *
this;
2828template <RunOn run_on>
2831 int srccomp,
int destcomp,
int numcomp)
noexcept
2837 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
2844 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
2856template <RunOn run_on,
typename T>
2859 const Box& srcbox,
const Box& destbox,
2860 int srccomp,
int destcomp,
int numcomp)
noexcept
2866 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
2869 T* p = d.
ptr(i,j,k,n+destcomp);
2874template <RunOn run_on,
typename T>
2875requires (!HasAtomicAdd<T>::value)
2876void basefab_atomic_add (BaseFab<T>& dfab,
const BaseFab<T>& sfab,
2877 const Box& srcbox,
const Box& destbox,
2878 int srccomp,
int destcomp,
int numcomp)
2887template <RunOn run_on>
2890 int srccomp,
int destcomp,
int numcomp)
noexcept
2896 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
2899 detail::basefab_atomic_add<run_on>(*
this, src, srcbox, destbox,
2900 srccomp, destcomp, numcomp);
2906template <RunOn run_on>
2909 int srccomp,
int destcomp,
int numcomp)
noexcept
2911#if defined(AMREX_USE_OMP) && (AMREX_SPACEDIM > 1)
2912#if defined(AMREX_USE_GPU)
2919 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
2928 Dim3 const offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
2944 for (
int ip = 0; ip < nplanes; ++ip) {
2949 int planes_left = nplanes;
2950 while (planes_left > 0) {
2952 auto const m = mm + plo;
2953 auto* lock = OpenMP::get_lock(m);
2954 if (omp_test_lock(lock))
2958 if (planedim == 1) {
2966 for (
int n = 0; n < numcomp; ++n) {
2967 for (
int k = lo.z; k <= hi.z; ++k) {
2968 for (
int j = lo.y; j <= hi.y; ++j) {
2969 auto *
pdst = d.
ptr(dlo.x,j ,k ,n+destcomp);
2972 for (
int ii = 0; ii < len.x; ++ii) {
2973 pdst[ii] += psrc[ii];
2981 omp_unset_lock(lock);
2982 if (planes_left == 0) {
break; }
2986 for (
int ip = 0; ip < nplanes; ++ip) {
2987 int new_mm = (mm+ip) % nplanes;
2988 if ( !
mask[new_mm] ) {
2999#if defined(AMREX_USE_GPU)
3001 return this->
template atomicAdd<run_on>(src, srcbox, destbox, srccomp, destcomp, numcomp);
3005 return this->
template atomicAdd<run_on>(src, srcbox, destbox, srccomp, destcomp, numcomp);
3010template <RunOn run_on>
3018template <RunOn run_on>
3026template <RunOn run_on>
3029 int srccomp,
int destcomp,
int numcomp)
noexcept
3035 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
3042 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
3052template <RunOn run_on>
3056 return this->mult<run_on>(r, this->domain,
DestComp{comp},
NumComps{numcomp});
3060template <RunOn run_on>
3068template <RunOn run_on>
3076template <RunOn run_on>
3084template <RunOn run_on>
3087 int srccomp,
int destcomp,
int numcomp)
noexcept
3093 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
3100 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
3110template <RunOn run_on>
3114 return this->divide<run_on>(r, this->domain,
DestComp{comp},
NumComps{numcomp});
3118template <RunOn run_on>
3126template <RunOn run_on>
3134template <RunOn run_on>
3142template <RunOn run_on>
3145 int srccomp,
int destcomp,
int numcomp)
noexcept
3151 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
3158 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
3168template <RunOn run_on>
3172 Box ovlp(this->domain);
3174 return ovlp.
ok() ? this->protected_divide<run_on>(src,ovlp,ovlp,0,0,this->nvar) : *
this;
3178template <RunOn run_on>
3182 Box ovlp(this->domain);
3184 return ovlp.
ok() ? this->protected_divide<run_on>(src,ovlp,ovlp,srccomp,destcomp,numcomp) : *
this;
3188template <RunOn run_on>
3191 int numcomp)
noexcept
3193 Box ovlp(this->domain);
3196 return ovlp.
ok() ? this->protected_divide<run_on>(src,ovlp,ovlp,srccomp,destcomp,numcomp) : *
this;
3200template <RunOn run_on>
3203 int srccomp,
int destcomp,
int numcomp)
noexcept
3209 BL_ASSERT(srccomp >= 0 && srccomp+numcomp <= src.nComp());
3216 const Dim3 offset{.
x = slo.x-dlo.x, .y = slo.y-dlo.y, .z = slo.z-dlo.z};
3220 d(i,j,k,n+destcomp) /= s(i+offset.x,j+offset.y,k+offset.z,n+srccomp);
3238template <RunOn run_on>
3243 const Box& b,
int comp,
int numcomp)
noexcept
3246 return copy<run_on>(f1,b1,comp1,b,comp,numcomp);
3248 return copy<run_on>(f2,b2,comp2,b,comp,numcomp);
3250 Real alpha = (t2-t)/(t2-t1);
3251 Real beta = (t-t1)/(t2-t1);
3252 return linComb<run_on>(f1,b1,comp1,f2,b2,comp2,alpha,beta,b,comp,numcomp);
3257template <RunOn run_on>
3262 const Box& b,
int comp,
int numcomp)
noexcept
3265 return copy<run_on>(f1,b,comp1,b,comp,numcomp);
3267 return copy<run_on>(f2,b,comp2,b,comp,numcomp);
3269 Real alpha = (t2-t)/(t2-t1);
3270 Real beta = (t-t1)/(t2-t1);
3271 return linComb<run_on>(f1,b,comp1,f2,b,comp2,alpha,beta,b,comp,numcomp);
3280template <RunOn run_on>
3284 this->setVal<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3288template <RunOn run_on>
3292 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3296 a(i,j,k,n+dcomp.i) = x;
3301template <RunOn run_on>
3309template <RunOn run_on>
3313 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3318 if (m(i,j,k)) { a(i,j,k,n+dcomp.i) = val; }
3323template <RunOn run_on>
3331template <RunOn run_on>
3335 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3340 if (!m(i,j,k)) { a(i,j,k,n+dcomp.i) = val; }
3345template <RunOn run_on>
3350 if (run_on == RunOn::Device && Gpu::inLaunchRegion()) {
3355 for (
int n = dcomp.i; n < ncomp.n+dcomp.i; ++n) {
3364 for (
auto const& b : b_lst) {
3365 this->setVal<RunOn::Host>(
x, b, dcomp, ncomp);
3371template <RunOn run_on>
3380template <RunOn run_on>
3386 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3387 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3395 d(i,j,k,n+dcomp.i) = s(i,j,k,n+scomp.i);
3402template <RunOn run_on>
3406 return this->plus<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3410template <RunOn run_on>
3414 return this->plus<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3418template <RunOn run_on>
3422 BL_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3427 a(i,j,k,n+dcomp.i) += val;
3434template <RunOn run_on>
3442template <RunOn run_on>
3450template <RunOn run_on>
3456 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3457 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3465 d(i,j,k,n+dcomp.i) += s(i,j,k,n+scomp.i);
3472template <RunOn run_on>
3476 return this->minus<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3480template <RunOn run_on>
3484 return this->minus<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3488template <RunOn run_on>
3492 BL_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3497 a(i,j,k,n+dcomp.i) -= val;
3504template <RunOn run_on>
3512template <RunOn run_on>
3520template <RunOn run_on>
3526 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3527 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3535 d(i,j,k,n+dcomp.i) -= s(i,j,k,n+scomp.i);
3542template <RunOn run_on>
3546 return this->mult<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3550template <RunOn run_on>
3554 return this->mult<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3558template <RunOn run_on>
3562 BL_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3567 a(i,j,k,n+dcomp.i) *= val;
3574template <RunOn run_on>
3582template <RunOn run_on>
3590template <RunOn run_on>
3596 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3597 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3605 d(i,j,k,n+dcomp.i) *= s(i,j,k,n+scomp.i);
3612template <RunOn run_on>
3616 return this->divide<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3620template <RunOn run_on>
3624 return this->divide<run_on>(val, this->domain,
DestComp{0},
NumComps{this->nvar});
3628template <RunOn run_on>
3632 BL_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3637 a(i,j,k,n+dcomp.i) /= val;
3644template <RunOn run_on>
3652template <RunOn run_on>
3660template <RunOn run_on>
3666 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3667 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3675 d(i,j,k,n+dcomp.i) /= s(i,j,k,n+scomp.i);
3682template <RunOn run_on>
3686 return this->negate<run_on>(this->domain,
DestComp{0},
NumComps{this->nvar});
3690template <RunOn run_on>
3694 BL_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3699 a(i,j,k,n+dcomp.i) = -a(i,j,k,n+dcomp.i);
3706template <RunOn run_on>
3710 return this->invert<run_on>(r, this->domain,
DestComp{0},
NumComps{this->nvar});
3714template <RunOn run_on>
3718 BL_ASSERT(dcomp.i >= 0 && dcomp.i + ncomp.n <= this->nvar);
3723 a(i,j,k,n+dcomp.i) = r / a(i,j,k,n+dcomp.i);
3730template <RunOn run_on>
3734 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3739 if (run_on == RunOn::Device && Gpu::inLaunchRegion()) {
3742 using ReduceTuple =
typename decltype(reduce_data)::Type;
3743 reduce_op.
eval(bx, reduce_data,
3747 for (
int n = 0; n < ncomp.n; ++n) {
3748 t += a(i,j,k,n+dcomp.i);
3752 ReduceTuple hv = reduce_data.
value(reduce_op);
3753 r = amrex::get<0>(hv);
3759 r += a(i,j,k,n+dcomp.i);
3767template <RunOn run_on>
3772 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3773 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3779 if (run_on == RunOn::Device && Gpu::inLaunchRegion()) {
3782 using ReduceTuple =
typename decltype(reduce_data)::Type;
3783 reduce_op.
eval(bx, reduce_data,
3787 for (
int n = 0; n < ncomp.n; ++n) {
3788 t += d(i,j,k,n+dcomp.i) * s(i,j,k,n+scomp.i);
3792 ReduceTuple hv = reduce_data.
value(reduce_op);
3793 r = amrex::get<0>(hv);
3799 r += d(i,j,k,n+dcomp.i) * s(i,j,k,n+scomp.i);
3807template <RunOn run_on>
3816template <RunOn run_on>
3820 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3825 if (run_on == RunOn::Device && Gpu::inLaunchRegion()) {
3828 using ReduceTuple =
typename decltype(reduce_data)::Type;
3829 reduce_op.
eval(bx, reduce_data,
3833 for (
int n = 0; n < ncomp.n; ++n) {
3834 t += a(i,j,k,n+dcomp.i)*a(i,j,k,n+dcomp.i);
3838 ReduceTuple hv = reduce_data.
value(reduce_op);
3839 r = amrex::get<0>(hv);
3845 r += a(i,j,k,n+dcomp.i)*a(i,j,k,n+dcomp.i);
3853template <RunOn run_on>
3860 AMREX_ASSERT(scomp.i >= 0 && scomp.i+ncomp.n <= src.nvar);
3861 AMREX_ASSERT(dcomp.i >= 0 && dcomp.i+ncomp.n <= this->nvar);
3868 if (run_on == RunOn::Device && Gpu::inLaunchRegion()) {
3871 using ReduceTuple =
typename decltype(reduce_data)::Type;
3872 reduce_op.
eval(bx, reduce_data,
3876 T mi =
static_cast<T
>(
static_cast<int>(
static_cast<bool>(m(i,j,k))));
3877 for (
int n = 0; n < ncomp.n; ++n) {
3878 t += d(i,j,k,n+dcomp.i)*s(i,j,k,n+scomp.i)*mi;
3882 ReduceTuple hv = reduce_data.
value(reduce_op);
3883 r = amrex::get<0>(hv);
3889 int mi = static_cast<int>(static_cast<bool>(m(i,j,k)));
3890 r += d(i,j,k,n+dcomp.i)*s(i,j,k,n+scomp.i)*mi;
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_DEFAULT_RUNON
Definition AMReX_GpuControl.H:69
#define AMREX_CUDA_SAFE_CALL(call)
Definition AMReX_GpuError.H:73
#define AMREX_HOST_DEVICE_FOR_1D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:105
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D_FLAG(where_to_run, box, nc, i, j, k, n, block)
Definition AMReX_GpuLaunch.nolint.H:75
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
int idir
Definition AMReX_HypreMLABecLap.cpp:1143
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
Real * pdst
Definition AMReX_HypreMLABecLap.cpp:1140
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
#define AMREX_LOOP_3D(bx, i, j, k, block)
Definition AMReX_Loop.nolint.H:4
#define AMREX_LOOP_4D(bx, ncomp, i, j, k, n, block)
Definition AMReX_Loop.nolint.H:16
void amrex_mempool_free(void *p)
Definition AMReX_MemPool.cpp:80
void * amrex_mempool_alloc(size_t nbytes)
Definition AMReX_MemPool.cpp:74
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A virtual base class for objects that manage their own dynamic memory allocation.
Definition AMReX_Arena.H:132
A FortranArrayBox(FAB)-like object.
Definition AMReX_BaseFab.H:194
int maskGE(BaseFab< int > &mask, T const &val, int comp=0) const noexcept
Same as above except mark cells with value greater than or equal to val.
Definition AMReX_BaseFab.H:2441
BaseFab< T > & saxpy(T a, const BaseFab< T > &x, const Box &srcbox, const Box &destbox, int srccomp, int destcomp, int numcomp=1) noexcept
FAB SAXPY (y[i] <- y[i] + a * x[i]), in place.
Definition AMReX_BaseFab.H:2497
Array4< T const > const_array() const noexcept
Definition AMReX_BaseFab.H:423
BaseFab< T > & copy(const BaseFab< T > &src, Box bx, SrcComp scomp, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx does not intersect with src fab.
Definition AMReX_BaseFab.H:3382
T sum(int comp, int numcomp=1) const noexcept
Returns sum of given component of FAB state vector.
Definition AMReX_BaseFab.H:2727
gpuStream_t alloc_stream
Definition AMReX_BaseFab.H:1176
Real norminfmask(const Box &subbox, const BaseFab< int > &mask, int scomp=0, int ncomp=1) const noexcept
Definition AMReX_BaseFab.H:1870
BaseFab< T > & divide(T const &val) noexcept
Scalar division on the whole domain and all components.
Definition AMReX_BaseFab.H:3614
const int * hiVect() const noexcept
Returns the upper corner of the domain.
Definition AMReX_BaseFab.H:334
BaseFab< T > & minus(const BaseFab< T > &src, Box bx, SrcComp scomp, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx does not intersect with src fab.
Definition AMReX_BaseFab.H:3522
BaseFab< T > & lockAdd(const BaseFab< T > &src, const Box &srcbox, const Box &destbox, int srccomp, int destcomp, int numcomp) noexcept
Atomically add srcbox region of src FAB to destbox region of this FAB. The srcbox and destbox must be...
Definition AMReX_BaseFab.H:2908
std::size_t copyToMem(const Box &srcbox, int srccomp, int numcomp, void *dst) const noexcept
Copy from the srcbox of this Fab to raw memory and return the number of bytes copied.
Definition AMReX_BaseFab.H:1751
std::size_t addFromMem(const Box &dstbox, int dstcomp, int numcomp, const void *src) noexcept
Add from raw memory to the dstbox of this Fab and return the number of bytes copied.
Definition AMReX_BaseFab.H:1806
std::size_t nBytesOwned() const noexcept
Definition AMReX_BaseFab.H:276
BaseFab< T > & copy(const BaseFab< T > &src) noexcept
Definition AMReX_BaseFab.H:3373
BaseFab< T > & addproduct(const Box &destbox, int destcomp, int numcomp, const BaseFab< T > &src1, int comp1, const BaseFab< T > &src2, int comp2) noexcept
y[i] <- y[i] + x1[i] * x2[i])
Definition AMReX_BaseFab.H:2562
BaseFab< T > & minus(T const &val) noexcept
Scalar subtraction on the whole domain and all components.
Definition AMReX_BaseFab.H:3474
int maskLT(BaseFab< int > &mask, T const &val, int comp=0) const noexcept
Compute mask array with value of 1 in cells where BaseFab has value less than val,...
Definition AMReX_BaseFab.H:2257
BaseFab< T > & plus(T const &val) noexcept
Scalar addition on the whole domain and all components.
Definition AMReX_BaseFab.H:3404
BaseFab< T > & mult(T const &val, Box const &bx, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3560
BaseFab< T > & mult(const BaseFab< T > &src, Box bx, SrcComp scomp, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx does not intersect with src fab.
Definition AMReX_BaseFab.H:3592
std::size_t nBytes(const Box &bx, int ncomps) const noexcept
Returns bytes used in the Box for those components.
Definition AMReX_BaseFab.H:281
void setPtr(T *p, Long sz) noexcept
Definition AMReX_BaseFab.H:381
BaseFab< T > & linComb(const BaseFab< T > &f1, const Box &b1, int comp1, const BaseFab< T > &f2, const Box &b2, int comp2, Real alpha, Real beta, const Box &b, int comp, int numcomp=1) noexcept
Linear combination. Result is alpha*f1 + beta*f2. Data is taken from b1 region of f1,...
Definition AMReX_BaseFab.H:2586
void define()
Allocates memory for the BaseFab<T>.
Definition AMReX_BaseFab.H:1462
BaseFab< T > & operator*=(T const &val) noexcept
Definition AMReX_BaseFab.H:3552
void resize(const Box &b, int N=1, Arena *ar=nullptr)
This function resizes a BaseFab so it covers the Box b with N components.
Definition AMReX_BaseFab.H:1633
void clear()
The function returns the BaseFab to the invalid state. The memory is freed.
Definition AMReX_BaseFab.H:1697
const IntVect & smallEnd() const noexcept
Returns the lower corner of the domain See class Box for analogue.
Definition AMReX_BaseFab.H:311
BaseFab< T > & mult(T const &r, int comp, int numcomp=1) noexcept
Scalar multiplication, except control which components are multiplied.
Definition AMReX_BaseFab.H:3054
BaseFab< T > & atomicAdd(const BaseFab< T > &x) noexcept
Atomic FAB addition (a[i] <- a[i] + b[i]).
Definition AMReX_BaseFab.H:2487
int maskEQ(BaseFab< int > &mask, T const &val, int comp=0) const noexcept
Same as above except mark cells with value equal to val.
Definition AMReX_BaseFab.H:2349
const int * loVect() const noexcept
Returns the lower corner of the domain.
Definition AMReX_BaseFab.H:324
bool contains(const BaseFab< T > &fab) const noexcept
Returns true if the domain of fab is totally contained within the domain of this BaseFab.
Definition AMReX_BaseFab.H:340
bool isAllocated() const noexcept
Returns true if the data for the FAB has been allocated.
Definition AMReX_BaseFab.H:441
std::unique_ptr< T, DataDeleter > release() noexcept
Release ownership of memory.
Definition AMReX_BaseFab.H:1733
void setVal(T const &x, Box const &bx, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3290
BaseFab< T > & operator-=(T const &val) noexcept
Definition AMReX_BaseFab.H:3482
const Box & box() const noexcept
Returns the domain (box) where the array is defined.
Definition AMReX_BaseFab.H:299
void setValIf(T const &val, Box const &bx, const BaseFab< int > &mask, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3311
Array4< T > array() noexcept
Definition AMReX_BaseFab.H:405
IntVect indexFromValue(const Box &subbox, int comp, T const &value) const noexcept
Definition AMReX_BaseFab.H:2166
BaseFab< T > & mult(const BaseFab< T > &src) noexcept
Definition AMReX_BaseFab.H:3576
bool shared_memory
Is the memory allocated in shared memory?
Definition AMReX_BaseFab.H:1174
int maskLE(BaseFab< int > &mask, T const &val, int comp=0) const noexcept
Same as above except mark cells with value less than or equal to val.
Definition AMReX_BaseFab.H:2303
void setValIf(T const &val, const BaseFab< int > &mask) noexcept
Definition AMReX_BaseFab.H:3303
BaseFab< T > & plus(const BaseFab< T > &src, Box bx, SrcComp scomp, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx does not intersect with src fab.
Definition AMReX_BaseFab.H:3452
void setValIfNot(T const &val, const BaseFab< int > &mask) noexcept
Definition AMReX_BaseFab.H:3325
BaseFab< T > & xpay(T a, const BaseFab< T > &x, const Box &srcbox, const Box &destbox, int srccomp, int destcomp, int numcomp=1) noexcept
FAB XPAY (y[i] <- x[i] + a * y[i])
Definition AMReX_BaseFab.H:2534
T & operator()(const IntVect &p, int N) noexcept
Returns a reference to the Nth component value defined at position p in the domain....
Definition AMReX_BaseFab.H:1273
std::size_t nBytes() const noexcept
Returns how many bytes used.
Definition AMReX_BaseFab.H:274
std::size_t copyFromMem(const Box &dstbox, int dstcomp, int numcomp, const void *src) noexcept
Copy from raw memory to the dstbox of this Fab and return the number of bytes copied.
Definition AMReX_BaseFab.H:1778
BaseFab< T > & negate() noexcept
on the whole domain and all components
Definition AMReX_BaseFab.H:3684
BaseFab< T > & minus(const BaseFab< T > &src) noexcept
Definition AMReX_BaseFab.H:3506
BaseFab< T > & minus(const BaseFab< T > &src, int srccomp, int destcomp, int numcomp=1) noexcept
Subtract src components (srccomp:srccomp+numcomp-1) to this FABs components (destcomp:destcomp+numcom...
Definition AMReX_BaseFab.H:3012
T value_type
Definition AMReX_BaseFab.H:199
void SetBoxType(const IndexType &typ) noexcept
Change the Box type without change the length.
Definition AMReX_BaseFab.H:986
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:387
T maxabs(int comp=0) const noexcept
Definition AMReX_BaseFab.H:2127
Elixir elixir() noexcept
Definition AMReX_BaseFab.H:1676
BaseFab< T > & operator+=(T const &val) noexcept
Definition AMReX_BaseFab.H:3412
void setComplement(T const &x, Box const &bx, DestComp dcomp, NumComps ncomp) noexcept
setVal on the complement of bx in the fab's domain
Definition AMReX_BaseFab.H:3347
BaseFab< T > & minus(T const &val, Box const &bx, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3490
Long truesize
nvar*numpts that was allocated on heap.
Definition AMReX_BaseFab.H:1172
void setVal(T const &val) noexcept
Set value on the whole domain and all components.
Definition AMReX_BaseFab.H:3282
const int * nCompPtr() const noexcept
for calls to fortran.
Definition AMReX_BaseFab.H:288
Array4< T const > const_array(int start_comp, int num_comps) const noexcept
Definition AMReX_BaseFab.H:435
Box domain
My index space.
Definition AMReX_BaseFab.H:1170
void fill_snan() noexcept
Definition AMReX_BaseFab.H:1380
bool contains(const Box &bx) const noexcept
Returns true if bx is totally contained within the domain of this BaseFab.
Definition AMReX_BaseFab.H:349
T * dptr
The data pointer.
Definition AMReX_BaseFab.H:1169
BaseFab< T > & shift(const IntVect &v) noexcept
Perform shifts upon the domain of the BaseFab. They are completely analogous to the corresponding Box...
Definition AMReX_BaseFab.H:1346
BaseFab< T > & divide(T const &val, Box const &bx, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3630
T max(int comp=0) const noexcept
Definition AMReX_BaseFab.H:2047
BaseFab< T > & copy(const BaseFab< T > &src, const Box &srcbox, int srccomp, const Box &destbox, int destcomp, int numcomp) noexcept
The copy functions copy the contents of one BaseFab into another. The destination BaseFab is always t...
Definition AMReX_BaseFab.H:1420
Array4< T > array(int start_comp, int num_comps) noexcept
Definition AMReX_BaseFab.H:417
int nvar
Number components.
Definition AMReX_BaseFab.H:1171
T dot(const Box &xbx, int xcomp, const BaseFab< T > &y, const Box &ybx, int ycomp, int numcomp=1) const noexcept
Dot product of x (i.e.,this) and y.
Definition AMReX_BaseFab.H:2623
BaseFab< T > & divide(T const &r, int comp, int numcomp=1) noexcept
As above except specify which components.
Definition AMReX_BaseFab.H:3112
Array4< T const > array(int start_comp) const noexcept
Definition AMReX_BaseFab.H:393
BaseFab< T > & operator/=(T const &val) noexcept
Definition AMReX_BaseFab.H:3622
std::pair< T, T > minmax(int comp=0) const noexcept
Definition AMReX_BaseFab.H:2085
Array4< T const > const_array(int start_comp) const noexcept
Definition AMReX_BaseFab.H:429
void setVal(T const &x, const Box &bx, int dcomp, int ncomp) noexcept
The setVal functions set sub-regions in the BaseFab to a constant value. This most general form speci...
Definition AMReX_BaseFab.H:1404
Long size() const noexcept
Returns the total number of points of all components.
Definition AMReX_BaseFab.H:296
BaseFab< T > & plus(T const &r, const Box &b, int comp=0, int numcomp=1) noexcept
Scalar addition (a[i] <- a[i] + r), most general.
Definition AMReX_BaseFab.H:2783
BaseFab< T > & operator=(const BaseFab< T > &rhs)=delete
void getVal(T *data, const IntVect &pos, int N, int numcomp) const noexcept
This function puts numcomp component values, starting at component N, from position pos in the domain...
Definition AMReX_BaseFab.H:1320
const IntVect & bigEnd() const noexcept
Returns the upper corner of the domain. See class Box for analogue.
Definition AMReX_BaseFab.H:314
Array4< T > array(int start_comp) noexcept
Definition AMReX_BaseFab.H:411
Long numPts() const noexcept
Returns the number of points.
Definition AMReX_BaseFab.H:293
const T * dataPtr(int n=0) const noexcept
Same as above except works on const FABs.
Definition AMReX_BaseFab.H:369
void setValIfNot(T const &val, Box const &bx, const BaseFab< int > &mask, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3333
BaseFab< T > & mult(T const &val) noexcept
Scalar multiplication on the whole domain and all components.
Definition AMReX_BaseFab.H:3544
BaseFab< T > & divide(const BaseFab< T > &src, Box bx, SrcComp scomp, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx does not intersect with src fab.
Definition AMReX_BaseFab.H:3662
void setValIfNot(T const &val, const Box &bx, const BaseFab< int > &mask, int nstart, int num) noexcept
Definition AMReX_BaseFab.H:1412
BaseFab< T > & shiftHalf(int dir, int n_cell) noexcept
Perform shifts upon the domain of the BaseFab. They are completely analogous to the corresponding Box...
Definition AMReX_BaseFab.H:1370
void setVal(T const &x, const Box &bx, int N=0) noexcept
Same as above, except the number of modified components is one. N is the component to be modified.
Definition AMReX_BaseFab.H:1388
void prefetchToDevice() const noexcept
Definition AMReX_BaseFab.H:1240
T * dataPtr(int n=0) noexcept
Returns a pointer to an object of type T that is the value of the Nth component associated with the c...
Definition AMReX_BaseFab.H:360
bool ptr_owner
Owner of T*?
Definition AMReX_BaseFab.H:1173
virtual ~BaseFab() noexcept
The destructor deletes the array memory.
Definition AMReX_BaseFab.H:1580
IntVect length() const noexcept
Returns a pointer to an array of SPACEDIM integers giving the length of the domain in each direction.
Definition AMReX_BaseFab.H:305
IntVect maxIndex(int comp=0) const noexcept
Definition AMReX_BaseFab.H:2231
BaseFab< T > & protected_divide(const BaseFab< T > &src) noexcept
Divide wherever "src" is "true" or "non-zero".
Definition AMReX_BaseFab.H:3170
friend class BaseFab
Definition AMReX_BaseFab.H:197
BaseFab< T > & invert(T const &r, const Box &b, int comp=0, int numcomp=1) noexcept
Most general version, specify subbox and which components.
Definition AMReX_BaseFab.H:2767
Array4< T const > array(int start_comp, int num_comps) const noexcept
Definition AMReX_BaseFab.H:399
int nComp() const noexcept
Returns the number of components.
Definition AMReX_BaseFab.H:285
int maskGT(BaseFab< int > &mask, T const &val, int comp=0) const noexcept
Same as above except mark cells with value greater than val.
Definition AMReX_BaseFab.H:2395
T dotmask(const BaseFab< int > &mask, const Box &xbx, int xcomp, const BaseFab< T > &y, const Box &ybx, int ycomp, int numcomp) const noexcept
Definition AMReX_BaseFab.H:2673
BaseFab< T > & plus(const BaseFab< T > &src) noexcept
Definition AMReX_BaseFab.H:3436
BaseFab() noexcept=default
Construct an empty BaseFab, which must be resized (see BaseFab::resize) before use.
BaseFab< T > & divide(const BaseFab< T > &src) noexcept
Definition AMReX_BaseFab.H:3646
void prefetchToHost() const noexcept
Definition AMReX_BaseFab.H:1208
T min(int comp=0) const noexcept
Definition AMReX_BaseFab.H:2009
void setComplement(T const &x, const Box &b, int ns, int num) noexcept
This function is analogous to the fourth form of setVal above, except that instead of setting values ...
Definition AMReX_BaseFab.H:1834
BaseFab< T > & linInterp(const BaseFab< T > &f1, const Box &b1, int comp1, const BaseFab< T > &f2, const Box &b2, int comp2, Real t1, Real t2, Real t, const Box &b, int comp, int numcomp=1) noexcept
Linear interpolation / extrapolation. Result is (t2-t)/(t2-t1)*f1 + (t-t1)/(t2-t1)*f2 Data is taken f...
Definition AMReX_BaseFab.H:3240
IntVect minIndex(int comp=0) const noexcept
Definition AMReX_BaseFab.H:2205
T * dataPtr(const IntVect &p, int n=0) noexcept
Definition AMReX_BaseFab.H:1183
void abs() noexcept
Compute absolute value for all components of this FAB.
Definition AMReX_BaseFab.H:1842
void getVal(T *data, const IntVect &pos) const noexcept
Same as above, except that starts at component 0 and copies all comps.
Definition AMReX_BaseFab.H:1338
BaseFab< T > & plus(T const &val, Box const &bx, DestComp dcomp, NumComps ncomp) noexcept
Do nothing if bx is empty.
Definition AMReX_BaseFab.H:3420
Real norm(int p, int scomp=0, int numcomp=1) const
Compute the Lp-norm of this FAB using components (scomp : scomp+ncomp-1). p < 0 -> ERROR p = 0 -> inf...
Definition AMReX_BaseFab.H:1914
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:124
__host__ __device__ const int * hiVect() const &noexcept
Return a constant pointer the array of high end coordinates. Useful for calls to FORTRAN.
Definition AMReX_Box.H:195
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:364
__host__ __device__ bool sameSize(const BoxND &b) const noexcept
Return true is Boxes same size, ie translates of each other,. It is an error if they have different t...
Definition AMReX_Box.H:295
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:155
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:216
__host__ __device__ const int * loVect() const &noexcept
Return a constant pointer the array of low end coordinates. Useful for calls to FORTRAN.
Definition AMReX_Box.H:190
__host__ __device__ BoxND & setType(const IndexTypeND< dim > &t) noexcept
Set indexing type.
Definition AMReX_Box.H:513
__host__ __device__ bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:212
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:112
GPU-compatible tuple.
Definition AMReX_Tuple.H:98
static int deviceId() noexcept
Definition AMReX_GpuDevice.cpp:692
static int devicePropMajor() noexcept
Definition AMReX_GpuDevice.H:203
Definition AMReX_GpuElixir.H:13
__host__ static __device__ constexpr IntVectND< dim > TheMinVector() noexcept
Definition AMReX_IntVect.H:819
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
T * data() noexcept
Definition AMReX_PODVector.H:666
Definition AMReX_Reduce.H:438
Type value()
Definition AMReX_Reduce.H:473
Definition AMReX_Reduce.H:597
void eval(MF const &mf, IntVect const &nghost, D &reduce_data, F &&f)
Definition AMReX_Reduce.H:731
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ Dim3 ubound(Array4< T > const &a) noexcept
Return the inclusive upper bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1359
__host__ __device__ Dim3 length(Array4< T > const &a) noexcept
Return the spatial extents of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1373
__host__ __device__ Dim3 lbound(Array4< T > const &a) noexcept
Return the inclusive lower bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1345
std::array< T, N > Array
Definition AMReX_Array.H:26
__host__ __device__ AMREX_FORCE_INLINE T Exch(T *address, T val) noexcept
Definition AMReX_GpuAtomic.H:487
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void dtoh_memcpy_async(void *p_h, const void *p_d, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:435
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:88
bool notInLaunchRegion() noexcept
Definition AMReX_GpuControl.H:89
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:421
gpuStream_t gpuStream() noexcept
Definition AMReX_GpuDevice.H:291
__host__ __device__ AMREX_FORCE_INLINE void Add(T *const sum, T const value) noexcept
Definition AMReX_GpuAtomic.H:640
Definition AMReX_Amr.cpp:50
MakeType
Definition AMReX_MakeType.H:7
@ make_deep_copy
Definition AMReX_MakeType.H:7
@ make_alias
Definition AMReX_MakeType.H:7
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2852
__host__ __device__ Array4< T > makeArray4(T *p, Box const &bx, int ncomp) noexcept
Definition AMReX_BaseFab.H:92
RunOn
Definition AMReX_GpuControl.H:65
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2018
cudaStream_t gpuStream_t
Definition AMReX_GpuControl.H:79
bool InitSNaN() noexcept
Definition AMReX.cpp:185
void ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:202
Long TotalBytesAllocatedInFabs() noexcept
Definition AMReX_BaseFab.cpp:66
void placementNew(T *const, Long)
Definition AMReX_BaseFab.H:100
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
void BaseFab_Initialize()
Definition AMReX_BaseFab.cpp:30
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:25
void BaseFab_Finalize()
Definition AMReX_BaseFab.cpp:59
void ResetTotalBytesAllocatedInFabsHWM() noexcept
Definition AMReX_BaseFab.cpp:134
BoxList boxDiff(const Box &b1in, const Box &b2)
Returns BoxList defining the compliment of b2 in b1in.
Definition AMReX_BoxList.cpp:599
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
Long TotalBytesAllocatedInFabsHWM() noexcept
Definition AMReX_BaseFab.cpp:83
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:45
Long TotalCellsAllocatedInFabsHWM() noexcept
Definition AMReX_BaseFab.cpp:117
Long TotalCellsAllocatedInFabs() noexcept
Definition AMReX_BaseFab.cpp:100
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:235
__host__ __device__ bool almostEqual(T x, T y, int ulp=2)
Definition AMReX_Algorithm.H:123
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
void placementDelete(T *const, Long)
Definition AMReX_BaseFab.H:127
void LoopOnCpu(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition AMReX_Loop.H:365
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2028
void update_fab_stats(Long n, Long s, size_t szt) noexcept
Definition AMReX_BaseFab.cpp:146
A multidimensional array accessor.
Definition AMReX_Array4.H:285
__host__ __device__ constexpr std::size_t size() const noexcept
Total number of elements in the ArrayND's index region.
Definition AMReX_Array4.H:639
__host__ __device__ T * ptr(idx... i) const noexcept
Multi-index ptr() for accessing pointer to element.
Definition AMReX_Array4.H:557
Definition AMReX_DataAllocator.H:9
void * alloc(std::size_t sz) const noexcept
Definition AMReX_DataAllocator.H:16
Arena * arena() const noexcept
Definition AMReX_DataAllocator.H:24
Definition AMReX_DataAllocator.H:29
Definition AMReX_BaseFab.H:77
int i
Definition AMReX_BaseFab.H:80
__host__ __device__ DestComp(int ai) noexcept
Definition AMReX_BaseFab.H:79
Definition AMReX_Dim3.H:13
int x
Definition AMReX_Dim3.H:13
Definition AMReX_TypeTraits.H:51
Definition AMReX_BaseFab.H:83
__host__ __device__ NumComps(int an) noexcept
Definition AMReX_BaseFab.H:85
int n
Definition AMReX_BaseFab.H:86
Definition AMReX_BaseFab.H:71
__host__ __device__ SrcComp(int ai) noexcept
Definition AMReX_BaseFab.H:73
int i
Definition AMReX_BaseFab.H:74