blob: 91bb802e81f4df81e01eb32029fdcd8ac7b70e47 (
plain) (
blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
|
#include <WarpX.H>
#include <AMReX_BoxIterator.H>
#include <array>
#include <algorithm>
using namespace amrex;
void
WarpX::ErrorEst (int lev, TagBoxArray& tags, Real time, int /*ngrow*/)
{
const Real* problo = Geom(lev).ProbLo();
const Real* dx = Geom(lev).CellSize();
#ifdef _OPENMP
#pragma omp parallel
#endif
for (MFIter mfi(tags); mfi.isValid(); ++mfi)
{
auto& fab = tags[mfi];
const Box& bx = fab.box();
for (BoxIterator bi(bx); bi.ok(); ++bi)
{
const IntVect& cell = bi();
RealVect pos {AMREX_D_DECL((cell[0]+0.5_rt)*dx[0]+problo[0],
(cell[1]+0.5_rt)*dx[1]+problo[1],
(cell[2]+0.5_rt)*dx[2]+problo[2])};
if (pos > fine_tag_lo && pos < fine_tag_hi) {
fab(cell) = TagBox::SET;
}
}
}
}
|