REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_PhysBCFunct.cpp
Go to the documentation of this file.
1#include "AMReX_PhysBCFunct.H"
4
5using namespace amrex;
6
7/**
8 * @param[inout] mf multifab to be filled
9 * @param[in ] msk land/sea mask for the variable
10 * @param[in ] icomp index into the MultiFab -- if cell-centered this can be any value from 0 to NCONS-1, if face-centered can be any value from 0 to 2 (inclusive)
11 * @param[in ] ncomp number of components -- if cell-centered this can be any value from 1 to NCONS as long as icomp+ncomp <= NCONS-1. If face-centered this must be 1
12 * @param[in ] nghost how many ghost cells to be filled
13 * @param[in ] time time at which the data should be filled
14 * @param[in ] bccomp index into both domain_bcs_type_bcr and bc_extdir_vals for icomp = 0 so this follows the BCVars enum
15 * @param[in ] n_not_fill halo size to not fill
16 * @param[in ] mf_calc data to use in calculation of RHS
17 * @param[in ] mf_msku land-sea mask at u-points
18 * @param[in ] mf_mskv land-sea mask at v-points
19*/
20void REMORAPhysBCFunct::operator() (MultiFab& mf, const MultiFab& msk, int icomp, int ncomp, IntVect const& nghost,
21 Real time, int bccomp,int n_not_fill, const MultiFab& mf_calc,
22 const MultiFab& mf_msku, const MultiFab& mf_mskv)
23{
24 if (m_geom.isAllPeriodic()) return;
25
26 BL_PROFILE("REMORAPhysBCFunct::()");
27
28 const auto& domain = m_geom.Domain();
29 const auto dxInv = m_geom.InvCellSizeArray();
30
31 // Create a grown domain box containing valid + periodic cells
32 Box gdomain = amrex::convert(domain, mf.boxArray().ixType());
33 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
34 if (m_geom.isPeriodic(i) and bccomp != BCVars::foextrap_bc) {
35 gdomain.grow(i, nghost[i]);
36 }
37 }
38 const bool null_mf_calc = (mf_calc.ok()) ? false : true;
39 const bool null_mf_msku = (mf_msku.ok()) ? false : true;
40 const bool null_mf_mskv = (mf_mskv.ok()) ? false : true;
41
42#ifdef AMREX_USE_OMP
43#pragma omp parallel if (Gpu::notInLaunchRegion())
44#endif
45 {
46 if (mf.boxArray()[0].ixType() == IndexType(IntVect(0,0,0))) {
47
48 // Cell-centered arrays only
49 for (MFIter mfi(mf); mfi.isValid(); ++mfi)
50 {
51 const Array4<Real>& dest_arr = mf.array(mfi);
52 const Array4<const Real>& msk_arr = msk.array(mfi);
53 //const Array4<const Real>& calc_arr = mf_calc.array(mfi);
54 const Array4<const Real>& calc_arr = (!null_mf_calc) ? mf_calc.const_array(mfi) : Array4<const Real>();
55 const Array4<const Real>& msku_arr = (!null_mf_msku) ? mf_msku.const_array(mfi) : Array4<const Real>();
56 const Array4<const Real>& mskv_arr = (!null_mf_mskv) ? mf_mskv.const_array(mfi) : Array4<const Real>();
57 Box bx = mfi.validbox(); bx.grow(nghost);
58 Box valid_bx = mfi.validbox();
59
60 if (!gdomain.contains(bx)) {
61 impose_cons_bcs(dest_arr,bx,valid_bx,domain,dxInv,msk_arr,
62 msku_arr,mskv_arr,calc_arr,icomp,ncomp,time,bccomp,n_not_fill);
63 }
64 } // mfi
65
66 } else {
67
68 // Face-based arrays only
69 for (MFIter mfi(mf); mfi.isValid(); ++mfi)
70 {
71 Box bx = mfi.validbox(); bx.grow(nghost);
72 const Array4<const Real>& msk_arr = msk.array(mfi);
73 Array4<const Real> calc_arr = Array4<const Real>();
74 if (!null_mf_calc) {
75 calc_arr = mf_calc.const_array(mfi);
76 }
77 if (!gdomain.contains(bx)) {
78 if(bx.ixType() == IndexType(IntVect(1,0,0))) {
79 const Array4<Real>& dest_arr = mf.array(mfi,icomp);
80 impose_xvel_bcs(dest_arr,bx,domain,dxInv,msk_arr,calc_arr,time,bccomp);
81 } else if (bx.ixType() == IndexType(IntVect(0,1,0))) {
82 const Array4<Real>& dest_arr = mf.array(mfi,icomp);
83 impose_yvel_bcs(dest_arr,bx,domain,dxInv,msk_arr,calc_arr,time,bccomp);
84 } else if (bx.ixType() == IndexType(IntVect(0,0,1))) {
85 const Array4<Real>& dest_arr = mf.array(mfi,icomp);
86 impose_zvel_bcs(dest_arr,bx,domain,dxInv,msk_arr,time,bccomp);
87 }
88 }
89 } // mfi
90 } // box type
91 } // OpenMP
92} // operator()
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
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
amrex::Geometry m_geom
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
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