What we want to achieve
When we simulate fluid flow, we have to cut a finite computational domain out of an entire flow region. For accurate simulation, we need to let fluid and sound wave flow smoothly out of the domain through the boundary.
The reflection at the boundary has a larger effect on the solution especially when we perform a compressible flow simulation as can be seen in the following movie (Upper: With reflection, Lower: Without reflection of sound wave).
In OpenFOAM, we can use two approximate non-reflecting boundary conditions:
They determine the boundary value by solving the following equation
\begin{align}
\frac{D \phi}{D t} = \frac{\partial \phi}{\partial t} + \boldsymbol{U} \cdot \nabla \phi = 0, \tag{1} \label{eq:advection}
\end{align}
where \(D/Dt\) is the material derivative and \(\boldsymbol{U}(\boldsymbol{x}, t)\) is the advection velocity.
We assume that the advection velocity \(\boldsymbol{U}\) is parallel to the boundary (face) normal direction and rewrite the eqn. \eqref{eq:advection} as
\begin{align}
\frac{D \phi}{Dt} \approx \frac{\partial \phi}{\partial t} + U_{n} \cdot \frac{\partial \phi}{\partial \boldsymbol{n}}= 0, \tag{2} \label{eq:advection2}
\end{align}
where \(\boldsymbol{n}\) is the outward-pointing unit normal vector.
These boundary conditions are different in how the advection speed (scalar quantity) \(U_{n}\) is calculated and it is calculated in advectionSpeed() member function.
advective B.C.
The advection speed is the component of the velocity normal to the boundary
\begin{align}
U_n = u_n. \tag{3} \label{eq:advectiveUn}
\end{align}
154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 |
template<class Type> Foam::tmp<Foam::scalarField> Foam::advectiveFvPatchField<Type>::advectionSpeed() const { const surfaceScalarField& phi = this->db().objectRegistry::template lookupObject<surfaceScalarField> (phiName_); fvsPatchField<scalar> phip = this->patch().template lookupPatchField<surfaceScalarField, scalar> ( phiName_ ); if (phi.dimensions() == dimDensity*dimVelocity*dimArea) { const fvPatchScalarField& rhop = this->patch().template lookupPatchField<volScalarField, scalar> ( rhoName_ ); return phip/(rhop*this->patch().magSf()); } else { return phip/this->patch().magSf(); } } |
waveTransmissive B.C.
The advection speed is the sum of the component of the velocity normal to the boundary and the speed of sound \(c\)
\begin{align}
U_n = u_n + c = u_n + \sqrt{\gamma/\psi}, \tag{4} \label{eq:waveTransmissiveUn}
\end{align}
where \(\gamma\) is the ratio of specific heats \(C_p/C_v\) and \(\psi\) is compressibility.
105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 |
template<class Type> Foam::tmp<Foam::scalarField> Foam::waveTransmissiveFvPatchField<Type>::advectionSpeed() const { // Lookup the velocity and compressibility of the patch const fvPatchField<scalar>& psip = this->patch().template lookupPatchField<volScalarField, scalar>(psiName_); const surfaceScalarField& phi = this->db().template lookupObject<surfaceScalarField>(this->phiName_); fvsPatchField<scalar> phip = this->patch().template lookupPatchField<surfaceScalarField, scalar>(this->phiName_); if (phi.dimensions() == dimDensity*dimVelocity*dimArea) { const fvPatchScalarField& rhop = this->patch().template lookupPatchField<volScalarField, scalar>(this->rhoName_); phip /= rhop; } // Calculate the speed of the field wave w // by summing the component of the velocity normal to the boundary // and the speed of sound (sqrt(gamma_/psi)). return phip/this->patch().magSf() + sqrt(gamma_/psip); } |
What do lInf and fieldInf mean?
Coming soon.