ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_AddThinBodySources.cpp File Reference
#include <AMReX_MultiFab.H>
#include <AMReX_ArrayLim.H>
#include <AMReX_BCRec.H>
#include <AMReX_TableData.H>
#include <AMReX_GpuContainers.H>
#include "ERF_TI_slow_headers.H"
#include "ERF_SrcHeaders.H"
#include "ERF_Utils.H"
Include dependency graph for ERF_AddThinBodySources.cpp:

Functions

void add_thin_body_sources (MultiFab &xmom_src, MultiFab &ymom_src, MultiFab &zmom_src, std::unique_ptr< iMultiFab > &xflux_imask_lev, std::unique_ptr< iMultiFab > &yflux_imask_lev, std::unique_ptr< iMultiFab > &zflux_imask_lev, std::unique_ptr< MultiFab > &thin_xforce_lev, std::unique_ptr< MultiFab > &thin_yforce_lev, std::unique_ptr< MultiFab > &thin_zforce_lev)
 

Function Documentation

◆ add_thin_body_sources()

void add_thin_body_sources ( MultiFab &  xmom_src,
MultiFab &  ymom_src,
MultiFab &  zmom_src,
std::unique_ptr< iMultiFab > &  xflux_imask_lev,
std::unique_ptr< iMultiFab > &  yflux_imask_lev,
std::unique_ptr< iMultiFab > &  zflux_imask_lev,
std::unique_ptr< MultiFab > &  thin_xforce_lev,
std::unique_ptr< MultiFab > &  thin_yforce_lev,
std::unique_ptr< MultiFab > &  thin_zforce_lev 
)

Function for computing the slow RHS for the evolution equations for the density, potential temperature and momentum.

Parameters
[in]xmom_srcsource terms for x-momentum
[in]ymom_srcsource terms for y-momentum
[in]zmom_srcsource terms for z-momentum
[in]xflux_imask_levthin-body mask on x-faces
[in]yflux_imask_levthin-body mask on y-faces
[in]zflux_imask_levthin-body mask on z-faces
[in]thin_xforce_levx-component of forces on thin immersed bodies
[in]thin_yforce_levy-component of forces on thin immersed bodies
[in]thin_zforce_levz-component of forces on thin immersed bodies
36 {
37  BL_PROFILE_REGION("erf_add_thin_body_sources()");
38 
39  const bool l_have_thin_xforce = (thin_xforce_lev != nullptr);
40  const bool l_have_thin_yforce = (thin_yforce_lev != nullptr);
41  const bool l_have_thin_zforce = (thin_zforce_lev != nullptr);
42 
43  // *****************************************************************************
44  // If a thin immersed body is present, add forcing terms
45  // *****************************************************************************
46  if (l_have_thin_xforce) {
47  MultiFab::Copy(*thin_xforce_lev, xmom_src, 0, 0, 1, 0);
48  thin_xforce_lev->mult(-1., 0, 1, 0);
49  ApplyInvertedMask(*thin_xforce_lev, *xflux_imask_lev, 0);
50  MultiFab::Add(xmom_src, *thin_xforce_lev, 0, 0, 1, 0);
51  }
52 
53  if (l_have_thin_yforce) {
54  MultiFab::Copy(*thin_yforce_lev, ymom_src, 0, 0, 1, 0);
55  thin_yforce_lev->mult(-1., 0, 1, 0);
56  ApplyInvertedMask(*thin_yforce_lev, *yflux_imask_lev, 0);
57  MultiFab::Add(ymom_src, *thin_yforce_lev, 0, 0, 1, 0);
58  }
59 
60  if (l_have_thin_zforce) {
61  MultiFab::Copy(*thin_zforce_lev, zmom_src, 0, 0, 1, 0);
62  thin_zforce_lev->mult(-1., 0, 1, 0);
63  ApplyInvertedMask(*thin_zforce_lev, *zflux_imask_lev, 0);
64  MultiFab::Add(zmom_src, *thin_zforce_lev, 0, 0, 1, 0);
65  }
66 
67 #if 0
68 #ifndef AMREX_USE_GPU
69  if (l_have_thin_xforce) {
70  // TODO: Implement particles to better track and output these data
71  if (nrk==2) {
72  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
73  {
74  const Box& tbx = mfi.nodaltilebox(0);
75  const Array4<const Real> & fx = thin_xforce_lev->const_array(mfi);
76  const Array4<const int> & mask = xflux_imask_lev->const_array(mfi);
77  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
78  if (mask(i,j,k)==0) {
79  amrex::AllPrint() << "thin body fx"<<IntVect(i,j,k)<<" = " << fx(i,j,k) << std::endl;
80  }
81  });
82  }
83  }
84  }
85 #endif
86 #endif
87 
88 #if 0
89 #ifndef AMREX_USE_GPU
90  if (l_have_thin_yforce) {
91  // TODO: Implement particles to better track and output these data
92  if (nrk==2) {
93  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
94  {
95  const Box& tby = mfi.nodaltilebox(1);
96  const Array4<const Real> & fy = thin_yforce_lev->const_array(mfi);
97  const Array4<const int> & mask = yflux_imask_lev->const_array(mfi);
98  ParallelFor(tby, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
99  if (mask(i,j,k)==0) {
100  amrex::AllPrint() << "thin body fy"<<IntVect(i,j,k)<<" = " << fy(i,j,k) << std::endl;
101  }
102  });
103  }
104  }
105  }
106 #endif
107 #endif
108 
109 #if 0
110 #ifndef AMREX_USE_GPU
111  if (l_have_thin_zforce) {
112  // TODO: Implement particles to better track and output these data
113  if (nrk==2) {
114  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
115  {
116  const Box& tbz = mfi.nodaltilebox(2);
117  const Array4<const Real> & fz = thin_zforce_lev->const_array(mfi);
118  const Array4<const int> & mask = zflux_imask_lev->const_array(mfi);
119  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
120  if (mask(i,j,k)==0) {
121  amrex::AllPrint() << "thin body fz"<<IntVect(i,j,k)<<" = " << fz(i,j,k) << std::endl;
122  }
123  });
124  }
125  }
126  }
127 #endif
128 #endif
129 }
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask(amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
Definition: ERF_Utils.H:409
@ cons
Definition: ERF_IndexDefines.H:158
Here is the call graph for this function: