ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
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_NumericalDiffusion.H>
#include <ERF_PlaneAverage.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
38 {
39  BL_PROFILE_REGION("erf_add_thin_body_sources()");
40 
41  const bool l_have_thin_xforce = (thin_xforce_lev != nullptr);
42  const bool l_have_thin_yforce = (thin_yforce_lev != nullptr);
43  const bool l_have_thin_zforce = (thin_zforce_lev != nullptr);
44 
45  // *****************************************************************************
46  // If a thin immersed body is present, add forcing terms
47  // *****************************************************************************
48  if (l_have_thin_xforce) {
49  MultiFab::Copy(*thin_xforce_lev, xmom_src, 0, 0, 1, 0);
50  thin_xforce_lev->mult(-1., 0, 1, 0);
51  ApplyInvertedMask(*thin_xforce_lev, *xflux_imask_lev, 0);
52  MultiFab::Add(xmom_src, *thin_xforce_lev, 0, 0, 1, 0);
53  }
54 
55  if (l_have_thin_yforce) {
56  MultiFab::Copy(*thin_yforce_lev, ymom_src, 0, 0, 1, 0);
57  thin_yforce_lev->mult(-1., 0, 1, 0);
58  ApplyInvertedMask(*thin_yforce_lev, *yflux_imask_lev, 0);
59  MultiFab::Add(ymom_src, *thin_yforce_lev, 0, 0, 1, 0);
60  }
61 
62  if (l_have_thin_zforce) {
63  MultiFab::Copy(*thin_zforce_lev, zmom_src, 0, 0, 1, 0);
64  thin_zforce_lev->mult(-1., 0, 1, 0);
65  ApplyInvertedMask(*thin_zforce_lev, *zflux_imask_lev, 0);
66  MultiFab::Add(zmom_src, *thin_zforce_lev, 0, 0, 1, 0);
67  }
68 
69 #if 0
70 #ifndef AMREX_USE_GPU
71  if (l_have_thin_xforce) {
72  // TODO: Implement particles to better track and output these data
73  if (nrk==2) {
74  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
75  {
76  const Box& tbx = mfi.nodaltilebox(0);
77  const Array4<const Real> & fx = thin_xforce_lev->const_array(mfi);
78  const Array4<const int> & mask = xflux_imask_lev->const_array(mfi);
79  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
80  if (mask(i,j,k)==0) {
81  amrex::AllPrint() << "thin body fx"<<IntVect(i,j,k)<<" = " << fx(i,j,k) << std::endl;
82  }
83  });
84  }
85  }
86  }
87 #endif
88 #endif
89 
90 #if 0
91 #ifndef AMREX_USE_GPU
92  if (l_have_thin_yforce) {
93  // TODO: Implement particles to better track and output these data
94  if (nrk==2) {
95  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
96  {
97  const Box& tby = mfi.nodaltilebox(1);
98  const Array4<const Real> & fy = thin_yforce_lev->const_array(mfi);
99  const Array4<const int> & mask = yflux_imask_lev->const_array(mfi);
100  ParallelFor(tby, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
101  if (mask(i,j,k)==0) {
102  amrex::AllPrint() << "thin body fy"<<IntVect(i,j,k)<<" = " << fy(i,j,k) << std::endl;
103  }
104  });
105  }
106  }
107  }
108 #endif
109 #endif
110 
111 #if 0
112 #ifndef AMREX_USE_GPU
113  if (l_have_thin_zforce) {
114  // TODO: Implement particles to better track and output these data
115  if (nrk==2) {
116  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
117  {
118  const Box& tbz = mfi.nodaltilebox(2);
119  const Array4<const Real> & fz = thin_zforce_lev->const_array(mfi);
120  const Array4<const int> & mask = zflux_imask_lev->const_array(mfi);
121  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
122  if (mask(i,j,k)==0) {
123  amrex::AllPrint() << "thin body fz"<<IntVect(i,j,k)<<" = " << fz(i,j,k) << std::endl;
124  }
125  });
126  }
127  }
128  }
129 #endif
130 #endif
131 }
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:383
@ cons
Definition: ERF_IndexDefines.H:139
Here is the call graph for this function: