1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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 * * * * * * * * * * * //
37 Foam::wordList Foam::errorEstimate<Type>::errorBCTypes() const
39 // Make the boundary condition type list
40 // Default types get over-ridden anyway
43 psi_.boundaryField().size(),
44 zeroGradientFvPatchField<Type>::typeName
47 forAll (psi_.boundaryField(), patchI)
49 if (psi_.boundaryField()[patchI].fixesValue())
51 ebct[patchI] = fixedValueFvPatchField<Type>::typeName;
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 // Construct from components
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
80 Foam::errorEstimate<Type>::errorEstimate(const Foam::errorEstimate<Type>& ee)
84 dimensions_(ee.dimensions_),
85 residual_(ee.residual_),
86 normFactor_(ee.normFactor_)
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 Foam::errorEstimate<Type>::~errorEstimate()
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
101 Foam::errorEstimate<Type>::residual() const
103 tmp<GeometricField<Type, fvPatchField, volMesh> > tres
105 new GeometricField<Type, fvPatchField, volMesh>
109 "residual" + psi_.name(),
110 psi_.mesh().time().timeName(),
116 psi_.dimensions()/dimTime,
121 GeometricField<Type, fvPatchField, volMesh>& res = tres();
123 res.internalField() = residual_;
124 res.boundaryField() == pTraits<Type>::zero;
126 res.correctBoundaryConditions();
133 Foam::tmp<Foam::volScalarField> Foam::errorEstimate<Type>::normFactor() const
135 tmp<volScalarField> tnormFactor
141 "normFactor" + psi_.name(),
142 psi_.mesh().time().timeName(),
153 volScalarField& normFactor = tnormFactor();
155 normFactor.internalField() = normFactor_;
156 normFactor.boundaryField() == pTraits<Type>::zero;
158 normFactor.correctBoundaryConditions();
164 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
165 Foam::errorEstimate<Type>::error() const
167 tmp<GeometricField<Type, fvPatchField, volMesh> > tresError
169 new GeometricField<Type, fvPatchField, volMesh>
173 "resError" + psi_.name(),
174 psi_.mesh().time().timeName(),
185 GeometricField<Type, fvPatchField, volMesh>& resError = tresError();
187 resError.internalField() = residual_/normFactor_;
188 resError.boundaryField() == pTraits<Type>::zero;
190 resError.correctBoundaryConditions();
195 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
198 void Foam::errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>& rhs)
200 // Check for assignment to self
205 "errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>&)"
206 ) << "Attempted assignment to self"
207 << abort(FatalError);
210 if (&psi_ != &(rhs.psi_))
214 "errorEstimate<Type>::operator=(const errorEstimate<Type>&)"
215 ) << "different fields"
216 << abort(FatalError);
219 residual_ = rhs.residual_;
220 normFactor_ = rhs.normFactor_;
225 void Foam::errorEstimate<Type>::operator=(const tmp<errorEstimate<Type> >& teev)
233 void Foam::errorEstimate<Type>::negate()
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_;
252 void Foam::errorEstimate<Type>::operator+=
254 const tmp<errorEstimate<Type> >& teev
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_;
274 void Foam::errorEstimate<Type>::operator-=(const tmp<errorEstimate<Type> >& teev)
282 void Foam::errorEstimate<Type>::operator+=
284 const GeometricField<Type, fvPatchField, volMesh>& su
287 checkMethod(*this, su, "+=");
288 residual_ -= su.internalField();
293 void Foam::errorEstimate<Type>::operator+=
295 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
304 void Foam::errorEstimate<Type>::operator-=
306 const GeometricField<Type, fvPatchField, volMesh>& su
309 checkMethod(*this, su, "-=");
310 residual_ += su.internalField();
315 void Foam::errorEstimate<Type>::operator-=
317 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
326 void Foam::errorEstimate<Type>::operator+=
328 const dimensioned<Type>& su
336 void Foam::errorEstimate<Type>::operator-=
338 const dimensioned<Type>& su
346 void Foam::errorEstimate<Type>::operator*=
348 const volScalarField& vsf
351 dimensions_ *= vsf.dimensions();
352 residual_ *= vsf.internalField();
353 normFactor_ *= vsf.internalField();
358 void Foam::errorEstimate<Type>::operator*=
360 const tmp<volScalarField>& tvsf
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 * * * * * * * * * * * * * //
383 void Foam::checkMethod
385 const errorEstimate<Type>& ee1,
386 const errorEstimate<Type>& ee2,
390 if (&ee1.psi() != &ee2.psi())
394 "checkMethod(const errorEstimate<Type>&, "
395 "const errorEstimate<Type>&)"
396 ) << "incompatible fields for operation "
398 << "[" << ee1.psi().name() << "] "
400 << " [" << ee2.psi().name() << "]"
401 << abort(FatalError);
404 if (dimensionSet::debug && ee1.dimensions() != ee2.dimensions())
408 "checkMethod(const errorEstimate<Type>&, "
409 "const errorEstimate<Type>&)"
410 ) << "incompatible dimensions for operation "
412 << "[" << ee1.psi().name() << ee1.dimensions()/dimVolume << " ] "
414 << " [" << ee2.psi().name() << ee2.dimensions()/dimVolume << " ]"
415 << abort(FatalError);
421 void Foam::checkMethod
423 const errorEstimate<Type>& ee,
424 const GeometricField<Type, fvPatchField, volMesh>& vf,
428 if (dimensionSet::debug && ee.dimensions()/dimVolume != vf.dimensions())
432 "checkMethod(const errorEstimate<Type>&, "
433 "const GeometricField<Type, fvPatchField, volMesh>&)"
434 ) << "incompatible dimensions for operation "
436 << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
438 << " [" << vf.name() << vf.dimensions() << " ]"
439 << abort(FatalError);
445 void Foam::checkMethod
447 const errorEstimate<Type>& ee,
448 const dimensioned<Type>& dt,
452 if (dimensionSet::debug && ee.dimensions()/dimVolume != dt.dimensions())
456 "checkMethod(const errorEstimate<Type>&, const dimensioned<Type>&)"
457 ) << "incompatible dimensions for operation "
459 << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
461 << " [" << dt.name() << dt.dimensions() << " ]"
462 << abort(FatalError);
467 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
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));
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());
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());
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());
530 tmp<errorEstimate<Type> > operator-
532 const errorEstimate<Type>& A
535 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
542 tmp<errorEstimate<Type> > operator-
544 const tmp<errorEstimate<Type> >& tA
547 tmp<errorEstimate<Type> > tC(tA.ptr());
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));
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());
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());
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());
612 tmp<errorEstimate<Type> > operator==
614 const errorEstimate<Type>& A,
615 const errorEstimate<Type>& B
618 checkMethod(A, B, "==");
624 tmp<errorEstimate<Type> > operator==
626 const tmp<errorEstimate<Type> >& tA,
627 const errorEstimate<Type>& B
630 checkMethod(tA(), B, "==");
636 tmp<errorEstimate<Type> > operator==
638 const errorEstimate<Type>& A,
639 const tmp<errorEstimate<Type> >& tB
642 checkMethod(A, tB(), "==");
648 tmp<errorEstimate<Type> > operator==
650 const tmp<errorEstimate<Type> >& tA,
651 const tmp<errorEstimate<Type> >& tB
654 checkMethod(tA(), tB(), "==");
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();
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();
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();
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();
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();
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();
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();
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();
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();
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();
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();
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();
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));
834 tC().res() -= su.internalField();
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());
849 tC().res() -= su.internalField();
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));
863 tC().res() -= tsu().internalField();
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());
879 tC().res() -= tsu().internalField();
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();
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();
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();
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();
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();
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();
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));
979 tC().res() -= su.value();
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());
994 tC().res() -= su.value();
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();
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();
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();
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();
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();
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();
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));
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));
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());
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());
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));
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());
1157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1159 } // End namespace Foam
1161 // ************************************************************************* //