Forward compatibility: flex
[foam-extend-3.2.git] / src / errorEstimation / errorEstimate / resErrorSup.C
blob73aeedeefbd8315d32b8a08ddaf0cdba9f47de8b
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 source operators.
27 \*---------------------------------------------------------------------------*/
29 #include "resErrorSup.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 namespace resError
41 template<class Type>
42 tmp<errorEstimate<Type> >
45     const volScalarField& sp,
46     GeometricField<Type, fvPatchField, volMesh>& vf
49     return tmp<errorEstimate<Type> >
50     (
51         new errorEstimate<Type>
52         (
53             vf,
54             sp.dimensions()*vf.dimensions()*dimVolume,
55             sp.internalField()*vf.internalField(),
56             scalarField(vf.internalField().size(), 0)
57         )
58     );
61 template<class Type>
62 tmp<errorEstimate<Type> >
65     const tmp<volScalarField>& tsp,
66     GeometricField<Type, fvPatchField, volMesh>& vf
69     tmp<errorEstimate<Type> > tee = resError::Sp(tsp(), vf);
70     tsp.clear();
71     return tee;
75 template<class Type>
76 tmp<errorEstimate<Type> >
79     const dimensionedScalar& sp,
80     GeometricField<Type, fvPatchField, volMesh>& vf
83     return tmp<errorEstimate<Type> >
84     (
85         new errorEstimate<Type>
86         (
87             vf,
88             sp.dimensions()*vf.dimensions()*dimVolume,
89             sp.value()*vf.internalField(),
90             scalarField(vf.internalField().size(), 0)
91         )
92     );
96 template<class Type>
97 tmp<errorEstimate<Type> >
98 SuSp
100     const volScalarField& sp,
101     GeometricField<Type, fvPatchField, volMesh>& vf
104     return Sp(sp, vf);
107 template<class Type>
108 tmp<errorEstimate<Type> >
109 SuSp
111     const tmp<volScalarField>& tsp,
112     GeometricField<Type, fvPatchField, volMesh>& vf
115     tmp<errorEstimate<Type> > tee = resError::SuSp(tsp(), vf);
116     tsp.clear();
117     return tee;
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 } // End namespace resError
125 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 } // End namespace Foam
129 // ************************************************************************* //
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //