1#ifndef REMORA_PhysBCFunct_H_
2#define REMORA_PhysBCFunct_H_
3#include <AMReX_Config.H>
5#include <AMReX_BCRec.H>
6#include <AMReX_Geometry.H>
7#include <AMReX_MultiFab.H>
8#include <AMReX_ArrayLim.H>
9#include <AMReX_FilCC_C.H>
10#include <AMReX_FilND_C.H>
11#include <AMReX_FilFC_C.H>
12#include "AMReX_TypeTraits.H"
13#include "AMReX_Orientation.H"
26 const amrex::Geometry& geom,
const amrex::Vector<amrex::BCRec>& domain_bcs_type,
27 const amrex::Gpu::DeviceVector<amrex::BCRec>& domain_bcs_type_d,
28 const amrex::Vector<amrex::GpuArray<amrex::Real,AMREX_SPACEDIM*2>>& bc_extdir_vals,
37 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
45 void operator() (amrex::MultiFab& mf,
const amrex::MultiFab& mask,
int icomp,
int ncomp, amrex::IntVect
const& nghost,
46 amrex::Real time,
int bccomp,
int n_not_fill=0,
47 const amrex::MultiFab& mf_calc = amrex::MultiFab(),
48 const amrex::MultiFab& mf_msku = amrex::MultiFab(),
49 const amrex::MultiFab& mf_mskv = amrex::MultiFab());
52 void impose_xvel_bcs (
const amrex::Array4<amrex::Real>& dest_arr,
const amrex::Box& bx,
const amrex::Box& domain,
53 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> dxInv,
const amrex::Array4<const amrex::Real>& msku,
54 const amrex::Array4<const amrex::Real>& calc_arr, amrex::Real time,
int bccomp);
57 void impose_yvel_bcs (
const amrex::Array4<amrex::Real>& dest_arr,
const amrex::Box& bx,
const amrex::Box& domain,
58 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> dxInv,
const amrex::Array4<const amrex::Real>& mskv,
59 const amrex::Array4<const amrex::Real>& calc_arr, amrex::Real time,
int bccomp);
64 const amrex::Box& bx,
const amrex::Box& domain,
65 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> dxInv,
const amrex::Array4<const amrex::Real>& mskr,
66 amrex::Real time,
int bccomp);
69 void impose_cons_bcs (
const amrex::Array4<amrex::Real>& mf,
const amrex::Box& bx,
const amrex::Box& valid_bx,
70 const amrex::Box& domain,
71 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> dxInv,
const amrex::Array4<const amrex::Real>& mskr,
72 const amrex::Array4<const amrex::Real>& msku,
const amrex::Array4<const amrex::Real>& mskv,
73 const amrex::Array4<const amrex::Real>& calc_arr,
int icomp,
int ncomp, amrex::Real time,
int bccomp,
int n_not_fill);
amrex::Vector< amrex::FArrayBox > PlaneVector
void impose_xvel_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 > &msku, const amrex::Array4< const amrex::Real > &calc_arr, amrex::Real time, int bccomp)
apply x-velocity type boundary conditions
amrex::Vector< amrex::GpuArray< amrex::Real, AMREX_SPACEDIM *2 > > m_bc_extdir_vals
amrex::Vector< amrex::BCRec > m_domain_bcs_type
void impose_cons_bcs(const amrex::Array4< amrex::Real > &mf, const amrex::Box &bx, const amrex::Box &valid_bx, const amrex::Box &domain, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, const amrex::Array4< const amrex::Real > &mskr, const amrex::Array4< const amrex::Real > &msku, const amrex::Array4< const amrex::Real > &mskv, const amrex::Array4< const amrex::Real > &calc_arr, int icomp, int ncomp, amrex::Real time, int bccomp, int n_not_fill)
apply scalar type boundary conditions
REMORAPhysBCFunct(const int lev, const amrex::Geometry &geom, const amrex::Vector< amrex::BCRec > &domain_bcs_type, const amrex::Gpu::DeviceVector< amrex::BCRec > &domain_bcs_type_d, const amrex::Vector< amrex::GpuArray< amrex::Real, AMREX_SPACEDIM *2 > > &bc_extdir_vals, int ncons)
Constructor for physical boundary condition functor.
amrex::Gpu::DeviceVector< amrex::BCRec > m_domain_bcs_type_d
void impose_yvel_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 > &mskv, const amrex::Array4< const amrex::Real > &calc_arr, amrex::Real time, int bccomp)
apply y-velocity type boundary conditions
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
void operator()(amrex::MultiFab &mf, const amrex::MultiFab &mask, int icomp, int ncomp, amrex::IntVect const &nghost, amrex::Real time, int bccomp, int n_not_fill=0, const amrex::MultiFab &mf_calc=amrex::MultiFab(), const amrex::MultiFab &mf_msku=amrex::MultiFab(), const amrex::MultiFab &mf_mskv=amrex::MultiFab())
apply boundary condition to mf