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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "errorEstimate.H"
27 #include "zeroGradientFvPatchField.H"
28 #include "fixedValueFvPatchField.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 Foam::wordList Foam::errorEstimate<Type>::errorBCTypes() const
38 // Make the boundary condition type list
39 // Default types get over-ridden anyway
42 psi_.boundaryField().size(),
43 zeroGradientFvPatchField<Type>::typeName
46 forAll (psi_.boundaryField(), patchI)
48 if (psi_.boundaryField()[patchI].fixesValue())
50 ebct[patchI] = fixedValueFvPatchField<Type>::typeName;
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 // Construct from components
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
79 Foam::errorEstimate<Type>::errorEstimate(const Foam::errorEstimate<Type>& ee)
83 dimensions_(ee.dimensions_),
84 residual_(ee.residual_),
85 normFactor_(ee.normFactor_)
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 Foam::errorEstimate<Type>::~errorEstimate()
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
100 Foam::errorEstimate<Type>::residual() const
102 tmp<GeometricField<Type, fvPatchField, volMesh> > tres
104 new GeometricField<Type, fvPatchField, volMesh>
108 "residual" + psi_.name(),
109 psi_.mesh().time().timeName(),
115 psi_.dimensions()/dimTime,
120 GeometricField<Type, fvPatchField, volMesh>& res = tres();
122 res.internalField() = residual_;
123 res.boundaryField() == pTraits<Type>::zero;
125 res.correctBoundaryConditions();
132 Foam::tmp<Foam::volScalarField> Foam::errorEstimate<Type>::normFactor() const
134 tmp<volScalarField> tnormFactor
140 "normFactor" + psi_.name(),
141 psi_.mesh().time().timeName(),
152 volScalarField& normFactor = tnormFactor();
154 normFactor.internalField() = normFactor_;
155 normFactor.boundaryField() == pTraits<Type>::zero;
157 normFactor.correctBoundaryConditions();
163 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
164 Foam::errorEstimate<Type>::error() const
166 tmp<GeometricField<Type, fvPatchField, volMesh> > tresError
168 new GeometricField<Type, fvPatchField, volMesh>
172 "resError" + psi_.name(),
173 psi_.mesh().time().timeName(),
184 GeometricField<Type, fvPatchField, volMesh>& resError = tresError();
186 resError.internalField() = residual_/normFactor_;
187 resError.boundaryField() == pTraits<Type>::zero;
189 resError.correctBoundaryConditions();
194 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
197 void Foam::errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>& rhs)
199 // Check for assignment to self
204 "errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>&)"
205 ) << "Attempted assignment to self"
206 << abort(FatalError);
209 if (&psi_ != &(rhs.psi_))
213 "errorEstimate<Type>::operator=(const errorEstimate<Type>&)"
214 ) << "different fields"
215 << abort(FatalError);
218 residual_ = rhs.residual_;
219 normFactor_ = rhs.normFactor_;
224 void Foam::errorEstimate<Type>::operator=(const tmp<errorEstimate<Type> >& teev)
232 void Foam::errorEstimate<Type>::negate()
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_;
251 void Foam::errorEstimate<Type>::operator+=
253 const tmp<errorEstimate<Type> >& teev
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_;
273 void Foam::errorEstimate<Type>::operator-=(const tmp<errorEstimate<Type> >& teev)
281 void Foam::errorEstimate<Type>::operator+=
283 const GeometricField<Type, fvPatchField, volMesh>& su
286 checkMethod(*this, su, "+=");
287 residual_ -= su.internalField();
292 void Foam::errorEstimate<Type>::operator+=
294 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
303 void Foam::errorEstimate<Type>::operator-=
305 const GeometricField<Type, fvPatchField, volMesh>& su
308 checkMethod(*this, su, "-=");
309 residual_ += su.internalField();
314 void Foam::errorEstimate<Type>::operator-=
316 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
325 void Foam::errorEstimate<Type>::operator+=
327 const dimensioned<Type>& su
335 void Foam::errorEstimate<Type>::operator-=
337 const dimensioned<Type>& su
345 void Foam::errorEstimate<Type>::operator*=
347 const volScalarField& vsf
350 dimensions_ *= vsf.dimensions();
351 residual_ *= vsf.internalField();
352 normFactor_ *= vsf.internalField();
357 void Foam::errorEstimate<Type>::operator*=
359 const tmp<volScalarField>& tvsf
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 * * * * * * * * * * * * * //
382 void Foam::checkMethod
384 const errorEstimate<Type>& ee1,
385 const errorEstimate<Type>& ee2,
389 if (&ee1.psi() != &ee2.psi())
393 "checkMethod(const errorEstimate<Type>&, "
394 "const errorEstimate<Type>&)"
395 ) << "incompatible fields for operation "
397 << "[" << ee1.psi().name() << "] "
399 << " [" << ee2.psi().name() << "]"
400 << abort(FatalError);
403 if (dimensionSet::debug && ee1.dimensions() != ee2.dimensions())
407 "checkMethod(const errorEstimate<Type>&, "
408 "const errorEstimate<Type>&)"
409 ) << "incompatible dimensions for operation "
411 << "[" << ee1.psi().name() << ee1.dimensions()/dimVolume << " ] "
413 << " [" << ee2.psi().name() << ee2.dimensions()/dimVolume << " ]"
414 << abort(FatalError);
420 void Foam::checkMethod
422 const errorEstimate<Type>& ee,
423 const GeometricField<Type, fvPatchField, volMesh>& vf,
427 if (dimensionSet::debug && ee.dimensions()/dimVolume != vf.dimensions())
431 "checkMethod(const errorEstimate<Type>&, "
432 "const GeometricField<Type, fvPatchField, volMesh>&)"
433 ) << "incompatible dimensions for operation "
435 << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
437 << " [" << vf.name() << vf.dimensions() << " ]"
438 << abort(FatalError);
444 void Foam::checkMethod
446 const errorEstimate<Type>& ee,
447 const dimensioned<Type>& dt,
451 if (dimensionSet::debug && ee.dimensions()/dimVolume != dt.dimensions())
455 "checkMethod(const errorEstimate<Type>&, const dimensioned<Type>&)"
456 ) << "incompatible dimensions for operation "
458 << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
460 << " [" << dt.name() << dt.dimensions() << " ]"
461 << abort(FatalError);
466 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
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));
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());
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());
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());
529 tmp<errorEstimate<Type> > operator-
531 const errorEstimate<Type>& A
534 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
541 tmp<errorEstimate<Type> > operator-
543 const tmp<errorEstimate<Type> >& tA
546 tmp<errorEstimate<Type> > tC(tA.ptr());
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));
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());
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());
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());
611 tmp<errorEstimate<Type> > operator==
613 const errorEstimate<Type>& A,
614 const errorEstimate<Type>& B
617 checkMethod(A, B, "==");
623 tmp<errorEstimate<Type> > operator==
625 const tmp<errorEstimate<Type> >& tA,
626 const errorEstimate<Type>& B
629 checkMethod(tA(), B, "==");
635 tmp<errorEstimate<Type> > operator==
637 const errorEstimate<Type>& A,
638 const tmp<errorEstimate<Type> >& tB
641 checkMethod(A, tB(), "==");
647 tmp<errorEstimate<Type> > operator==
649 const tmp<errorEstimate<Type> >& tA,
650 const tmp<errorEstimate<Type> >& tB
653 checkMethod(tA(), tB(), "==");
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();
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();
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();
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();
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();
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();
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();
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();
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();
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();
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();
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();
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));
833 tC().res() -= su.internalField();
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());
848 tC().res() -= su.internalField();
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));
862 tC().res() -= tsu().internalField();
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());
878 tC().res() -= tsu().internalField();
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();
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();
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();
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();
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();
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();
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));
978 tC().res() -= su.value();
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());
993 tC().res() -= su.value();
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();
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();
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();
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();
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();
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();
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));
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));
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());
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());
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));
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());
1156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1158 } // End namespace Foam
1160 // ************************************************************************* //