Merge branch 'master' of ssh://opencfd@repo.or.cz/srv/git/OpenFOAM-1.5.x
[OpenFOAM-1.5.x.git] / src / errorEstimation / errorEstimate / resErrorDiv.C
blob67e8e6060a3355b3d15c2e48e79781450753f3b8
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "resErrorDiv.H"
28 #include "fvc.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace resError
40 template<class Type>
41 tmp<errorEstimate<Type> >
42 div
44     const surfaceScalarField& flux,
45     const GeometricField<Type, fvPatchField, volMesh>& vf
48     const fvMesh& mesh = vf.mesh();
50     const scalarField& vols = mesh.V();
51     const surfaceVectorField& faceCentres = mesh.Cf();
52     const volVectorField& cellCentres = mesh.C();
53     const fvPatchList& patches = mesh.boundary();
54     const unallocLabelList& owner = mesh.owner();
55     const unallocLabelList& neighbour = mesh.neighbour();
57     Field<Type> res(vols.size(), pTraits<Type>::zero);
58     scalarField aNorm(vols.size(), 0.0);
60     // Get sign of flux
61     const surfaceScalarField signF = pos(flux);
63     // Calculate gradient of the solution
64     GeometricField
65     <
66         typename outerProduct<vector, Type>::type, fvPatchField, volMesh
67     > gradVf = fvc::grad(vf);
69     // Internal faces
70     forAll (owner, faceI)
71     {
72         // Calculate the centre of the face
73         const vector& curFaceCentre = faceCentres[faceI];
75         // Owner
76         vector ownD = curFaceCentre - cellCentres[owner[faceI]];
78         // Subtract convection
79         res[owner[faceI]] -=
80         (
81             vf[owner[faceI]]
82           + (ownD & gradVf[owner[faceI]])
83         )*flux[faceI];
85         aNorm[owner[faceI]] += signF[faceI]*flux[faceI];
87         // Neighbour
88         vector neiD = curFaceCentre - cellCentres[neighbour[faceI]];
90         // Subtract convection
91         res[neighbour[faceI]] +=
92         (
93             vf[neighbour[faceI]]
94           + (neiD & gradVf[neighbour[faceI]])
95         )*flux[faceI];
97         aNorm[neighbour[faceI]] -= (1.0 - signF[faceI])*flux[faceI];
98     }
100     forAll (patches, patchI)
101     {
102         const vectorField& patchFaceCentres =
103             faceCentres.boundaryField()[patchI];
105         const scalarField& patchFlux = flux.boundaryField()[patchI];
106         const scalarField& patchSignFlux = signF.boundaryField()[patchI];
108         const labelList& fCells = patches[patchI].faceCells();
110         forAll (fCells, faceI)
111         {
112             vector d =
113                 patchFaceCentres[faceI] - cellCentres[fCells[faceI]];
115             // Subtract convection
116             res[fCells[faceI]] -=
117             (
118                 vf[fCells[faceI]]
119               + (d & gradVf[fCells[faceI]])
120             )*patchFlux[faceI];
122             aNorm[fCells[faceI]] += patchSignFlux[faceI]*patchFlux[faceI];
123         }
124     }
126     res /= vols;
127     aNorm /= vols;
129     return tmp<errorEstimate<Type> >
130     (
131         new errorEstimate<Type>
132         (
133             vf,
134             flux.dimensions()*vf.dimensions(),
135             res,
136             aNorm
137         )
138     );
142 template<class Type>
143 tmp<errorEstimate<Type> >
146     const tmp<surfaceScalarField>& tflux,
147     const GeometricField<Type, fvPatchField, volMesh>& vf
150     tmp<errorEstimate<Type> > Div(resError::div(tflux(), vf));
151     tflux.clear();
152     return Div;
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 } // End namespace resError
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 } // End namespace Foam
164 // ************************************************************************* //
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //