constrainHbyA |
The boundary conditions of HbyA are now constrained by the new “constrainHbyA” function which applies the velocity boundary values for patches for which the velocity cannot be modified by assignment and pressure extrapolation is not specified via the new
“fixedFluxExtrapolatedPressureFvPatchScalarField”.
32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 |
Foam::tmp<Foam::volVectorField> Foam::constrainHbyA ( const tmp<volVectorField>& tHbyA, const volVectorField& U, const volScalarField& p ) { tmp<volVectorField> tHbyANew; if (tHbyA.isTmp()) { tHbyANew = tHbyA; tHbyANew.ref().rename("HbyA"); } else { tHbyANew = new volVectorField("HbyA", tHbyA); } volVectorField& HbyA = tHbyANew.ref(); volVectorField::Boundary& HbyAbf = HbyA.boundaryFieldRef(); forAll(U.boundaryField(), patchi) { if ( !U.boundaryField()[patchi].assignable() && !isA<fixedFluxExtrapolatedPressureFvPatchScalarField> ( p.boundaryField()[patchi] ) ) { HbyAbf[patchi] = U.boundaryField()[patchi]; } } return tHbyANew; } |
The member function assignable() is defined in the abstract base class of boundary conditions (fvPatchField):
320 321 322 323 324 325 |
//- Return true if the value of the patch field // is altered by assignment (the default) virtual bool assignable() const { return true; } |
It is overridden in the basic boundary condition classes (fixedValue, mixed, directionMixed and sliced):
163 164 165 166 167 |
//- Return false: this patch field is not altered by assignment virtual bool assignable() const { return false; } |
fvc::flux and dotInterpolate |
surfaceInterpolationScheme: Added dotInterpolate member-function
dotInterpolate interpolates the field and “dots” the resulting face-values with the vector field provided which removes the need to create a temporary field for the interpolate. This reduces the peak storage of OpenFOAM caused by the divergence of the gradient of vector fields, improves memory management and under some conditions decreases run-time.
This development is based on a patch contributed by Paul Edwards, Intel.
dotInterpolate function avoids storing an interpolated field in flux calculation in #OpenFOAM-dev, contrib @Intel, https://t.co/A3IjVo9SA7
— CFD Direct OpenFOAM (@CFDdirect) 2016年4月8日
31 32 33 34 35 36 37 38 39 40 41 |
Foam::tmp<Foam::surfaceScalarField> Foam::fvc::flux ( const volVectorField& vvf ) { return scheme<vector> ( vvf.mesh(), "flux(" + vvf.name() + ')' )().dotInterpolate(vvf.mesh().Sf(), vvf); } |
constrainPressure |
In the case of compressible flow, the following relation \eqref{eq:1} holds for the velocity on the face:
\begin{equation}
\rho \boldsymbol{U}_f \cdot \boldsymbol{S}_f = \left(\frac{\rho \boldsymbol{H}}{A_p}\right)_f \cdot \boldsymbol{S}_f\;-\;\left(\frac{\rho}{A_p}\right)_f (\nabla p)_f \cdot \boldsymbol{S}_f. \tag{1} \label{eq:1}
\end{equation}
Thus
\begin{equation}
(\nabla p)_f \cdot \boldsymbol{n}_f = \frac{1}{|\boldsymbol{S}_f|\left(\frac{\rho}{A_p}\right)_f} \left(\left(\frac{\rho \boldsymbol{H}}{A_p}\right)_f \cdot \boldsymbol{S}_f\;-\;\rho \boldsymbol{U}_f \cdot \boldsymbol{S}_f \right). \tag{2} \label{eq:2}
\end{equation}
The fixedFluxPressure boundary condition calculates the boundary normal gradient of the pressure field using the Eq. \eqref{eq:2} so that the boundary fluxes satisfy the Eq. \eqref{eq:1}.
34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 |
template<class RhoType, class RAUType, class MRFType> void Foam::constrainPressure ( volScalarField& p, const RhoType& rho, const volVectorField& U, const surfaceScalarField& phiHbyA, const RAUType& rhorAU, const MRFType& MRF ) { const fvMesh& mesh = p.mesh(); volScalarField::Boundary& pBf = p.boundaryFieldRef(); const volVectorField::Boundary& UBf = U.boundaryField(); const surfaceScalarField::Boundary& phiHbyABf = phiHbyA.boundaryField(); const typename RAUType::Boundary& rhorAUBf = rhorAU.boundaryField(); const surfaceVectorField::Boundary& SfBf = mesh.Sf().boundaryField(); const surfaceScalarField::Boundary& magSfBf = mesh.magSf().boundaryField(); forAll(pBf, patchi) { if (isA<fixedFluxPressureFvPatchScalarField>(pBf[patchi])) { refCast<fixedFluxPressureFvPatchScalarField> ( pBf[patchi] ).updateCoeffs ( ( phiHbyABf[patchi] - rho.boundaryField()[patchi] *MRF.relative(SfBf[patchi] & UBf[patchi], patchi) ) /(magSfBf[patchi]*rhorAUBf[patchi]) ); } } } |