1#ifdef REMORA_USE_PARTICLES
6#include <AMReX_TracerParticle_mod_K.H>
16void REMORAPC::EvolveParticles (
int a_lev,
18 Vector<MultiFab const*>& a_flow_vel,
19 const Vector<std::unique_ptr<MultiFab>>& a_z_phys_nd )
21 BL_PROFILE(
"REMORAPCPC::EvolveParticles()");
23 if (m_advect_w_flow) {
24 AdvectWithFlow( a_lev, a_dt_lev, a_flow_vel, a_z_phys_nd[a_lev] );
37void REMORAPC::AdvectWithFlow (
int a_lev,
39 Vector<MultiFab const*>& a_umac,
40 const std::unique_ptr<MultiFab>& a_z_height )
42 BL_PROFILE(
"REMORAPCPC::AdvectWithUmac()");
43 AMREX_ASSERT(OK(a_lev, a_lev, a_umac[0]->nGrow()-1));
44 AMREX_ASSERT(a_lev >= 0 && a_lev < GetParticles().size());
46 AMREX_D_TERM(AMREX_ASSERT(a_umac[0]->nGrow() >= 1);,
47 AMREX_ASSERT(a_umac[1]->nGrow() >= 1);,
48 AMREX_ASSERT(a_umac[2]->nGrow() >= 1););
50 AMREX_D_TERM(AMREX_ASSERT(!a_umac[0]->contains_nan());,
51 AMREX_ASSERT(!a_umac[1]->contains_nan());,
52 AMREX_ASSERT(!a_umac[2]->contains_nan()););
54 const auto strttime = amrex::second();
55 const Geometry& geom = m_gdb->Geom(a_lev);
56 const auto plo = geom.ProbLoArray();
57 const auto dxi = geom.InvCellSizeArray();
59 for (
int ipass = 0; ipass < 2; ipass++)
62#pragma omp parallel if (Gpu::notInLaunchRegion())
64 for (ParIterType pti(*
this, a_lev); pti.isValid(); ++pti)
66 int grid = pti.index();
67 auto& ptile = ParticlesAt(a_lev, pti);
68 auto& aos = ptile.GetArrayOfStructs();
69 auto& soa = ptile.GetStructOfArrays();
70 const int n = aos.numParticles();
71 auto *p_pbox = aos().data();
73 Array<ParticleReal*,AMREX_SPACEDIM> v_ptr;
74 v_ptr[0] = soa.GetRealData(REMORAParticlesRealIdxSoA::vx).data();
75 v_ptr[1] = soa.GetRealData(REMORAParticlesRealIdxSoA::vy).data();
76 v_ptr[2] = soa.GetRealData(REMORAParticlesRealIdxSoA::vz).data();
78 const FArrayBox* fab[AMREX_SPACEDIM] = { AMREX_D_DECL(&((*a_umac[0])[grid]),
79 &((*a_umac[1])[grid]),
80 &((*a_umac[2])[grid])) };
83 GpuArray<Array4<const Real>, AMREX_SPACEDIM>
84 const umacarr {fab[0]->array(), fab[1]->array(), fab[2]->array()};
86 bool use_terrain = (a_z_height !=
nullptr);
87 auto zheight = use_terrain ? (*a_z_height)[grid].array() : Array4<Real>{};
89 ParallelFor(n, [=] AMREX_GPU_DEVICE (
int i)
91 ParticleType& p = p_pbox[i];
92 if (p.id() <= 0) { return; }
94 ParticleReal v[AMREX_SPACEDIM];
96 mac_interpolate_mapped_z(p, plo, dxi, umacarr, zheight, v);
98 mac_interpolate(p, plo, dxi, umacarr, v);
102 for (
int dim=0; dim < AMREX_SPACEDIM; dim++)
104 v_ptr[dim][i] = p.pos(dim);
105 p.pos(dim) +=
static_cast<ParticleReal
>(ParticleReal(0.5)*a_dt*
v[dim]);
108 update_location_idata(p,plo,dxi,zheight);
110 for (
int dim=0; dim < AMREX_SPACEDIM; dim++)
112 p.pos(dim) = v_ptr[dim][i] +
static_cast<ParticleReal
>(a_dt*
v[dim]);
113 v_ptr[dim][i] =
v[dim];
116 update_location_idata(p,plo,dxi,zheight);
124 auto stoptime = amrex::second() - strttime;
127 Lazy::QueueReduction( [=] ()
mutable {
129 ParallelReduce::Max(stoptime, ParallelContext::IOProcessorNumberSub(),
130 ParallelContext::CommunicatorSub());
132 Print() <<
"REMORAPC::AdvectWithFlow() time: " << stoptime <<
'\n';