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 "extendedLeastSquaresGrad.H"
00027 #include "extendedLeastSquaresVectors.H"
00028 #include <finiteVolume/gaussGrad.H>
00029 #include <finiteVolume/fvMesh.H>
00030 #include <finiteVolume/volMesh.H>
00031 #include <finiteVolume/surfaceMesh.H>
00032 #include <OpenFOAM/GeometricField.H>
00033 #include <finiteVolume/zeroGradientFvPatchField.H>
00034
00035
00036
00037 namespace Foam
00038 {
00039
00040
00041
00042 namespace fv
00043 {
00044
00045
00046
00047 template<class Type>
00048 tmp
00049 <
00050 GeometricField
00051 <
00052 typename outerProduct<vector, Type>::type, fvPatchField, volMesh
00053 >
00054 >
00055 extendedLeastSquaresGrad<Type>::grad
00056 (
00057 const GeometricField<Type, fvPatchField, volMesh>& vsf
00058 ) const
00059 {
00060 typedef typename outerProduct<vector, Type>::type GradType;
00061
00062 const fvMesh& mesh = vsf.mesh();
00063
00064 tmp<GeometricField<GradType, fvPatchField, volMesh> > tlsGrad
00065 (
00066 new GeometricField<GradType, fvPatchField, volMesh>
00067 (
00068 IOobject
00069 (
00070 "grad("+vsf.name()+')',
00071 vsf.instance(),
00072 mesh,
00073 IOobject::NO_READ,
00074 IOobject::NO_WRITE
00075 ),
00076 mesh,
00077 dimensioned<GradType>
00078 (
00079 "zero",
00080 vsf.dimensions()/dimLength,
00081 pTraits<GradType>::zero
00082 ),
00083 zeroGradientFvPatchField<GradType>::typeName
00084 )
00085 );
00086 GeometricField<GradType, fvPatchField, volMesh>& lsGrad = tlsGrad();
00087
00088
00089 const extendedLeastSquaresVectors& lsv = extendedLeastSquaresVectors::New
00090 (
00091 mesh,
00092 minDet_
00093 );
00094
00095 const surfaceVectorField& ownLs = lsv.pVectors();
00096 const surfaceVectorField& neiLs = lsv.nVectors();
00097
00098 const unallocLabelList& owner = mesh.owner();
00099 const unallocLabelList& neighbour = mesh.neighbour();
00100
00101 forAll(owner, facei)
00102 {
00103 label own = owner[facei];
00104 label nei = neighbour[facei];
00105
00106 Type deltaVsf = vsf[nei] - vsf[own];
00107
00108 lsGrad[own] += ownLs[facei]*deltaVsf;
00109 lsGrad[nei] -= neiLs[facei]*deltaVsf;
00110 }
00111
00112
00113 forAll(vsf.boundaryField(), patchi)
00114 {
00115 const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
00116
00117 const unallocLabelList& faceCells =
00118 lsGrad.boundaryField()[patchi].patch().faceCells();
00119
00120 if (vsf.boundaryField()[patchi].coupled())
00121 {
00122 Field<Type> neiVsf =
00123 vsf.boundaryField()[patchi].patchNeighbourField();
00124
00125 forAll(neiVsf, patchFaceI)
00126 {
00127 lsGrad[faceCells[patchFaceI]] +=
00128 patchOwnLs[patchFaceI]
00129 *(neiVsf[patchFaceI] - vsf[faceCells[patchFaceI]]);
00130 }
00131 }
00132 else
00133 {
00134 const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
00135
00136 forAll(patchVsf, patchFaceI)
00137 {
00138 lsGrad[faceCells[patchFaceI]] +=
00139 patchOwnLs[patchFaceI]
00140 *(patchVsf[patchFaceI] - vsf[faceCells[patchFaceI]]);
00141 }
00142 }
00143 }
00144
00145
00146 const List<labelPair>& additionalCells = lsv.additionalCells();
00147 const vectorField& additionalVectors = lsv.additionalVectors();
00148
00149 forAll(additionalCells, i)
00150 {
00151 lsGrad[additionalCells[i][0]] +=
00152 additionalVectors[i]
00153 *(vsf[additionalCells[i][1]] - vsf[additionalCells[i][0]]);
00154 }
00155
00156
00157 lsGrad.correctBoundaryConditions();
00158 gaussGrad<Type>::correctBoundaryConditions(vsf, lsGrad);
00159
00160 return tlsGrad;
00161 }
00162
00163
00164
00165
00166 }
00167
00168
00169
00170 }
00171
00172