REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_BoundaryConditions_zvel.cpp
Go to the documentation of this file.
1#include "AMReX_PhysBCFunct.H"
4
5using namespace amrex;
6
7/**
8 * @param[inout] dest_arr data on which to apply BCs
9 * @param[in ] bx box to update on
10 * @param[in ] domain domain box
11 * @param[in ] dxInv pm or pn
12 * @param[in ] mskr land-sea mask on rho-points
13 * @param[in ] calc_arr data to use in the RHS of calculations
14 * @param[in ] time current time
15 * @param[in ] bccomp index into both domain_bcs_type_bcr and bc_extdir_vals for icomp=0
16 */
17void REMORAPhysBCFunct::impose_zvel_bcs (const Array4<Real>& dest_arr, const Box& bx, const Box& domain,
18 const GpuArray<Real,AMREX_SPACEDIM> /*dxInv*/,const Array4<const Real>& /*mskr*/,
19 Real /*time*/, int bccomp)
20{
21 const auto& dom_lo = amrex::lbound(domain);
22 const auto& dom_hi = amrex::ubound(domain);
23
24 // Based on BCRec for the domain, we need to make BCRec for this Box
25 // bccomp is used as starting index for m_domain_bcs_type
26 // 0 is used as starting index for bcrs
27 int ncomp = 1;
28 Vector<BCRec> bcrs(ncomp);
29 amrex::setBC(bx, domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
30
31 // xlo: ori = 0
32 // ylo: ori = 1
33 // zlo: ori = 2
34 // xhi: ori = 3
35 // yhi: ori = 4
36 // zhi: ori = 5
37
38 //! Based on BCRec for the domain, we need to make BCRec for this Box
39 // bccomp is used as starting index for m_domain_bcs_type
40 // 0 is used as starting index for bcrs
41 amrex::setBC(bx, domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
42
43 amrex::Gpu::DeviceVector<BCRec> bcrs_d(ncomp);
44#ifdef AMREX_USE_GPU
45 Gpu::htod_memcpy_async(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
46#else
47 std::memcpy(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
48#endif
49 const amrex::BCRec* bc_ptr = bcrs_d.data();
50 const auto* bc_extdir_vals_ptr = m_bc_extdir_vals_d.data();
51
52 GeometryData const& geomdata = m_geom.data();
53 bool is_periodic_in_x = geomdata.isPeriodic(0);
54 bool is_periodic_in_y = geomdata.isPeriodic(1);
55
56 // First do all ext_dir bcs
57 if (!is_periodic_in_x or bccomp == BCVars::foextrap_bc(m_ncons))
58 {
59 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
60 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
61 ParallelFor(
62 bx_xlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
63 int iflip = dom_lo.x - 1 - i;
64 if (bc_ptr[n].lo(0) == REMORABCType::ext_dir) {
65 dest_arr(i,j,k) = bc_extdir_vals_ptr[bccomp+n][0];
66 } else if (bc_ptr[n].lo(0) == REMORABCType::foextrap || bc_ptr[n].lo(0) == REMORABCType::clamped) {
67 dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
68 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_even) {
69 dest_arr(i,j,k) = dest_arr(iflip,j,k);
70 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_odd) {
71 dest_arr(i,j,k) = -dest_arr(iflip,j,k);
72 }
73 },
74 bx_xhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
75 int iflip = 2*dom_hi.x + 1 - i;
76 if (bc_ptr[n].hi(0) == REMORABCType::ext_dir) {
77 dest_arr(i,j,k) = bc_extdir_vals_ptr[bccomp+n][3];
78 } else if (bc_ptr[n].hi(0) == REMORABCType::foextrap || bc_ptr[n].hi(0) == REMORABCType::clamped) {
79 dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
80 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_even) {
81 dest_arr(i,j,k) = dest_arr(iflip,j,k);
82 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_odd) {
83 dest_arr(i,j,k) = -dest_arr(iflip,j,k);
84 }
85 }
86 );
87 }
88
89 // First do all ext_dir bcs
90 if (!is_periodic_in_y or bccomp == BCVars::foextrap_bc(m_ncons))
91 {
92 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
93 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
94 ParallelFor(bx_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
95 int jflip = dom_lo.y - 1 - j;
96 if (bc_ptr[n].lo(1) == REMORABCType::ext_dir) {
97 dest_arr(i,j,k) = bc_extdir_vals_ptr[bccomp+n][1];
98 } else if (bc_ptr[n].lo(1) == REMORABCType::foextrap || bc_ptr[n].lo(1) == REMORABCType::clamped) {
99 dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
100 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_even) {
101 dest_arr(i,j,k) = dest_arr(i,jflip,k);
102 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_odd) {
103 dest_arr(i,j,k) = -dest_arr(i,jflip,k);
104 }
105 },
106 bx_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
107 int jflip = 2*dom_hi.y + 1 - j;
108 if (bc_ptr[n].hi(1) == REMORABCType::ext_dir) {
109 dest_arr(i,j,k) = bc_extdir_vals_ptr[bccomp+n][4];
110 } else if (bc_ptr[n].hi(1) == REMORABCType::foextrap || bc_ptr[n].hi(1) == REMORABCType::clamped) {
111 dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
112 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_even) {
113 dest_arr(i,j,k) = dest_arr(i,jflip,k);
114 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_odd) {
115 dest_arr(i,j,k) = -dest_arr(i,jflip,k);
116 }
117 });
118 }
119
120 // Never need to apply a boundary condition to the top
121
122 Gpu::streamSynchronize();
123}
amrex::Vector< amrex::BCRec > m_domain_bcs_type
amrex::Geometry m_geom
amrex::Gpu::DeviceVector< amrex::GpuArray< amrex::Real, AMREX_SPACEDIM *2 > > m_bc_extdir_vals_d
void impose_zvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, const amrex::Array4< const amrex::Real > &mskr, amrex::Real time, int bccomp)
apply z-velocity type boundary conditions
int foextrap_bc(int ncons) noexcept