Forward compatibility: flex
[foam-extend-3.2.git] / src / errorEstimation / errorEstimate / resErrorDiv.C
blobb1d14b5f66942a66a1de116cec50ea644d7d2755
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Description
25     Residual error estimate for the fv convection operators.
27 \*---------------------------------------------------------------------------*/
29 #include "resErrorDiv.H"
30 #include "fvc.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 namespace resError
42 template<class Type>
43 tmp<errorEstimate<Type> >
44 div
46     const surfaceScalarField& flux,
47     const GeometricField<Type, fvPatchField, volMesh>& vf
50     const fvMesh& mesh = vf.mesh();
52     const scalarField& vols = mesh.V();
53     const surfaceVectorField& faceCentres = mesh.Cf();
54     const volVectorField& cellCentres = mesh.C();
55     const fvPatchList& patches = mesh.boundary();
56     const unallocLabelList& owner = mesh.owner();
57     const unallocLabelList& neighbour = mesh.neighbour();
59     Field<Type> res(vols.size(), pTraits<Type>::zero);
60     scalarField aNorm(vols.size(), 0.0);
62     // Get sign of flux
63     const surfaceScalarField signF = pos(flux);
65     // Calculate gradient of the solution
66     GeometricField
67     <
68         typename outerProduct<vector, Type>::type, fvPatchField, volMesh
69     > gradVf = fvc::grad(vf);
71     // Internal faces
72     forAll (owner, faceI)
73     {
74         // Calculate the centre of the face
75         const vector& curFaceCentre = faceCentres[faceI];
77         // Owner
78         vector ownD = curFaceCentre - cellCentres[owner[faceI]];
80         // Subtract convection
81         res[owner[faceI]] -=
82         (
83             vf[owner[faceI]]
84           + (ownD & gradVf[owner[faceI]])
85         )*flux[faceI];
87         aNorm[owner[faceI]] += signF[faceI]*flux[faceI];
89         // Neighbour
90         vector neiD = curFaceCentre - cellCentres[neighbour[faceI]];
92         // Subtract convection
93         res[neighbour[faceI]] +=
94         (
95             vf[neighbour[faceI]]
96           + (neiD & gradVf[neighbour[faceI]])
97         )*flux[faceI];
99         aNorm[neighbour[faceI]] -= (1.0 - signF[faceI])*flux[faceI];
100     }
102     forAll (patches, patchI)
103     {
104         const vectorField& patchFaceCentres =
105             faceCentres.boundaryField()[patchI];
107         const scalarField& patchFlux = flux.boundaryField()[patchI];
108         const scalarField& patchSignFlux = signF.boundaryField()[patchI];
110         const labelList& fCells = patches[patchI].faceCells();
112         forAll (fCells, faceI)
113         {
114             vector d =
115                 patchFaceCentres[faceI] - cellCentres[fCells[faceI]];
117             // Subtract convection
118             res[fCells[faceI]] -=
119             (
120                 vf[fCells[faceI]]
121               + (d & gradVf[fCells[faceI]])
122             )*patchFlux[faceI];
124             aNorm[fCells[faceI]] += patchSignFlux[faceI]*patchFlux[faceI];
125         }
126     }
128     res /= vols;
129     aNorm /= vols;
131     return tmp<errorEstimate<Type> >
132     (
133         new errorEstimate<Type>
134         (
135             vf,
136             flux.dimensions()*vf.dimensions(),
137             res,
138             aNorm
139         )
140     );
144 template<class Type>
145 tmp<errorEstimate<Type> >
148     const tmp<surfaceScalarField>& tflux,
149     const GeometricField<Type, fvPatchField, volMesh>& vf
152     tmp<errorEstimate<Type> > Div(resError::div(tflux(), vf));
153     tflux.clear();
154     return Div;
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 } // End namespace resError
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 } // End namespace Foam
166 // ************************************************************************* //
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //