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 "leastSquaresGrad.H"
00027 #include "leastSquaresVectors.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 leastSquaresGrad<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 leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
00090
00091 const surfaceVectorField& ownLs = lsv.pVectors();
00092 const surfaceVectorField& neiLs = lsv.nVectors();
00093
00094 const unallocLabelList& own = mesh.owner();
00095 const unallocLabelList& nei = mesh.neighbour();
00096
00097 forAll(own, facei)
00098 {
00099 register label ownFaceI = own[facei];
00100 register label neiFaceI = nei[facei];
00101
00102 Type deltaVsf = vsf[neiFaceI] - vsf[ownFaceI];
00103
00104 lsGrad[ownFaceI] += ownLs[facei]*deltaVsf;
00105 lsGrad[neiFaceI] -= neiLs[facei]*deltaVsf;
00106 }
00107
00108
00109 forAll(vsf.boundaryField(), patchi)
00110 {
00111 const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
00112
00113 const unallocLabelList& faceCells =
00114 lsGrad.boundaryField()[patchi].patch().faceCells();
00115
00116 if (vsf.boundaryField()[patchi].coupled())
00117 {
00118 Field<Type> neiVsf =
00119 vsf.boundaryField()[patchi].patchNeighbourField();
00120
00121 forAll(neiVsf, patchFaceI)
00122 {
00123 lsGrad[faceCells[patchFaceI]] +=
00124 patchOwnLs[patchFaceI]
00125 *(neiVsf[patchFaceI] - vsf[faceCells[patchFaceI]]);
00126 }
00127 }
00128 else
00129 {
00130 const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
00131
00132 forAll(patchVsf, patchFaceI)
00133 {
00134 lsGrad[faceCells[patchFaceI]] +=
00135 patchOwnLs[patchFaceI]
00136 *(patchVsf[patchFaceI] - vsf[faceCells[patchFaceI]]);
00137 }
00138 }
00139 }
00140
00141
00142 lsGrad.correctBoundaryConditions();
00143 gaussGrad<Type>::correctBoundaryConditions(vsf, lsGrad);
00144
00145 return tlsGrad;
00146 }
00147
00148
00149
00150
00151 }
00152
00153
00154
00155 }
00156
00157