1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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/>.
28 Residual error estimation
33 \*---------------------------------------------------------------------------*/
35 #ifndef errorEstimate_H
36 #define errorEstimate_H
38 #include "volFields.H"
39 #include "surfaceFields.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 /*---------------------------------------------------------------------------*\
47 Class errorEstimate Declaration
48 \*---------------------------------------------------------------------------*/
57 // Reference to GeometricField<Type, fvPatchField, volMesh>
58 const GeometricField<Type, fvPatchField, volMesh>& psi_;
61 dimensionSet dimensions_;
63 //- Cell residual pointer
64 Field<Type> residual_;
66 //- Normalisation factor
67 scalarField normFactor_;
70 // Private Member Functions
72 //- Return boundary condition types for the error field
73 wordList errorBCTypes() const;
77 // Static data members
79 ClassName("errorEstimate");
84 //- Construct from components
87 const GeometricField<Type, fvPatchField, volMesh>& psi,
88 const dimensionSet& ds,
89 const Field<Type>& res,
90 const scalarField& norm
94 errorEstimate(const errorEstimate<Type>&);
107 const GeometricField<Type, fvPatchField, volMesh>& psi() const
112 //- Return residual dimensions
113 const dimensionSet& dimensions() const
118 // Raw residual (for calculus)
125 const Field<Type>& res() const
133 //- Cell residual (volume intensive)
134 tmp<GeometricField<Type, fvPatchField, volMesh> > residual() const;
136 //- Normalisation factor
137 tmp<volScalarField> normFactor() const;
140 tmp<GeometricField<Type, fvPatchField, volMesh> > error() const;
145 void operator=(const errorEstimate<Type>&);
146 void operator=(const tmp<errorEstimate<Type> >&);
150 void operator+=(const errorEstimate<Type>&);
151 void operator+=(const tmp<errorEstimate<Type> >&);
153 void operator-=(const errorEstimate<Type>&);
154 void operator-=(const tmp<errorEstimate<Type> >&);
156 void operator+=(const GeometricField<Type,fvPatchField,volMesh>&);
157 void operator+=(const tmp<GeometricField<Type,fvPatchField,volMesh> >&);
159 void operator-=(const GeometricField<Type,fvPatchField,volMesh>&);
160 void operator-=(const tmp<GeometricField<Type,fvPatchField,volMesh> >&);
162 void operator+=(const dimensioned<Type>&);
163 void operator-=(const dimensioned<Type>&);
165 void operator*=(const volScalarField&);
166 void operator*=(const tmp<volScalarField>&);
168 void operator*=(const dimensioned<scalar>&);
177 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
182 const errorEstimate<Type>&,
183 const errorEstimate<Type>&,
190 const errorEstimate<Type>&,
191 const GeometricField<Type, fvPatchField, volMesh>&,
198 const errorEstimate<Type>&,
199 const dimensioned<Type>&,
204 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
207 tmp<errorEstimate<Type> > operator-
209 const errorEstimate<Type>&
213 tmp<errorEstimate<Type> > operator-
215 const tmp<errorEstimate<Type> >&
219 tmp<errorEstimate<Type> > operator+
221 const errorEstimate<Type>&,
222 const errorEstimate<Type>&
226 tmp<errorEstimate<Type> > operator+
228 const tmp<errorEstimate<Type> >&,
229 const errorEstimate<Type>&
233 tmp<errorEstimate<Type> > operator+
235 const errorEstimate<Type>&,
236 const tmp<errorEstimate<Type> >&
240 tmp<errorEstimate<Type> > operator+
242 const tmp<errorEstimate<Type> >&,
243 const tmp<errorEstimate<Type> >&
247 tmp<errorEstimate<Type> > operator-
249 const errorEstimate<Type>&,
250 const errorEstimate<Type>&
254 tmp<errorEstimate<Type> > operator-
256 const tmp<errorEstimate<Type> >&,
257 const errorEstimate<Type>&
261 tmp<errorEstimate<Type> > operator-
263 const errorEstimate<Type>&,
264 const tmp<errorEstimate<Type> >&
268 tmp<errorEstimate<Type> > operator-
270 const tmp<errorEstimate<Type> >&,
271 const tmp<errorEstimate<Type> >&
275 tmp<errorEstimate<Type> > operator==
277 const errorEstimate<Type>&,
278 const errorEstimate<Type>&
282 tmp<errorEstimate<Type> > operator==
284 const tmp<errorEstimate<Type> >&,
285 const errorEstimate<Type>&
289 tmp<errorEstimate<Type> > operator==
291 const errorEstimate<Type>&,
292 const tmp<errorEstimate<Type> >&
296 tmp<errorEstimate<Type> > operator==
298 const tmp<errorEstimate<Type> >&,
299 const tmp<errorEstimate<Type> >&
303 tmp<errorEstimate<Type> > operator+
305 const errorEstimate<Type>&,
306 const GeometricField<Type, fvPatchField, volMesh>&
310 tmp<errorEstimate<Type> > operator+
312 const tmp<errorEstimate<Type> >&,
313 const GeometricField<Type, fvPatchField, volMesh>&
317 tmp<errorEstimate<Type> > operator+
319 const errorEstimate<Type>&,
320 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
324 tmp<errorEstimate<Type> > operator+
326 const tmp<errorEstimate<Type> >&,
327 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
331 tmp<errorEstimate<Type> > operator+
333 const GeometricField<Type, fvPatchField, volMesh>&,
334 const errorEstimate<Type>&
338 tmp<errorEstimate<Type> > operator+
340 const GeometricField<Type, fvPatchField, volMesh>&,
341 const tmp<errorEstimate<Type> >&
345 tmp<errorEstimate<Type> > operator+
347 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
348 const errorEstimate<Type>&
352 tmp<errorEstimate<Type> > operator+
354 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
355 const tmp<errorEstimate<Type> >&
359 tmp<errorEstimate<Type> > operator-
361 const errorEstimate<Type>&,
362 const GeometricField<Type, fvPatchField, volMesh>&
366 tmp<errorEstimate<Type> > operator-
368 const tmp<errorEstimate<Type> >&,
369 const GeometricField<Type, fvPatchField, volMesh>&
373 tmp<errorEstimate<Type> > operator-
375 const errorEstimate<Type>&,
376 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
380 tmp<errorEstimate<Type> > operator-
382 const tmp<errorEstimate<Type> >&,
383 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
387 tmp<errorEstimate<Type> > operator-
389 const GeometricField<Type, fvPatchField, volMesh>&,
390 const errorEstimate<Type>&
394 tmp<errorEstimate<Type> > operator-
396 const GeometricField<Type, fvPatchField, volMesh>&,
397 const tmp<errorEstimate<Type> >&
401 tmp<errorEstimate<Type> > operator-
403 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
404 const errorEstimate<Type>&
408 tmp<errorEstimate<Type> > operator-
410 const tmp<GeometricField<Type, fvPatchField, volMesh> >&,
411 const tmp<errorEstimate<Type> >&
415 tmp<errorEstimate<Type> > operator+
417 const tmp<errorEstimate<Type> >&,
418 const dimensioned<Type>&
422 tmp<errorEstimate<Type> > operator+
424 const dimensioned<Type>&,
425 const tmp<errorEstimate<Type> >&
429 tmp<errorEstimate<Type> > operator-
431 const tmp<errorEstimate<Type> >&,
432 const dimensioned<Type>&
436 tmp<errorEstimate<Type> > operator-
438 const dimensioned<Type>&,
439 const tmp<errorEstimate<Type> >&
443 tmp<errorEstimate<Type> > operator==
445 const errorEstimate<Type>&,
446 const GeometricField<Type, fvPatchField, volMesh>&
450 tmp<errorEstimate<Type> > operator==
452 const tmp<errorEstimate<Type> >&,
453 const GeometricField<Type, fvPatchField, volMesh>&
457 tmp<errorEstimate<Type> > operator==
459 const errorEstimate<Type>&,
460 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
464 tmp<errorEstimate<Type> > operator==
466 const tmp<errorEstimate<Type> >&,
467 const tmp<GeometricField<Type, fvPatchField, volMesh> >&
471 tmp<errorEstimate<Type> > operator==
473 const errorEstimate<Type>&,
474 const dimensioned<Type>&
478 tmp<errorEstimate<Type> > operator==
480 const tmp<errorEstimate<Type> >&,
481 const dimensioned<Type>&
486 tmp<errorEstimate<Type> > operator*
488 const volScalarField&,
489 const errorEstimate<Type>&
493 tmp<errorEstimate<Type> > operator*
495 const volScalarField&,
496 const tmp<errorEstimate<Type> >&
500 tmp<errorEstimate<Type> > operator*
502 const tmp<volScalarField>&,
503 const errorEstimate<Type>&
507 tmp<errorEstimate<Type> > operator*
509 const tmp<volScalarField>&,
510 const tmp<errorEstimate<Type> >&
515 tmp<errorEstimate<Type> > operator*
517 const dimensioned<scalar>&,
518 const errorEstimate<Type>&
522 tmp<errorEstimate<Type> > operator*
524 const dimensioned<scalar>&,
525 const tmp<errorEstimate<Type> >&
529 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
531 } // End namespace Foam
533 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
536 # include "errorEstimate.C"
539 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
543 // ************************************************************************* //