Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / errorEstimation / errorEstimate / errorEstimate.C
blob94b01e5dbe1e3dc3295d82a1b3b3d51e581b28e8
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 \*---------------------------------------------------------------------------*/
26 #include "errorEstimate.H"
27 #include "zeroGradientFvPatchField.H"
28 #include "fixedValueFvPatchField.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
35 template<class Type>
36 Foam::wordList Foam::errorEstimate<Type>::errorBCTypes() const
38     // Make the boundary condition type list
39     // Default types get over-ridden anyway
40     wordList ebct
41     (
42         psi_.boundaryField().size(),
43         zeroGradientFvPatchField<Type>::typeName
44     );
46     forAll (psi_.boundaryField(), patchI)
47     {
48         if (psi_.boundaryField()[patchI].fixesValue())
49         {
50             ebct[patchI] = fixedValueFvPatchField<Type>::typeName;
51         }
52     }
54     return ebct;
58 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
60 // Construct from components
61 template<class Type>
62 Foam::errorEstimate<Type>::errorEstimate
64     const GeometricField<Type, fvPatchField, volMesh>& psi,
65     const dimensionSet& ds,
66     const Field<Type>& res,
67     const scalarField& norm
70     psi_(psi),
71     dimensions_(ds),
72     residual_(res),
73     normFactor_(norm)
77 // Construct as copy
78 template<class Type>
79 Foam::errorEstimate<Type>::errorEstimate(const Foam::errorEstimate<Type>& ee)
81     refCount(),
82     psi_(ee.psi_),
83     dimensions_(ee.dimensions_),
84     residual_(ee.residual_),
85     normFactor_(ee.normFactor_)
89 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
91 template<class Type>
92 Foam::errorEstimate<Type>::~errorEstimate()
96 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
98 template<class Type>
99 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
100 Foam::errorEstimate<Type>::residual() const
102     tmp<GeometricField<Type, fvPatchField, volMesh> > tres
103     (
104         new GeometricField<Type, fvPatchField, volMesh>
105         (
106             IOobject
107             (
108                 "residual" + psi_.name(),
109                 psi_.mesh().time().timeName(),
110                 psi_.db(),
111                 IOobject::NO_READ,
112                 IOobject::NO_WRITE
113             ),
114             psi_.mesh(),
115             psi_.dimensions()/dimTime,
116             errorBCTypes()
117         )
118     );
120     GeometricField<Type, fvPatchField, volMesh>& res = tres();
122     res.internalField() = residual_;
123     res.boundaryField() == pTraits<Type>::zero;
125     res.correctBoundaryConditions();
127     return tres;
131 template<class Type>
132 Foam::tmp<Foam::volScalarField> Foam::errorEstimate<Type>::normFactor() const
134     tmp<volScalarField> tnormFactor
135     (
136         new volScalarField
137         (
138             IOobject
139             (
140                 "normFactor" + psi_.name(),
141                 psi_.mesh().time().timeName(),
142                 psi_.db(),
143                 IOobject::NO_READ,
144                 IOobject::NO_WRITE
145             ),
146             psi_.mesh(),
147             dimless/dimTime,
148             errorBCTypes()
149         )
150     );
152     volScalarField& normFactor = tnormFactor();
154     normFactor.internalField() = normFactor_;
155     normFactor.boundaryField() == pTraits<Type>::zero;
157     normFactor.correctBoundaryConditions();
159     return tnormFactor;
162 template<class Type>
163 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
164 Foam::errorEstimate<Type>::error() const
166     tmp<GeometricField<Type, fvPatchField, volMesh> > tresError
167     (
168         new GeometricField<Type, fvPatchField, volMesh>
169         (
170             IOobject
171             (
172                 "resError" + psi_.name(),
173                 psi_.mesh().time().timeName(),
174                 psi_.db(),
175                 IOobject::NO_READ,
176                 IOobject::NO_WRITE
177             ),
178             psi_.mesh(),
179             psi_.dimensions(),
180             errorBCTypes()
181         )
182     );
184     GeometricField<Type, fvPatchField, volMesh>& resError = tresError();
186     resError.internalField() = residual_/normFactor_;
187     resError.boundaryField() == pTraits<Type>::zero;
189     resError.correctBoundaryConditions();
191     return tresError;
194 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
196 template<class Type>
197 void Foam::errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>& rhs)
199     // Check for assignment to self
200     if (this == &rhs)
201     {
202         FatalErrorIn
203         (
204             "errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>&)"
205         )   << "Attempted assignment to self"
206             << abort(FatalError);
207     }
209     if (&psi_ != &(rhs.psi_))
210     {
211         FatalErrorIn
212         (
213             "errorEstimate<Type>::operator=(const errorEstimate<Type>&)"
214         )   << "different fields"
215             << abort(FatalError);
216     }
218     residual_ = rhs.residual_;
219     normFactor_ = rhs.normFactor_;
223 template<class Type>
224 void Foam::errorEstimate<Type>::operator=(const tmp<errorEstimate<Type> >& teev)
226     operator=(teev());
227     teev.clear();
231 template<class Type>
232 void Foam::errorEstimate<Type>::negate()
234     residual_.negate();
238 template<class Type>
239 void Foam::errorEstimate<Type>::operator+=(const errorEstimate<Type>& eev)
241     checkMethod(*this, eev, "+=");
243     dimensions_ += eev.dimensions_;
245     residual_ += eev.residual_;
246     normFactor_ += eev.normFactor_;
250 template<class Type>
251 void Foam::errorEstimate<Type>::operator+=
253     const tmp<errorEstimate<Type> >& teev
256     operator+=(teev());
257     teev.clear();
261 template<class Type>
262 void Foam::errorEstimate<Type>::operator-=(const errorEstimate<Type>& eev)
264     checkMethod(*this, eev, "+=");
266     dimensions_ -= eev.dimensions_;
267     residual_ -= eev.residual_;
268     normFactor_ += eev.normFactor_;
272 template<class Type>
273 void Foam::errorEstimate<Type>::operator-=(const tmp<errorEstimate<Type> >& teev)
275     operator-=(teev());
276     teev.clear();
280 template<class Type>
281 void Foam::errorEstimate<Type>::operator+=
283     const GeometricField<Type, fvPatchField, volMesh>& su
286     checkMethod(*this, su, "+=");
287     residual_ -= su.internalField();
291 template<class Type>
292 void Foam::errorEstimate<Type>::operator+=
294     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
297     operator+=(tsu());
298     tsu.clear();
302 template<class Type>
303 void Foam::errorEstimate<Type>::operator-=
305     const GeometricField<Type, fvPatchField, volMesh>& su
308     checkMethod(*this, su, "-=");
309     residual_ += su.internalField();
313 template<class Type>
314 void Foam::errorEstimate<Type>::operator-=
316     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
319     operator-=(tsu());
320     tsu.clear();
324 template<class Type>
325 void Foam::errorEstimate<Type>::operator+=
327     const dimensioned<Type>& su
330     residual_ -= su;
334 template<class Type>
335 void Foam::errorEstimate<Type>::operator-=
337     const dimensioned<Type>& su
340     residual_ += su;
344 template<class Type>
345 void Foam::errorEstimate<Type>::operator*=
347     const volScalarField& vsf
350     dimensions_ *= vsf.dimensions();
351     residual_ *= vsf.internalField();
352     normFactor_ *= vsf.internalField();
356 template<class Type>
357 void Foam::errorEstimate<Type>::operator*=
359     const tmp<volScalarField>& tvsf
362     operator*=(tvsf());
363     tvsf.clear();
367 template<class Type>
368 void Foam::errorEstimate<Type>::operator*=
370     const dimensioned<scalar>& ds
373     dimensions_ *= ds.dimensions();
374     residual_ *= ds.value();
375     normFactor_ *= ds.value();
379 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
381 template<class Type>
382 void Foam::checkMethod
384     const errorEstimate<Type>& ee1,
385     const errorEstimate<Type>& ee2,
386     const char* op
389     if (&ee1.psi() != &ee2.psi())
390     {
391         FatalErrorIn
392         (
393             "checkMethod(const errorEstimate<Type>&, "
394             "const errorEstimate<Type>&)"
395         )   << "incompatible fields for operation "
396             << endl << "    "
397             << "[" << ee1.psi().name() << "] "
398             << op
399             << " [" << ee2.psi().name() << "]"
400             << abort(FatalError);
401     }
403     if (dimensionSet::debug && ee1.dimensions() != ee2.dimensions())
404     {
405         FatalErrorIn
406         (
407             "checkMethod(const errorEstimate<Type>&, "
408             "const errorEstimate<Type>&)"
409         )   << "incompatible dimensions for operation "
410             << endl << "    "
411             << "[" << ee1.psi().name() << ee1.dimensions()/dimVolume << " ] "
412             << op
413             << " [" << ee2.psi().name() << ee2.dimensions()/dimVolume << " ]"
414             << abort(FatalError);
415     }
419 template<class Type>
420 void Foam::checkMethod
422     const errorEstimate<Type>& ee,
423     const GeometricField<Type, fvPatchField, volMesh>& vf,
424     const char* op
427     if (dimensionSet::debug && ee.dimensions()/dimVolume != vf.dimensions())
428     {
429         FatalErrorIn
430         (
431             "checkMethod(const errorEstimate<Type>&, "
432             "const GeometricField<Type, fvPatchField, volMesh>&)"
433         )   <<  "incompatible dimensions for operation "
434             << endl << "    "
435             << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
436             << op
437             << " [" << vf.name() << vf.dimensions() << " ]"
438             << abort(FatalError);
439     }
443 template<class Type>
444 void Foam::checkMethod
446     const errorEstimate<Type>& ee,
447     const dimensioned<Type>& dt,
448     const char* op
451     if (dimensionSet::debug && ee.dimensions()/dimVolume != dt.dimensions())
452     {
453         FatalErrorIn
454         (
455             "checkMethod(const errorEstimate<Type>&, const dimensioned<Type>&)"
456         )   << "incompatible dimensions for operation "
457             << endl << "    "
458             << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
459             << op
460             << " [" << dt.name() << dt.dimensions() << " ]"
461             << abort(FatalError);
462     }
466 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
468 namespace Foam
471 template<class Type>
472 tmp<errorEstimate<Type> > operator+
474     const errorEstimate<Type>& A,
475     const errorEstimate<Type>& B
478     checkMethod(A, B, "+");
479     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
480     tC() += B;
481     return tC;
485 template<class Type>
486 tmp<errorEstimate<Type> > operator+
488     const tmp<errorEstimate<Type> >& tA,
489     const errorEstimate<Type>& B
492     checkMethod(tA(), B, "+");
493     tmp<errorEstimate<Type> > tC(tA.ptr());
494     tC() += B;
495     return tC;
499 template<class Type>
500 tmp<errorEstimate<Type> > operator+
502     const errorEstimate<Type>& A,
503     const tmp<errorEstimate<Type> >& tB
506     checkMethod(A, tB(), "+");
507     tmp<errorEstimate<Type> > tC(tB.ptr());
508     tC() += A;
509     return tC;
513 template<class Type>
514 tmp<errorEstimate<Type> > operator+
516     const tmp<errorEstimate<Type> >& tA,
517     const tmp<errorEstimate<Type> >& tB
520     checkMethod(tA(), tB(), "+");
521     tmp<errorEstimate<Type> > tC(tA.ptr());
522     tC() += tB();
523     tB.clear();
524     return tC;
528 template<class Type>
529 tmp<errorEstimate<Type> > operator-
531     const errorEstimate<Type>& A
534     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
535     tC().negate();
536     return tC;
540 template<class Type>
541 tmp<errorEstimate<Type> > operator-
543     const tmp<errorEstimate<Type> >& tA
546     tmp<errorEstimate<Type> > tC(tA.ptr());
547     tC().negate();
548     return tC;
552 template<class Type>
553 tmp<errorEstimate<Type> > operator-
555     const errorEstimate<Type>& A,
556     const errorEstimate<Type>& B
559     checkMethod(A, B, "-");
560     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
561     tC() -= B;
562     return tC;
566 template<class Type>
567 tmp<errorEstimate<Type> > operator-
569     const tmp<errorEstimate<Type> >& tA,
570     const errorEstimate<Type>& B
573     checkMethod(tA(), B, "-");
574     tmp<errorEstimate<Type> > tC(tA.ptr());
575     tC() -= B;
576     return tC;
580 template<class Type>
581 tmp<errorEstimate<Type> > operator-
583     const errorEstimate<Type>& A,
584     const tmp<errorEstimate<Type> >& tB
587     checkMethod(A, tB(), "-");
588     tmp<errorEstimate<Type> > tC(tB.ptr());
589     tC() -= A;
590     tC().negate();
591     return tC;
595 template<class Type>
596 tmp<errorEstimate<Type> > operator-
598     const tmp<errorEstimate<Type> >& tA,
599     const tmp<errorEstimate<Type> >& tB
602     checkMethod(tA(), tB(), "-");
603     tmp<errorEstimate<Type> > tC(tA.ptr());
604     tC() -= tB();
605     tB.clear();
606     return tC;
610 template<class Type>
611 tmp<errorEstimate<Type> > operator==
613     const errorEstimate<Type>& A,
614     const errorEstimate<Type>& B
617     checkMethod(A, B, "==");
618     return (A - B);
622 template<class Type>
623 tmp<errorEstimate<Type> > operator==
625     const tmp<errorEstimate<Type> >& tA,
626     const errorEstimate<Type>& B
629     checkMethod(tA(), B, "==");
630     return (tA - B);
634 template<class Type>
635 tmp<errorEstimate<Type> > operator==
637     const errorEstimate<Type>& A,
638     const tmp<errorEstimate<Type> >& tB
641     checkMethod(A, tB(), "==");
642     return (A - tB);
646 template<class Type>
647 tmp<errorEstimate<Type> > operator==
649     const tmp<errorEstimate<Type> >& tA,
650     const tmp<errorEstimate<Type> >& tB
653     checkMethod(tA(), tB(), "==");
654     return (tA - tB);
658 template<class Type>
659 tmp<errorEstimate<Type> > operator+
661     const errorEstimate<Type>& A,
662     const GeometricField<Type, fvPatchField, volMesh>& su
665     checkMethod(A, su, "+");
666     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
667     tC().res() -= su.internalField();
668     return tC;
671 template<class Type>
672 tmp<errorEstimate<Type> > operator+
674     const tmp<errorEstimate<Type> >& tA,
675     const GeometricField<Type, fvPatchField, volMesh>& su
678     checkMethod(tA(), su, "+");
679     tmp<errorEstimate<Type> > tC(tA.ptr());
680     tC().res() -= su.internalField();
681     return tC;
684 template<class Type>
685 tmp<errorEstimate<Type> > operator+
687     const errorEstimate<Type>& A,
688     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
691     checkMethod(A, tsu(), "+");
692     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
693     tC().res() -= tsu().internalField();
694     tsu.clear();
695     return tC;
699 template<class Type>
700 tmp<errorEstimate<Type> > operator+
702     const tmp<errorEstimate<Type> >& tA,
703     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
706     checkMethod(tA(), tsu(), "+");
707     tmp<errorEstimate<Type> > tC(tA.ptr());
708     tC().res() -= tsu().internalField();
709     tsu.clear();
710     return tC;
713 template<class Type>
714 tmp<errorEstimate<Type> > operator+
716     const GeometricField<Type, fvPatchField, volMesh>& su,
717     const errorEstimate<Type>& A
720     checkMethod(A, su, "+");
721     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
722     tC().res() -= su.internalField();
723     return tC;
726 template<class Type>
727 tmp<errorEstimate<Type> > operator+
729     const GeometricField<Type, fvPatchField, volMesh>& su,
730     const tmp<errorEstimate<Type> >& tA
733     checkMethod(tA(), su, "+");
734     tmp<errorEstimate<Type> > tC(tA.ptr());
735     tC().res() -= su.internalField();
736     return tC;
739 template<class Type>
740 tmp<errorEstimate<Type> > operator+
742     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
743     const errorEstimate<Type>& A
746     checkMethod(A, tsu(), "+");
747     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
748     tC().res() -= tsu().internalField();
749     tsu.clear();
750     return tC;
753 template<class Type>
754 tmp<errorEstimate<Type> > operator+
756     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
757     const tmp<errorEstimate<Type> >& tA
760     checkMethod(tA(), tsu(), "+");
761     tmp<errorEstimate<Type> > tC(tA.ptr());
762     tC().res() -= tsu().internalField();
763     tsu.clear();
764     return tC;
768 template<class Type>
769 tmp<errorEstimate<Type> > operator-
771     const errorEstimate<Type>& A,
772     const GeometricField<Type, fvPatchField, volMesh>& su
775     checkMethod(A, su, "-");
776     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
777     tC().res() += su.internalField();
778     return tC;
781 template<class Type>
782 tmp<errorEstimate<Type> > operator-
784     const tmp<errorEstimate<Type> >& tA,
785     const GeometricField<Type, fvPatchField, volMesh>& su
788     checkMethod(tA(), su, "-");
789     tmp<errorEstimate<Type> > tC(tA.ptr());
790     tC().res() += su.internalField();
791     return tC;
794 template<class Type>
795 tmp<errorEstimate<Type> > operator-
797     const errorEstimate<Type>& A,
798     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
801     checkMethod(A, tsu(), "-");
802     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
803     tC().res() += tsu().internalField();
804     tsu.clear();
805     return tC;
808 template<class Type>
809 tmp<errorEstimate<Type> > operator-
811     const tmp<errorEstimate<Type> >& tA,
812     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
815     checkMethod(tA(), tsu(), "-");
816     tmp<errorEstimate<Type> > tC(tA.ptr());
817     tC().res() += tsu().internalField();
818     tsu.clear();
819     return tC;
823 template<class Type>
824 tmp<errorEstimate<Type> > operator-
826     const GeometricField<Type, fvPatchField, volMesh>& su,
827     const errorEstimate<Type>& A
830     checkMethod(A, su, "-");
831     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
832     tC().negate();
833     tC().res() -= su.internalField();
834     return tC;
838 template<class Type>
839 tmp<errorEstimate<Type> > operator-
841     const GeometricField<Type, fvPatchField, volMesh>& su,
842     const tmp<errorEstimate<Type> >& tA
845     checkMethod(tA(), su, "-");
846     tmp<errorEstimate<Type> > tC(tA.ptr());
847     tC().negate();
848     tC().res() -= su.internalField();
849     return tC;
852 template<class Type>
853 tmp<errorEstimate<Type> > operator-
855     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
856     const errorEstimate<Type>& A
859     checkMethod(A, tsu(), "-");
860     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
861     tC().negate();
862     tC().res() -= tsu().internalField();
863     tsu.clear();
864     return tC;
868 template<class Type>
869 tmp<errorEstimate<Type> > operator-
871     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
872     const tmp<errorEstimate<Type> >& tA
875     checkMethod(tA(), tsu(), "-");
876     tmp<errorEstimate<Type> > tC(tA.ptr());
877     tC().negate();
878     tC().res() -= tsu().internalField();
879     tsu.clear();
880     return tC;
884 template<class Type>
885 tmp<errorEstimate<Type> > operator+
887     const errorEstimate<Type>& A,
888     const dimensioned<Type>& su
891     checkMethod(A, su, "+");
892     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
893     tC().res() -= su.value();
894     return tC;
898 template<class Type>
899 tmp<errorEstimate<Type> > operator+
901     const tmp<errorEstimate<Type> >& tA,
902     const dimensioned<Type>& su
905     checkMethod(tA(), su, "+");
906     tmp<errorEstimate<Type> > tC(tA.ptr());
907     tC().res() -= su.value();
908     return tC;
912 template<class Type>
913 tmp<errorEstimate<Type> > operator+
915     const dimensioned<Type>& su,
916     const errorEstimate<Type>& A
919     checkMethod(A, su, "+");
920     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
921     tC().res() -= su.value();
922     return tC;
926 template<class Type>
927 tmp<errorEstimate<Type> > operator+
929     const dimensioned<Type>& su,
930     const tmp<errorEstimate<Type> >& tA
933     checkMethod(tA(), su, "+");
934     tmp<errorEstimate<Type> > tC(tA.ptr());
935     tC().res() -= su.value();
936     return tC;
940 template<class Type>
941 tmp<errorEstimate<Type> > operator-
943     const errorEstimate<Type>& A,
944     const dimensioned<Type>& su
947     checkMethod(A, su, "-");
948     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
949     tC().res() += su.value();
950     return tC;
954 template<class Type>
955 tmp<errorEstimate<Type> > operator-
957     const tmp<errorEstimate<Type> >& tA,
958     const dimensioned<Type>& su
961     checkMethod(tA(), su, "-");
962     tmp<errorEstimate<Type> > tC(tA.ptr());
963     tC().res() += su.value();
964     return tC;
968 template<class Type>
969 tmp<errorEstimate<Type> > operator-
971     const dimensioned<Type>& su,
972     const errorEstimate<Type>& A
975     checkMethod(A, su, "-");
976     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
977     tC().negate();
978     tC().res() -= su.value();
979     return tC;
983 template<class Type>
984 tmp<errorEstimate<Type> > operator-
986     const dimensioned<Type>& su,
987     const tmp<errorEstimate<Type> >& tA
990     checkMethod(tA(), su, "-");
991     tmp<errorEstimate<Type> > tC(tA.ptr());
992     tC().negate();
993     tC().res() -= su.value();
994     return tC;
998 template<class Type>
999 tmp<errorEstimate<Type> > operator==
1001     const errorEstimate<Type>& A,
1002     const GeometricField<Type, fvPatchField, volMesh>& su
1005     checkMethod(A, su, "==");
1006     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1007     tC().res() += su.internalField();
1008     return tC;
1011 template<class Type>
1012 tmp<errorEstimate<Type> > operator==
1014     const tmp<errorEstimate<Type> >& tA,
1015     const GeometricField<Type, fvPatchField, volMesh>& su
1018     checkMethod(tA(), su, "==");
1019     tmp<errorEstimate<Type> > tC(tA.ptr());
1020     tC().res() += su.internalField();
1021     return tC;
1024 template<class Type>
1025 tmp<errorEstimate<Type> > operator==
1027     const errorEstimate<Type>& A,
1028     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1031     checkMethod(A, tsu(), "==");
1032     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1033     tC().res() += tsu().internalField();
1034     tsu.clear();
1035     return tC;
1038 template<class Type>
1039 tmp<errorEstimate<Type> > operator==
1041     const tmp<errorEstimate<Type> >& tA,
1042     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1045     checkMethod(tA(), tsu(), "==");
1046     tmp<errorEstimate<Type> > tC(tA.ptr());
1047     tC().res() += tsu().internalField();
1048     tsu.clear();
1049     return tC;
1053 template<class Type>
1054 tmp<errorEstimate<Type> > operator==
1056     const errorEstimate<Type>& A,
1057     const dimensioned<Type>& su
1060     checkMethod(A, su, "==");
1061     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1062     tC().res() += su.value();
1063     return tC;
1067 template<class Type>
1068 tmp<errorEstimate<Type> > operator==
1070     const tmp<errorEstimate<Type> >& tA,
1071     const dimensioned<Type>& su
1074     checkMethod(tA(), su, "==");
1075     tmp<errorEstimate<Type> > tC(tA.ptr());
1076     tC().res() += su.value();
1077     return tC;
1081 template<class Type>
1082 tmp<errorEstimate<Type> > operator*
1084     const volScalarField& vsf,
1085     const errorEstimate<Type>& A
1088     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1089     tC() *= vsf;
1090     return tC;
1093 template<class Type>
1094 tmp<errorEstimate<Type> > operator*
1096     const tmp<volScalarField>& tvsf,
1097     const errorEstimate<Type>& A
1100     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1101     tC() *= tvsf;
1102     return tC;
1105 template<class Type>
1106 tmp<errorEstimate<Type> > operator*
1108     const volScalarField& vsf,
1109     const tmp<errorEstimate<Type> >& tA
1112     tmp<errorEstimate<Type> > tC(tA.ptr());
1113     tC() *= vsf;
1114     return tC;
1117 template<class Type>
1118 tmp<errorEstimate<Type> > operator*
1120     const tmp<volScalarField>& tvsf,
1121     const tmp<errorEstimate<Type> >& tA
1124     tmp<errorEstimate<Type> > tC(tA.ptr());
1125     tC() *= tvsf;
1126     return tC;
1130 template<class Type>
1131 tmp<errorEstimate<Type> > operator*
1133     const dimensioned<scalar>& ds,
1134     const errorEstimate<Type>& A
1137     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1138     tC() *= ds;
1139     return tC;
1143 template<class Type>
1144 tmp<errorEstimate<Type> > operator*
1146     const dimensioned<scalar>& ds,
1147     const tmp<errorEstimate<Type> >& tA
1150     tmp<errorEstimate<Type> > tC(tA.ptr());
1151     tC() *= ds;
1152     return tC;
1156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1158 } // End namespace Foam
1160 // ************************************************************************* //