Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "linearUpwind.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/zeroGradientFvPatchField.H>
00029
00030
00031
00032 template<class Type>
00033 Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh> >
00034 Foam::linearUpwind<Type>::correction
00035 (
00036 const GeometricField<Type, fvPatchField, volMesh>& vf
00037 ) const
00038 {
00039 const fvMesh& mesh = this->mesh();
00040
00041 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsfCorr
00042 (
00043 new GeometricField<Type, fvsPatchField, surfaceMesh>
00044 (
00045 IOobject
00046 (
00047 "linearUpwindCorrection(" + vf.name() + ')',
00048 mesh.time().timeName(),
00049 mesh,
00050 IOobject::NO_READ,
00051 IOobject::NO_WRITE,
00052 false
00053 ),
00054 mesh,
00055 dimensioned<Type>(vf.name(), vf.dimensions(), pTraits<Type>::zero)
00056 )
00057 );
00058
00059 GeometricField<Type, fvsPatchField, surfaceMesh>& sfCorr = tsfCorr();
00060
00061 const surfaceScalarField& faceFlux = this->faceFlux_;
00062
00063 const labelList& owner = mesh.owner();
00064 const labelList& neighbour = mesh.neighbour();
00065
00066 const volVectorField& C = mesh.C();
00067 const surfaceVectorField& Cf = mesh.Cf();
00068
00069 GeometricField
00070 <typename outerProduct<vector, Type>::type, fvPatchField, volMesh>
00071 gradVf = gradScheme_().grad(vf);
00072
00073 forAll(faceFlux, facei)
00074 {
00075 if (faceFlux[facei] > 0)
00076 {
00077 label own = owner[facei];
00078 sfCorr[facei] = (Cf[facei] - C[own]) & gradVf[own];
00079 }
00080 else
00081 {
00082 label nei = neighbour[facei];
00083 sfCorr[facei] = (Cf[facei] - C[nei]) & gradVf[nei];
00084 }
00085 }
00086
00087
00088 typename GeometricField<Type, fvsPatchField, surfaceMesh>::
00089 GeometricBoundaryField& bSfCorr = sfCorr.boundaryField();
00090
00091 forAll(bSfCorr, patchi)
00092 {
00093 fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
00094
00095 if (pSfCorr.coupled())
00096 {
00097 const unallocLabelList& pOwner =
00098 mesh.boundary()[patchi].faceCells();
00099
00100 const vectorField& pCf = Cf.boundaryField()[patchi];
00101
00102 const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
00103
00104 Field<typename outerProduct<vector, Type>::type> pGradVfNei =
00105 gradVf.boundaryField()[patchi].patchNeighbourField();
00106
00107
00108 vectorField pd =
00109 mesh.Sf().boundaryField()[patchi]
00110 /(
00111 mesh.magSf().boundaryField()[patchi]
00112 *mesh.deltaCoeffs().boundaryField()[patchi]
00113 );
00114
00115 if (!mesh.orthogonal())
00116 {
00117 pd -= mesh.correctionVectors().boundaryField()[patchi]
00118 /mesh.deltaCoeffs().boundaryField()[patchi];
00119 }
00120
00121 forAll(pOwner, facei)
00122 {
00123 label own = pOwner[facei];
00124
00125 if (pFaceFlux[facei] > 0)
00126 {
00127 pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
00128 }
00129 else
00130 {
00131 pSfCorr[facei] =
00132 (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
00133 }
00134 }
00135 }
00136 }
00137
00138 return tsfCorr;
00139 }
00140
00141
00142 namespace Foam
00143 {
00144
00145 makelimitedSurfaceInterpolationTypeScheme(linearUpwind, scalar)
00146 makelimitedSurfaceInterpolationTypeScheme(linearUpwind, vector)
00147 }
00148
00149