ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ComputeStress_EB.cpp File Reference
#include <ERF_Diffusion.H>
Include dependency graph for ERF_ComputeStress_EB.cpp:

Functions

void ComputeStressConsVisc_EB (Box bxcc, Box tbxxy, Box tbxxz, Box tbxyz, Real mu_eff, const Array4< const Real > &cell_data, Array4< Real > &tau11, Array4< Real > &tau22, Array4< Real > &tau33, Array4< Real > &tau12, Array4< Real > &tau13, Array4< Real > &tau23, const Array4< const Real > &er_arr, Array4< const Real > &vfrac, Array4< Real > &tau13i, Array4< Real > &tau23i, Array4< Real > &tau33i)
 

Function Documentation

◆ ComputeStressConsVisc_EB()

void ComputeStressConsVisc_EB ( Box  bxcc,
Box  tbxxy,
Box  tbxxz,
Box  tbxyz,
Real  mu_eff,
const Array4< const Real > &  cell_data,
Array4< Real > &  tau11,
Array4< Real > &  tau22,
Array4< Real > &  tau33,
Array4< Real > &  tau12,
Array4< Real > &  tau13,
Array4< Real > &  tau23,
const Array4< const Real > &  er_arr,
Array4< const Real > &  vfrac,
Array4< Real > &  tau13i,
Array4< Real > &  tau23i,
Array4< Real > &  tau33i 
)

Function for computing the stress with constant viscosity for EB.

Parameters
[in]bxcccell center box for tau_ii
[in]tbxxynodal xy box for tau_12
[in]tbxxznodal xz box for tau_13
[in]tbxyznodal yz box for tau_23
[in]mu_effconstant molecular viscosity
[in]cell_datato access rho if ConstantAlpha
[in,out]tau1111 strain -> stress
[in,out]tau2222 strain -> stress
[in,out]tau3333 strain -> stress
[in,out]tau1212 strain -> stress
[in,out]tau1313 strain -> stress
[in,out]tau2323 strain -> stress
[in]er_arrexpansion rate
[in,out]tau13icontribution to stress from du/dz
[in,out]tau23icontribution to stress from dv/dz
[in,out]tau33icontribution to stress from dw/dz
35 {
36  Real OneThird = (1./3.);
37 
38  // NOTE: mu_eff includes factor of 2
39 
40  if (cell_data)
41  // constant alpha (stored in mu_eff)
42  {
43  // Cell centered strains
44  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
45  {
46  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
47  if (tau33i) tau33i(i,j,k) = -rhoAlpha * tau33(i,j,k);
48  tau11(i,j,k) = -rhoAlpha * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
49  tau22(i,j,k) = -rhoAlpha * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
50  tau33(i,j,k) = -rhoAlpha * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
51  });
52 
53  // Off-diagonal strains
54  ParallelFor(tbxxy,tbxxz,tbxyz,
55  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
56  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j-1,k) + vfrac(i-1,j-1,k);
57  Real rho_bar = 0.0;
58  if (vol_sum > 1.e-16) {
59  rho_bar = ( vfrac(i-1,j,k) * cell_data(i-1, j , k, Rho_comp)
60  + vfrac(i,j,k) * cell_data(i, j , k, Rho_comp)
61  + vfrac(i-1,j-1,k) * cell_data(i-1, j-1, k, Rho_comp)
62  + vfrac(i,j-1,k) * cell_data(i, j-1, k, Rho_comp) ) / vol_sum;
63  } else {
64  rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + vfrac(i,j,k) * cell_data(i, j , k, Rho_comp)
65  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
66  }
67  tau12(i,j,k) *= -rho_bar * mu_eff;
68  },
69  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
70  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j,k-1) + vfrac(i-1,j,k-1);
71  Real rho_bar = 0.0;
72  if (vol_sum > 1.e-16) {
73  rho_bar = ( vfrac(i-1,j,k) * cell_data(i-1, j, k , Rho_comp)
74  + vfrac(i,j,k) * cell_data(i, j, k , Rho_comp)
75  + vfrac(i-1,j,k-1) * cell_data(i-1, j, k-1, Rho_comp)
76  + vfrac(i,j,k-1) * cell_data(i, j, k-1, Rho_comp) )/ vol_sum;
77  } else {
78  rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
79  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
80  }
81  tau13(i,j,k) *= -rho_bar * mu_eff;
82 
83  if (tau13i) tau13i(i,j,k) *= -rho_bar * mu_eff;
84  },
85  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
86  Real vol_sum = vfrac(i,j,k) + vfrac(i,j-1,k) + vfrac(i,j,k-1) + vfrac(i,j-1,k-1);
87  Real rho_bar = 0.0;
88  if (vol_sum > 1.e-16) {
89  rho_bar = ( vfrac(i,j-1,k) * cell_data(i, j-1, k , Rho_comp)
90  + vfrac(i,j,k) * cell_data(i, j, k , Rho_comp)
91  + vfrac(i,j-1,k-1) * cell_data(i, j-1, k-1, Rho_comp)
92  + vfrac(i,j,k-1) * cell_data(i, j, k-1, Rho_comp) ) / vol_sum;
93  } else {
94  rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
95  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
96  }
97  tau23(i,j,k) *= -rho_bar * mu_eff;
98 
99  if (tau23i) tau23i(i,j,k) *= -rho_bar * mu_eff;
100  });
101  }
102  else
103  // constant mu_eff
104  {
105  // Cell centered strains
106  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
107  {
108  if (tau33i) tau33i(i,j,k) = -mu_eff * tau33(i,j,k);
109  tau11(i,j,k) = -mu_eff * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
110  tau22(i,j,k) = -mu_eff * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
111  tau33(i,j,k) = -mu_eff * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
112  });
113 
114  // Off-diagonal strains
115  ParallelFor(tbxxy,tbxxz,tbxyz,
116  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
117  tau12(i,j,k) *= -mu_eff;
118  },
119  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
120  tau13(i,j,k) *= -mu_eff;
121 
122  if (tau13i) tau13i(i,j,k) *= -mu_eff;
123  },
124  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
125  tau23(i,j,k) *= -mu_eff;
126 
127  if (tau23i) tau23i(i,j,k) *= -mu_eff;
128  });
129  }
130 }
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
#define Rho_comp
Definition: ERF_IndexDefines.H:36
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by erf_make_tau_terms().

Here is the caller graph for this function: