Fix tutorials: coupled/conjugateHeatFoam/conjugateCavity: fix Allrun file
[OpenFOAM-1.6-ext.git] / src / errorEstimation / errorEstimate / errorEstimate.C
blob4c338d0f6a1b60cda3a7a94db5463e3b41f27c96
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "errorEstimate.H"
28 #include "zeroGradientFvPatchField.H"
29 #include "fixedValueFvPatchField.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
36 template<class Type>
37 Foam::wordList Foam::errorEstimate<Type>::errorBCTypes() const
39     // Make the boundary condition type list
40     // Default types get over-ridden anyway
41     wordList ebct
42     (
43         psi_.boundaryField().size(),
44         zeroGradientFvPatchField<Type>::typeName
45     );
47     forAll (psi_.boundaryField(), patchI)
48     {
49         if (psi_.boundaryField()[patchI].fixesValue())
50         {
51             ebct[patchI] = fixedValueFvPatchField<Type>::typeName;
52         }
53     }
55     return ebct;
59 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
61 // Construct from components
62 template<class Type>
63 Foam::errorEstimate<Type>::errorEstimate
65     const GeometricField<Type, fvPatchField, volMesh>& psi,
66     const dimensionSet& ds,
67     const Field<Type>& res,
68     const scalarField& norm
71     psi_(psi),
72     dimensions_(ds),
73     residual_(res),
74     normFactor_(norm)
78 // Construct as copy
79 template<class Type>
80 Foam::errorEstimate<Type>::errorEstimate(const Foam::errorEstimate<Type>& ee)
82     refCount(),
83     psi_(ee.psi_),
84     dimensions_(ee.dimensions_),
85     residual_(ee.residual_),
86     normFactor_(ee.normFactor_)
90 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
92 template<class Type>
93 Foam::errorEstimate<Type>::~errorEstimate()
97 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
99 template<class Type>
100 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
101 Foam::errorEstimate<Type>::residual() const
103     tmp<GeometricField<Type, fvPatchField, volMesh> > tres
104     (
105         new GeometricField<Type, fvPatchField, volMesh>
106         (
107             IOobject
108             (
109                 "residual" + psi_.name(),
110                 psi_.mesh().time().timeName(),
111                 psi_.db(),
112                 IOobject::NO_READ,
113                 IOobject::NO_WRITE
114             ),
115             psi_.mesh(),
116             psi_.dimensions()/dimTime,
117             errorBCTypes()
118         )
119     );
121     GeometricField<Type, fvPatchField, volMesh>& res = tres();
123     res.internalField() = residual_;
124     res.boundaryField() == pTraits<Type>::zero;
126     res.correctBoundaryConditions();
128     return tres;
132 template<class Type>
133 Foam::tmp<Foam::volScalarField> Foam::errorEstimate<Type>::normFactor() const
135     tmp<volScalarField> tnormFactor
136     (
137         new volScalarField
138         (
139             IOobject
140             (
141                 "normFactor" + psi_.name(),
142                 psi_.mesh().time().timeName(),
143                 psi_.db(),
144                 IOobject::NO_READ,
145                 IOobject::NO_WRITE
146             ),
147             psi_.mesh(),
148             dimless/dimTime,
149             errorBCTypes()
150         )
151     );
153     volScalarField& normFactor = tnormFactor();
155     normFactor.internalField() = normFactor_;
156     normFactor.boundaryField() == pTraits<Type>::zero;
158     normFactor.correctBoundaryConditions();
160     return tnormFactor;
163 template<class Type>
164 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
165 Foam::errorEstimate<Type>::error() const
167     tmp<GeometricField<Type, fvPatchField, volMesh> > tresError
168     (
169         new GeometricField<Type, fvPatchField, volMesh>
170         (
171             IOobject
172             (
173                 "resError" + psi_.name(),
174                 psi_.mesh().time().timeName(),
175                 psi_.db(),
176                 IOobject::NO_READ,
177                 IOobject::NO_WRITE
178             ),
179             psi_.mesh(),
180             psi_.dimensions(),
181             errorBCTypes()
182         )
183     );
185     GeometricField<Type, fvPatchField, volMesh>& resError = tresError();
187     resError.internalField() = residual_/normFactor_;
188     resError.boundaryField() == pTraits<Type>::zero;
190     resError.correctBoundaryConditions();
192     return tresError;
195 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
197 template<class Type>
198 void Foam::errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>& rhs)
200     // Check for assignment to self
201     if (this == &rhs)
202     {
203         FatalErrorIn
204         (
205             "errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>&)"
206         )   << "Attempted assignment to self"
207             << abort(FatalError);
208     }
210     if (&psi_ != &(rhs.psi_))
211     {
212         FatalErrorIn
213         (
214             "errorEstimate<Type>::operator=(const errorEstimate<Type>&)"
215         )   << "different fields"
216             << abort(FatalError);
217     }
219     residual_ = rhs.residual_;
220     normFactor_ = rhs.normFactor_;
224 template<class Type>
225 void Foam::errorEstimate<Type>::operator=(const tmp<errorEstimate<Type> >& teev)
227     operator=(teev());
228     teev.clear();
232 template<class Type>
233 void Foam::errorEstimate<Type>::negate()
235     residual_.negate();
239 template<class Type>
240 void Foam::errorEstimate<Type>::operator+=(const errorEstimate<Type>& eev)
242     checkMethod(*this, eev, "+=");
244     dimensions_ += eev.dimensions_;
246     residual_ += eev.residual_;
247     normFactor_ += eev.normFactor_;
251 template<class Type>
252 void Foam::errorEstimate<Type>::operator+=
254     const tmp<errorEstimate<Type> >& teev
257     operator+=(teev());
258     teev.clear();
262 template<class Type>
263 void Foam::errorEstimate<Type>::operator-=(const errorEstimate<Type>& eev)
265     checkMethod(*this, eev, "+=");
267     dimensions_ -= eev.dimensions_;
268     residual_ -= eev.residual_;
269     normFactor_ += eev.normFactor_;
273 template<class Type>
274 void Foam::errorEstimate<Type>::operator-=(const tmp<errorEstimate<Type> >& teev)
276     operator-=(teev());
277     teev.clear();
281 template<class Type>
282 void Foam::errorEstimate<Type>::operator+=
284     const GeometricField<Type, fvPatchField, volMesh>& su
287     checkMethod(*this, su, "+=");
288     residual_ -= su.internalField();
292 template<class Type>
293 void Foam::errorEstimate<Type>::operator+=
295     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
298     operator+=(tsu());
299     tsu.clear();
303 template<class Type>
304 void Foam::errorEstimate<Type>::operator-=
306     const GeometricField<Type, fvPatchField, volMesh>& su
309     checkMethod(*this, su, "-=");
310     residual_ += su.internalField();
314 template<class Type>
315 void Foam::errorEstimate<Type>::operator-=
317     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
320     operator-=(tsu());
321     tsu.clear();
325 template<class Type>
326 void Foam::errorEstimate<Type>::operator+=
328     const dimensioned<Type>& su
331     residual_ -= su;
335 template<class Type>
336 void Foam::errorEstimate<Type>::operator-=
338     const dimensioned<Type>& su
341     residual_ += su;
345 template<class Type>
346 void Foam::errorEstimate<Type>::operator*=
348     const volScalarField& vsf
351     dimensions_ *= vsf.dimensions();
352     residual_ *= vsf.internalField();
353     normFactor_ *= vsf.internalField();
357 template<class Type>
358 void Foam::errorEstimate<Type>::operator*=
360     const tmp<volScalarField>& tvsf
363     operator*=(tvsf());
364     tvsf.clear();
368 template<class Type>
369 void Foam::errorEstimate<Type>::operator*=
371     const dimensioned<scalar>& ds
374     dimensions_ *= ds.dimensions();
375     residual_ *= ds.value();
376     normFactor_ *= ds.value();
380 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
382 template<class Type>
383 void Foam::checkMethod
385     const errorEstimate<Type>& ee1,
386     const errorEstimate<Type>& ee2,
387     const char* op
390     if (&ee1.psi() != &ee2.psi())
391     {
392         FatalErrorIn
393         (
394             "checkMethod(const errorEstimate<Type>&, "
395             "const errorEstimate<Type>&)"
396         )   << "incompatible fields for operation "
397             << endl << "    "
398             << "[" << ee1.psi().name() << "] "
399             << op
400             << " [" << ee2.psi().name() << "]"
401             << abort(FatalError);
402     }
404     if (dimensionSet::debug && ee1.dimensions() != ee2.dimensions())
405     {
406         FatalErrorIn
407         (
408             "checkMethod(const errorEstimate<Type>&, "
409             "const errorEstimate<Type>&)"
410         )   << "incompatible dimensions for operation "
411             << endl << "    "
412             << "[" << ee1.psi().name() << ee1.dimensions()/dimVolume << " ] "
413             << op
414             << " [" << ee2.psi().name() << ee2.dimensions()/dimVolume << " ]"
415             << abort(FatalError);
416     }
420 template<class Type>
421 void Foam::checkMethod
423     const errorEstimate<Type>& ee,
424     const GeometricField<Type, fvPatchField, volMesh>& vf,
425     const char* op
428     if (dimensionSet::debug && ee.dimensions()/dimVolume != vf.dimensions())
429     {
430         FatalErrorIn
431         (
432             "checkMethod(const errorEstimate<Type>&, "
433             "const GeometricField<Type, fvPatchField, volMesh>&)"
434         )   <<  "incompatible dimensions for operation "
435             << endl << "    "
436             << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
437             << op
438             << " [" << vf.name() << vf.dimensions() << " ]"
439             << abort(FatalError);
440     }
444 template<class Type>
445 void Foam::checkMethod
447     const errorEstimate<Type>& ee,
448     const dimensioned<Type>& dt,
449     const char* op
452     if (dimensionSet::debug && ee.dimensions()/dimVolume != dt.dimensions())
453     {
454         FatalErrorIn
455         (
456             "checkMethod(const errorEstimate<Type>&, const dimensioned<Type>&)"
457         )   << "incompatible dimensions for operation "
458             << endl << "    "
459             << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
460             << op
461             << " [" << dt.name() << dt.dimensions() << " ]"
462             << abort(FatalError);
463     }
467 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
469 namespace Foam
472 template<class Type>
473 tmp<errorEstimate<Type> > operator+
475     const errorEstimate<Type>& A,
476     const errorEstimate<Type>& B
479     checkMethod(A, B, "+");
480     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
481     tC() += B;
482     return tC;
486 template<class Type>
487 tmp<errorEstimate<Type> > operator+
489     const tmp<errorEstimate<Type> >& tA,
490     const errorEstimate<Type>& B
493     checkMethod(tA(), B, "+");
494     tmp<errorEstimate<Type> > tC(tA.ptr());
495     tC() += B;
496     return tC;
500 template<class Type>
501 tmp<errorEstimate<Type> > operator+
503     const errorEstimate<Type>& A,
504     const tmp<errorEstimate<Type> >& tB
507     checkMethod(A, tB(), "+");
508     tmp<errorEstimate<Type> > tC(tB.ptr());
509     tC() += A;
510     return tC;
514 template<class Type>
515 tmp<errorEstimate<Type> > operator+
517     const tmp<errorEstimate<Type> >& tA,
518     const tmp<errorEstimate<Type> >& tB
521     checkMethod(tA(), tB(), "+");
522     tmp<errorEstimate<Type> > tC(tA.ptr());
523     tC() += tB();
524     tB.clear();
525     return tC;
529 template<class Type>
530 tmp<errorEstimate<Type> > operator-
532     const errorEstimate<Type>& A
535     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
536     tC().negate();
537     return tC;
541 template<class Type>
542 tmp<errorEstimate<Type> > operator-
544     const tmp<errorEstimate<Type> >& tA
547     tmp<errorEstimate<Type> > tC(tA.ptr());
548     tC().negate();
549     return tC;
553 template<class Type>
554 tmp<errorEstimate<Type> > operator-
556     const errorEstimate<Type>& A,
557     const errorEstimate<Type>& B
560     checkMethod(A, B, "-");
561     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
562     tC() -= B;
563     return tC;
567 template<class Type>
568 tmp<errorEstimate<Type> > operator-
570     const tmp<errorEstimate<Type> >& tA,
571     const errorEstimate<Type>& B
574     checkMethod(tA(), B, "-");
575     tmp<errorEstimate<Type> > tC(tA.ptr());
576     tC() -= B;
577     return tC;
581 template<class Type>
582 tmp<errorEstimate<Type> > operator-
584     const errorEstimate<Type>& A,
585     const tmp<errorEstimate<Type> >& tB
588     checkMethod(A, tB(), "-");
589     tmp<errorEstimate<Type> > tC(tB.ptr());
590     tC() -= A;
591     tC().negate();
592     return tC;
596 template<class Type>
597 tmp<errorEstimate<Type> > operator-
599     const tmp<errorEstimate<Type> >& tA,
600     const tmp<errorEstimate<Type> >& tB
603     checkMethod(tA(), tB(), "-");
604     tmp<errorEstimate<Type> > tC(tA.ptr());
605     tC() -= tB();
606     tB.clear();
607     return tC;
611 template<class Type>
612 tmp<errorEstimate<Type> > operator==
614     const errorEstimate<Type>& A,
615     const errorEstimate<Type>& B
618     checkMethod(A, B, "==");
619     return (A - B);
623 template<class Type>
624 tmp<errorEstimate<Type> > operator==
626     const tmp<errorEstimate<Type> >& tA,
627     const errorEstimate<Type>& B
630     checkMethod(tA(), B, "==");
631     return (tA - B);
635 template<class Type>
636 tmp<errorEstimate<Type> > operator==
638     const errorEstimate<Type>& A,
639     const tmp<errorEstimate<Type> >& tB
642     checkMethod(A, tB(), "==");
643     return (A - tB);
647 template<class Type>
648 tmp<errorEstimate<Type> > operator==
650     const tmp<errorEstimate<Type> >& tA,
651     const tmp<errorEstimate<Type> >& tB
654     checkMethod(tA(), tB(), "==");
655     return (tA - tB);
659 template<class Type>
660 tmp<errorEstimate<Type> > operator+
662     const errorEstimate<Type>& A,
663     const GeometricField<Type, fvPatchField, volMesh>& su
666     checkMethod(A, su, "+");
667     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
668     tC().res() -= su.internalField();
669     return tC;
672 template<class Type>
673 tmp<errorEstimate<Type> > operator+
675     const tmp<errorEstimate<Type> >& tA,
676     const GeometricField<Type, fvPatchField, volMesh>& su
679     checkMethod(tA(), su, "+");
680     tmp<errorEstimate<Type> > tC(tA.ptr());
681     tC().res() -= su.internalField();
682     return tC;
685 template<class Type>
686 tmp<errorEstimate<Type> > operator+
688     const errorEstimate<Type>& A,
689     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
692     checkMethod(A, tsu(), "+");
693     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
694     tC().res() -= tsu().internalField();
695     tsu.clear();
696     return tC;
700 template<class Type>
701 tmp<errorEstimate<Type> > operator+
703     const tmp<errorEstimate<Type> >& tA,
704     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
707     checkMethod(tA(), tsu(), "+");
708     tmp<errorEstimate<Type> > tC(tA.ptr());
709     tC().res() -= tsu().internalField();
710     tsu.clear();
711     return tC;
714 template<class Type>
715 tmp<errorEstimate<Type> > operator+
717     const GeometricField<Type, fvPatchField, volMesh>& su,
718     const errorEstimate<Type>& A
721     checkMethod(A, su, "+");
722     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
723     tC().res() -= su.internalField();
724     return tC;
727 template<class Type>
728 tmp<errorEstimate<Type> > operator+
730     const GeometricField<Type, fvPatchField, volMesh>& su,
731     const tmp<errorEstimate<Type> >& tA
734     checkMethod(tA(), su, "+");
735     tmp<errorEstimate<Type> > tC(tA.ptr());
736     tC().res() -= su.internalField();
737     return tC;
740 template<class Type>
741 tmp<errorEstimate<Type> > operator+
743     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
744     const errorEstimate<Type>& A
747     checkMethod(A, tsu(), "+");
748     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
749     tC().res() -= tsu().internalField();
750     tsu.clear();
751     return tC;
754 template<class Type>
755 tmp<errorEstimate<Type> > operator+
757     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
758     const tmp<errorEstimate<Type> >& tA
761     checkMethod(tA(), tsu(), "+");
762     tmp<errorEstimate<Type> > tC(tA.ptr());
763     tC().res() -= tsu().internalField();
764     tsu.clear();
765     return tC;
769 template<class Type>
770 tmp<errorEstimate<Type> > operator-
772     const errorEstimate<Type>& A,
773     const GeometricField<Type, fvPatchField, volMesh>& su
776     checkMethod(A, su, "-");
777     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
778     tC().res() += su.internalField();
779     return tC;
782 template<class Type>
783 tmp<errorEstimate<Type> > operator-
785     const tmp<errorEstimate<Type> >& tA,
786     const GeometricField<Type, fvPatchField, volMesh>& su
789     checkMethod(tA(), su, "-");
790     tmp<errorEstimate<Type> > tC(tA.ptr());
791     tC().res() += su.internalField();
792     return tC;
795 template<class Type>
796 tmp<errorEstimate<Type> > operator-
798     const errorEstimate<Type>& A,
799     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
802     checkMethod(A, tsu(), "-");
803     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
804     tC().res() += tsu().internalField();
805     tsu.clear();
806     return tC;
809 template<class Type>
810 tmp<errorEstimate<Type> > operator-
812     const tmp<errorEstimate<Type> >& tA,
813     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
816     checkMethod(tA(), tsu(), "-");
817     tmp<errorEstimate<Type> > tC(tA.ptr());
818     tC().res() += tsu().internalField();
819     tsu.clear();
820     return tC;
824 template<class Type>
825 tmp<errorEstimate<Type> > operator-
827     const GeometricField<Type, fvPatchField, volMesh>& su,
828     const errorEstimate<Type>& A
831     checkMethod(A, su, "-");
832     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
833     tC().negate();
834     tC().res() -= su.internalField();
835     return tC;
839 template<class Type>
840 tmp<errorEstimate<Type> > operator-
842     const GeometricField<Type, fvPatchField, volMesh>& su,
843     const tmp<errorEstimate<Type> >& tA
846     checkMethod(tA(), su, "-");
847     tmp<errorEstimate<Type> > tC(tA.ptr());
848     tC().negate();
849     tC().res() -= su.internalField();
850     return tC;
853 template<class Type>
854 tmp<errorEstimate<Type> > operator-
856     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
857     const errorEstimate<Type>& A
860     checkMethod(A, tsu(), "-");
861     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
862     tC().negate();
863     tC().res() -= tsu().internalField();
864     tsu.clear();
865     return tC;
869 template<class Type>
870 tmp<errorEstimate<Type> > operator-
872     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
873     const tmp<errorEstimate<Type> >& tA
876     checkMethod(tA(), tsu(), "-");
877     tmp<errorEstimate<Type> > tC(tA.ptr());
878     tC().negate();
879     tC().res() -= tsu().internalField();
880     tsu.clear();
881     return tC;
885 template<class Type>
886 tmp<errorEstimate<Type> > operator+
888     const errorEstimate<Type>& A,
889     const dimensioned<Type>& su
892     checkMethod(A, su, "+");
893     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
894     tC().res() -= su.value();
895     return tC;
899 template<class Type>
900 tmp<errorEstimate<Type> > operator+
902     const tmp<errorEstimate<Type> >& tA,
903     const dimensioned<Type>& su
906     checkMethod(tA(), su, "+");
907     tmp<errorEstimate<Type> > tC(tA.ptr());
908     tC().res() -= su.value();
909     return tC;
913 template<class Type>
914 tmp<errorEstimate<Type> > operator+
916     const dimensioned<Type>& su,
917     const errorEstimate<Type>& A
920     checkMethod(A, su, "+");
921     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
922     tC().res() -= su.value();
923     return tC;
927 template<class Type>
928 tmp<errorEstimate<Type> > operator+
930     const dimensioned<Type>& su,
931     const tmp<errorEstimate<Type> >& tA
934     checkMethod(tA(), su, "+");
935     tmp<errorEstimate<Type> > tC(tA.ptr());
936     tC().res() -= su.value();
937     return tC;
941 template<class Type>
942 tmp<errorEstimate<Type> > operator-
944     const errorEstimate<Type>& A,
945     const dimensioned<Type>& su
948     checkMethod(A, su, "-");
949     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
950     tC().res() += su.value();
951     return tC;
955 template<class Type>
956 tmp<errorEstimate<Type> > operator-
958     const tmp<errorEstimate<Type> >& tA,
959     const dimensioned<Type>& su
962     checkMethod(tA(), su, "-");
963     tmp<errorEstimate<Type> > tC(tA.ptr());
964     tC().res() += su.value();
965     return tC;
969 template<class Type>
970 tmp<errorEstimate<Type> > operator-
972     const dimensioned<Type>& su,
973     const errorEstimate<Type>& A
976     checkMethod(A, su, "-");
977     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
978     tC().negate();
979     tC().res() -= su.value();
980     return tC;
984 template<class Type>
985 tmp<errorEstimate<Type> > operator-
987     const dimensioned<Type>& su,
988     const tmp<errorEstimate<Type> >& tA
991     checkMethod(tA(), su, "-");
992     tmp<errorEstimate<Type> > tC(tA.ptr());
993     tC().negate();
994     tC().res() -= su.value();
995     return tC;
999 template<class Type>
1000 tmp<errorEstimate<Type> > operator==
1002     const errorEstimate<Type>& A,
1003     const GeometricField<Type, fvPatchField, volMesh>& su
1006     checkMethod(A, su, "==");
1007     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1008     tC().res() += su.internalField();
1009     return tC;
1012 template<class Type>
1013 tmp<errorEstimate<Type> > operator==
1015     const tmp<errorEstimate<Type> >& tA,
1016     const GeometricField<Type, fvPatchField, volMesh>& su
1019     checkMethod(tA(), su, "==");
1020     tmp<errorEstimate<Type> > tC(tA.ptr());
1021     tC().res() += su.internalField();
1022     return tC;
1025 template<class Type>
1026 tmp<errorEstimate<Type> > operator==
1028     const errorEstimate<Type>& A,
1029     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1032     checkMethod(A, tsu(), "==");
1033     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1034     tC().res() += tsu().internalField();
1035     tsu.clear();
1036     return tC;
1039 template<class Type>
1040 tmp<errorEstimate<Type> > operator==
1042     const tmp<errorEstimate<Type> >& tA,
1043     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1046     checkMethod(tA(), tsu(), "==");
1047     tmp<errorEstimate<Type> > tC(tA.ptr());
1048     tC().res() += tsu().internalField();
1049     tsu.clear();
1050     return tC;
1054 template<class Type>
1055 tmp<errorEstimate<Type> > operator==
1057     const errorEstimate<Type>& A,
1058     const dimensioned<Type>& su
1061     checkMethod(A, su, "==");
1062     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1063     tC().res() += su.value();
1064     return tC;
1068 template<class Type>
1069 tmp<errorEstimate<Type> > operator==
1071     const tmp<errorEstimate<Type> >& tA,
1072     const dimensioned<Type>& su
1075     checkMethod(tA(), su, "==");
1076     tmp<errorEstimate<Type> > tC(tA.ptr());
1077     tC().res() += su.value();
1078     return tC;
1082 template<class Type>
1083 tmp<errorEstimate<Type> > operator*
1085     const volScalarField& vsf,
1086     const errorEstimate<Type>& A
1089     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1090     tC() *= vsf;
1091     return tC;
1094 template<class Type>
1095 tmp<errorEstimate<Type> > operator*
1097     const tmp<volScalarField>& tvsf,
1098     const errorEstimate<Type>& A
1101     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1102     tC() *= tvsf;
1103     return tC;
1106 template<class Type>
1107 tmp<errorEstimate<Type> > operator*
1109     const volScalarField& vsf,
1110     const tmp<errorEstimate<Type> >& tA
1113     tmp<errorEstimate<Type> > tC(tA.ptr());
1114     tC() *= vsf;
1115     return tC;
1118 template<class Type>
1119 tmp<errorEstimate<Type> > operator*
1121     const tmp<volScalarField>& tvsf,
1122     const tmp<errorEstimate<Type> >& tA
1125     tmp<errorEstimate<Type> > tC(tA.ptr());
1126     tC() *= tvsf;
1127     return tC;
1131 template<class Type>
1132 tmp<errorEstimate<Type> > operator*
1134     const dimensioned<scalar>& ds,
1135     const errorEstimate<Type>& A
1138     tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
1139     tC() *= ds;
1140     return tC;
1144 template<class Type>
1145 tmp<errorEstimate<Type> > operator*
1147     const dimensioned<scalar>& ds,
1148     const tmp<errorEstimate<Type> >& tA
1151     tmp<errorEstimate<Type> > tC(tA.ptr());
1152     tC() *= ds;
1153     return tC;
1157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1159 } // End namespace Foam
1161 // ************************************************************************* //