00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "errorEstimate.H"
00027 #include <finiteVolume/zeroGradientFvPatchField.H>
00028 #include <finiteVolume/fixedValueFvPatchField.H>
00029
00030
00031
00032
00033
00034
00035 template<class Type>
00036 Foam::wordList Foam::errorEstimate<Type>::errorBCTypes() const
00037 {
00038
00039
00040 wordList ebct
00041 (
00042 psi_.boundaryField().size(),
00043 zeroGradientFvPatchField<Type>::typeName
00044 );
00045
00046 forAll (psi_.boundaryField(), patchI)
00047 {
00048 if (psi_.boundaryField()[patchI].fixesValue())
00049 {
00050 ebct[patchI] = fixedValueFvPatchField<Type>::typeName;
00051 }
00052 }
00053
00054 return ebct;
00055 }
00056
00057
00058
00059
00060
00061 template<class Type>
00062 Foam::errorEstimate<Type>::errorEstimate
00063 (
00064 const GeometricField<Type, fvPatchField, volMesh>& psi,
00065 const dimensionSet& ds,
00066 const Field<Type>& res,
00067 const scalarField& norm
00068 )
00069 :
00070 psi_(psi),
00071 dimensions_(ds),
00072 residual_(res),
00073 normFactor_(norm)
00074 {}
00075
00076
00077
00078 template<class Type>
00079 Foam::errorEstimate<Type>::errorEstimate(const Foam::errorEstimate<Type>& ee)
00080 :
00081 refCount(),
00082 psi_(ee.psi_),
00083 dimensions_(ee.dimensions_),
00084 residual_(ee.residual_),
00085 normFactor_(ee.normFactor_)
00086 {}
00087
00088
00089
00090
00091 template<class Type>
00092 Foam::errorEstimate<Type>::~errorEstimate()
00093 {}
00094
00095
00096
00097
00098 template<class Type>
00099 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
00100 Foam::errorEstimate<Type>::residual() const
00101 {
00102 tmp<GeometricField<Type, fvPatchField, volMesh> > tres
00103 (
00104 new GeometricField<Type, fvPatchField, volMesh>
00105 (
00106 IOobject
00107 (
00108 "residual" + psi_.name(),
00109 psi_.mesh().time().timeName(),
00110 psi_.db(),
00111 IOobject::NO_READ,
00112 IOobject::NO_WRITE
00113 ),
00114 psi_.mesh(),
00115 psi_.dimensions()/dimTime,
00116 errorBCTypes()
00117 )
00118 );
00119
00120 GeometricField<Type, fvPatchField, volMesh>& res = tres();
00121
00122 res.internalField() = residual_;
00123 res.boundaryField() == pTraits<Type>::zero;
00124
00125 res.correctBoundaryConditions();
00126
00127 return tres;
00128 }
00129
00130
00131 template<class Type>
00132 Foam::tmp<Foam::volScalarField> Foam::errorEstimate<Type>::normFactor() const
00133 {
00134 tmp<volScalarField> tnormFactor
00135 (
00136 new volScalarField
00137 (
00138 IOobject
00139 (
00140 "normFactor" + psi_.name(),
00141 psi_.mesh().time().timeName(),
00142 psi_.db(),
00143 IOobject::NO_READ,
00144 IOobject::NO_WRITE
00145 ),
00146 psi_.mesh(),
00147 dimless/dimTime,
00148 errorBCTypes()
00149 )
00150 );
00151
00152 volScalarField& normFactor = tnormFactor();
00153
00154 normFactor.internalField() = normFactor_;
00155 normFactor.boundaryField() == pTraits<Type>::zero;
00156
00157 normFactor.correctBoundaryConditions();
00158
00159 return tnormFactor;
00160 }
00161
00162 template<class Type>
00163 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
00164 Foam::errorEstimate<Type>::error() const
00165 {
00166 tmp<GeometricField<Type, fvPatchField, volMesh> > tresError
00167 (
00168 new GeometricField<Type, fvPatchField, volMesh>
00169 (
00170 IOobject
00171 (
00172 "resError" + psi_.name(),
00173 psi_.mesh().time().timeName(),
00174 psi_.db(),
00175 IOobject::NO_READ,
00176 IOobject::NO_WRITE
00177 ),
00178 psi_.mesh(),
00179 psi_.dimensions(),
00180 errorBCTypes()
00181 )
00182 );
00183
00184 GeometricField<Type, fvPatchField, volMesh>& resError = tresError();
00185
00186 resError.internalField() = residual_/normFactor_;
00187 resError.boundaryField() == pTraits<Type>::zero;
00188
00189 resError.correctBoundaryConditions();
00190
00191 return tresError;
00192 }
00193
00194
00195
00196 template<class Type>
00197 void Foam::errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>& rhs)
00198 {
00199
00200 if (this == &rhs)
00201 {
00202 FatalErrorIn
00203 (
00204 "errorEstimate<Type>::operator=(const Foam::errorEstimate<Type>&)"
00205 ) << "Attempted assignment to self"
00206 << abort(FatalError);
00207 }
00208
00209 if (&psi_ != &(rhs.psi_))
00210 {
00211 FatalErrorIn
00212 (
00213 "errorEstimate<Type>::operator=(const errorEstimate<Type>&)"
00214 ) << "different fields"
00215 << abort(FatalError);
00216 }
00217
00218 residual_ = rhs.residual_;
00219 normFactor_ = rhs.normFactor_;
00220 }
00221
00222
00223 template<class Type>
00224 void Foam::errorEstimate<Type>::operator=(const tmp<errorEstimate<Type> >& teev)
00225 {
00226 operator=(teev());
00227 teev.clear();
00228 }
00229
00230
00231 template<class Type>
00232 void Foam::errorEstimate<Type>::negate()
00233 {
00234 residual_.negate();
00235 }
00236
00237
00238 template<class Type>
00239 void Foam::errorEstimate<Type>::operator+=(const errorEstimate<Type>& eev)
00240 {
00241 checkMethod(*this, eev, "+=");
00242
00243 dimensions_ += eev.dimensions_;
00244
00245 residual_ += eev.residual_;
00246 normFactor_ += eev.normFactor_;
00247 }
00248
00249
00250 template<class Type>
00251 void Foam::errorEstimate<Type>::operator+=
00252 (
00253 const tmp<errorEstimate<Type> >& teev
00254 )
00255 {
00256 operator+=(teev());
00257 teev.clear();
00258 }
00259
00260
00261 template<class Type>
00262 void Foam::errorEstimate<Type>::operator-=(const errorEstimate<Type>& eev)
00263 {
00264 checkMethod(*this, eev, "+=");
00265
00266 dimensions_ -= eev.dimensions_;
00267 residual_ -= eev.residual_;
00268 normFactor_ += eev.normFactor_;
00269 }
00270
00271
00272 template<class Type>
00273 void Foam::errorEstimate<Type>::operator-=(const tmp<errorEstimate<Type> >& teev)
00274 {
00275 operator-=(teev());
00276 teev.clear();
00277 }
00278
00279
00280 template<class Type>
00281 void Foam::errorEstimate<Type>::operator+=
00282 (
00283 const GeometricField<Type, fvPatchField, volMesh>& su
00284 )
00285 {
00286 checkMethod(*this, su, "+=");
00287 residual_ -= su.internalField();
00288 }
00289
00290
00291 template<class Type>
00292 void Foam::errorEstimate<Type>::operator+=
00293 (
00294 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
00295 )
00296 {
00297 operator+=(tsu());
00298 tsu.clear();
00299 }
00300
00301
00302 template<class Type>
00303 void Foam::errorEstimate<Type>::operator-=
00304 (
00305 const GeometricField<Type, fvPatchField, volMesh>& su
00306 )
00307 {
00308 checkMethod(*this, su, "-=");
00309 residual_ += su.internalField();
00310 }
00311
00312
00313 template<class Type>
00314 void Foam::errorEstimate<Type>::operator-=
00315 (
00316 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
00317 )
00318 {
00319 operator-=(tsu());
00320 tsu.clear();
00321 }
00322
00323
00324 template<class Type>
00325 void Foam::errorEstimate<Type>::operator+=
00326 (
00327 const dimensioned<Type>& su
00328 )
00329 {
00330 residual_ -= su;
00331 }
00332
00333
00334 template<class Type>
00335 void Foam::errorEstimate<Type>::operator-=
00336 (
00337 const dimensioned<Type>& su
00338 )
00339 {
00340 residual_ += su;
00341 }
00342
00343
00344 template<class Type>
00345 void Foam::errorEstimate<Type>::operator*=
00346 (
00347 const volScalarField& vsf
00348 )
00349 {
00350 dimensions_ *= vsf.dimensions();
00351 residual_ *= vsf.internalField();
00352 normFactor_ *= vsf.internalField();
00353 }
00354
00355
00356 template<class Type>
00357 void Foam::errorEstimate<Type>::operator*=
00358 (
00359 const tmp<volScalarField>& tvsf
00360 )
00361 {
00362 operator*=(tvsf());
00363 tvsf.clear();
00364 }
00365
00366
00367 template<class Type>
00368 void Foam::errorEstimate<Type>::operator*=
00369 (
00370 const dimensioned<scalar>& ds
00371 )
00372 {
00373 dimensions_ *= ds.dimensions();
00374 residual_ *= ds.value();
00375 normFactor_ *= ds.value();
00376 }
00377
00378
00379
00380
00381 template<class Type>
00382 void Foam::checkMethod
00383 (
00384 const errorEstimate<Type>& ee1,
00385 const errorEstimate<Type>& ee2,
00386 const char* op
00387 )
00388 {
00389 if (&ee1.psi() != &ee2.psi())
00390 {
00391 FatalErrorIn
00392 (
00393 "checkMethod(const errorEstimate<Type>&, "
00394 "const errorEstimate<Type>&)"
00395 ) << "incompatible fields for operation "
00396 << endl << " "
00397 << "[" << ee1.psi().name() << "] "
00398 << op
00399 << " [" << ee2.psi().name() << "]"
00400 << abort(FatalError);
00401 }
00402
00403 if (dimensionSet::debug && ee1.dimensions() != ee2.dimensions())
00404 {
00405 FatalErrorIn
00406 (
00407 "checkMethod(const errorEstimate<Type>&, "
00408 "const errorEstimate<Type>&)"
00409 ) << "incompatible dimensions for operation "
00410 << endl << " "
00411 << "[" << ee1.psi().name() << ee1.dimensions()/dimVolume << " ] "
00412 << op
00413 << " [" << ee2.psi().name() << ee2.dimensions()/dimVolume << " ]"
00414 << abort(FatalError);
00415 }
00416 }
00417
00418
00419 template<class Type>
00420 void Foam::checkMethod
00421 (
00422 const errorEstimate<Type>& ee,
00423 const GeometricField<Type, fvPatchField, volMesh>& vf,
00424 const char* op
00425 )
00426 {
00427 if (dimensionSet::debug && ee.dimensions()/dimVolume != vf.dimensions())
00428 {
00429 FatalErrorIn
00430 (
00431 "checkMethod(const errorEstimate<Type>&, "
00432 "const GeometricField<Type, fvPatchField, volMesh>&)"
00433 ) << "incompatible dimensions for operation "
00434 << endl << " "
00435 << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
00436 << op
00437 << " [" << vf.name() << vf.dimensions() << " ]"
00438 << abort(FatalError);
00439 }
00440 }
00441
00442
00443 template<class Type>
00444 void Foam::checkMethod
00445 (
00446 const errorEstimate<Type>& ee,
00447 const dimensioned<Type>& dt,
00448 const char* op
00449 )
00450 {
00451 if (dimensionSet::debug && ee.dimensions()/dimVolume != dt.dimensions())
00452 {
00453 FatalErrorIn
00454 (
00455 "checkMethod(const errorEstimate<Type>&, const dimensioned<Type>&)"
00456 ) << "incompatible dimensions for operation "
00457 << endl << " "
00458 << "[" << ee.psi().name() << ee.dimensions()/dimVolume << " ] "
00459 << op
00460 << " [" << dt.name() << dt.dimensions() << " ]"
00461 << abort(FatalError);
00462 }
00463 }
00464
00465
00466
00467
00468 namespace Foam
00469 {
00470
00471 template<class Type>
00472 tmp<errorEstimate<Type> > operator+
00473 (
00474 const errorEstimate<Type>& A,
00475 const errorEstimate<Type>& B
00476 )
00477 {
00478 checkMethod(A, B, "+");
00479 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00480 tC() += B;
00481 return tC;
00482 }
00483
00484
00485 template<class Type>
00486 tmp<errorEstimate<Type> > operator+
00487 (
00488 const tmp<errorEstimate<Type> >& tA,
00489 const errorEstimate<Type>& B
00490 )
00491 {
00492 checkMethod(tA(), B, "+");
00493 tmp<errorEstimate<Type> > tC(tA.ptr());
00494 tC() += B;
00495 return tC;
00496 }
00497
00498
00499 template<class Type>
00500 tmp<errorEstimate<Type> > operator+
00501 (
00502 const errorEstimate<Type>& A,
00503 const tmp<errorEstimate<Type> >& tB
00504 )
00505 {
00506 checkMethod(A, tB(), "+");
00507 tmp<errorEstimate<Type> > tC(tB.ptr());
00508 tC() += A;
00509 return tC;
00510 }
00511
00512
00513 template<class Type>
00514 tmp<errorEstimate<Type> > operator+
00515 (
00516 const tmp<errorEstimate<Type> >& tA,
00517 const tmp<errorEstimate<Type> >& tB
00518 )
00519 {
00520 checkMethod(tA(), tB(), "+");
00521 tmp<errorEstimate<Type> > tC(tA.ptr());
00522 tC() += tB();
00523 tB.clear();
00524 return tC;
00525 }
00526
00527
00528 template<class Type>
00529 tmp<errorEstimate<Type> > operator-
00530 (
00531 const errorEstimate<Type>& A
00532 )
00533 {
00534 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00535 tC().negate();
00536 return tC;
00537 }
00538
00539
00540 template<class Type>
00541 tmp<errorEstimate<Type> > operator-
00542 (
00543 const tmp<errorEstimate<Type> >& tA
00544 )
00545 {
00546 tmp<errorEstimate<Type> > tC(tA.ptr());
00547 tC().negate();
00548 return tC;
00549 }
00550
00551
00552 template<class Type>
00553 tmp<errorEstimate<Type> > operator-
00554 (
00555 const errorEstimate<Type>& A,
00556 const errorEstimate<Type>& B
00557 )
00558 {
00559 checkMethod(A, B, "-");
00560 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00561 tC() -= B;
00562 return tC;
00563 }
00564
00565
00566 template<class Type>
00567 tmp<errorEstimate<Type> > operator-
00568 (
00569 const tmp<errorEstimate<Type> >& tA,
00570 const errorEstimate<Type>& B
00571 )
00572 {
00573 checkMethod(tA(), B, "-");
00574 tmp<errorEstimate<Type> > tC(tA.ptr());
00575 tC() -= B;
00576 return tC;
00577 }
00578
00579
00580 template<class Type>
00581 tmp<errorEstimate<Type> > operator-
00582 (
00583 const errorEstimate<Type>& A,
00584 const tmp<errorEstimate<Type> >& tB
00585 )
00586 {
00587 checkMethod(A, tB(), "-");
00588 tmp<errorEstimate<Type> > tC(tB.ptr());
00589 tC() -= A;
00590 tC().negate();
00591 return tC;
00592 }
00593
00594
00595 template<class Type>
00596 tmp<errorEstimate<Type> > operator-
00597 (
00598 const tmp<errorEstimate<Type> >& tA,
00599 const tmp<errorEstimate<Type> >& tB
00600 )
00601 {
00602 checkMethod(tA(), tB(), "-");
00603 tmp<errorEstimate<Type> > tC(tA.ptr());
00604 tC() -= tB();
00605 tB.clear();
00606 return tC;
00607 }
00608
00609
00610 template<class Type>
00611 tmp<errorEstimate<Type> > operator==
00612 (
00613 const errorEstimate<Type>& A,
00614 const errorEstimate<Type>& B
00615 )
00616 {
00617 checkMethod(A, B, "==");
00618 return (A - B);
00619 }
00620
00621
00622 template<class Type>
00623 tmp<errorEstimate<Type> > operator==
00624 (
00625 const tmp<errorEstimate<Type> >& tA,
00626 const errorEstimate<Type>& B
00627 )
00628 {
00629 checkMethod(tA(), B, "==");
00630 return (tA - B);
00631 }
00632
00633
00634 template<class Type>
00635 tmp<errorEstimate<Type> > operator==
00636 (
00637 const errorEstimate<Type>& A,
00638 const tmp<errorEstimate<Type> >& tB
00639 )
00640 {
00641 checkMethod(A, tB(), "==");
00642 return (A - tB);
00643 }
00644
00645
00646 template<class Type>
00647 tmp<errorEstimate<Type> > operator==
00648 (
00649 const tmp<errorEstimate<Type> >& tA,
00650 const tmp<errorEstimate<Type> >& tB
00651 )
00652 {
00653 checkMethod(tA(), tB(), "==");
00654 return (tA - tB);
00655 }
00656
00657
00658 template<class Type>
00659 tmp<errorEstimate<Type> > operator+
00660 (
00661 const errorEstimate<Type>& A,
00662 const GeometricField<Type, fvPatchField, volMesh>& su
00663 )
00664 {
00665 checkMethod(A, su, "+");
00666 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00667 tC().res() -= su.internalField();
00668 return tC;
00669 }
00670
00671 template<class Type>
00672 tmp<errorEstimate<Type> > operator+
00673 (
00674 const tmp<errorEstimate<Type> >& tA,
00675 const GeometricField<Type, fvPatchField, volMesh>& su
00676 )
00677 {
00678 checkMethod(tA(), su, "+");
00679 tmp<errorEstimate<Type> > tC(tA.ptr());
00680 tC().res() -= su.internalField();
00681 return tC;
00682 }
00683
00684 template<class Type>
00685 tmp<errorEstimate<Type> > operator+
00686 (
00687 const errorEstimate<Type>& A,
00688 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
00689 )
00690 {
00691 checkMethod(A, tsu(), "+");
00692 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00693 tC().res() -= tsu().internalField();
00694 tsu.clear();
00695 return tC;
00696 }
00697
00698
00699 template<class Type>
00700 tmp<errorEstimate<Type> > operator+
00701 (
00702 const tmp<errorEstimate<Type> >& tA,
00703 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
00704 )
00705 {
00706 checkMethod(tA(), tsu(), "+");
00707 tmp<errorEstimate<Type> > tC(tA.ptr());
00708 tC().res() -= tsu().internalField();
00709 tsu.clear();
00710 return tC;
00711 }
00712
00713 template<class Type>
00714 tmp<errorEstimate<Type> > operator+
00715 (
00716 const GeometricField<Type, fvPatchField, volMesh>& su,
00717 const errorEstimate<Type>& A
00718 )
00719 {
00720 checkMethod(A, su, "+");
00721 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00722 tC().res() -= su.internalField();
00723 return tC;
00724 }
00725
00726 template<class Type>
00727 tmp<errorEstimate<Type> > operator+
00728 (
00729 const GeometricField<Type, fvPatchField, volMesh>& su,
00730 const tmp<errorEstimate<Type> >& tA
00731 )
00732 {
00733 checkMethod(tA(), su, "+");
00734 tmp<errorEstimate<Type> > tC(tA.ptr());
00735 tC().res() -= su.internalField();
00736 return tC;
00737 }
00738
00739 template<class Type>
00740 tmp<errorEstimate<Type> > operator+
00741 (
00742 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
00743 const errorEstimate<Type>& A
00744 )
00745 {
00746 checkMethod(A, tsu(), "+");
00747 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00748 tC().res() -= tsu().internalField();
00749 tsu.clear();
00750 return tC;
00751 }
00752
00753 template<class Type>
00754 tmp<errorEstimate<Type> > operator+
00755 (
00756 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
00757 const tmp<errorEstimate<Type> >& tA
00758 )
00759 {
00760 checkMethod(tA(), tsu(), "+");
00761 tmp<errorEstimate<Type> > tC(tA.ptr());
00762 tC().res() -= tsu().internalField();
00763 tsu.clear();
00764 return tC;
00765 }
00766
00767
00768 template<class Type>
00769 tmp<errorEstimate<Type> > operator-
00770 (
00771 const errorEstimate<Type>& A,
00772 const GeometricField<Type, fvPatchField, volMesh>& su
00773 )
00774 {
00775 checkMethod(A, su, "-");
00776 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00777 tC().res() += su.internalField();
00778 return tC;
00779 }
00780
00781 template<class Type>
00782 tmp<errorEstimate<Type> > operator-
00783 (
00784 const tmp<errorEstimate<Type> >& tA,
00785 const GeometricField<Type, fvPatchField, volMesh>& su
00786 )
00787 {
00788 checkMethod(tA(), su, "-");
00789 tmp<errorEstimate<Type> > tC(tA.ptr());
00790 tC().res() += su.internalField();
00791 return tC;
00792 }
00793
00794 template<class Type>
00795 tmp<errorEstimate<Type> > operator-
00796 (
00797 const errorEstimate<Type>& A,
00798 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
00799 )
00800 {
00801 checkMethod(A, tsu(), "-");
00802 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00803 tC().res() += tsu().internalField();
00804 tsu.clear();
00805 return tC;
00806 }
00807
00808 template<class Type>
00809 tmp<errorEstimate<Type> > operator-
00810 (
00811 const tmp<errorEstimate<Type> >& tA,
00812 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
00813 )
00814 {
00815 checkMethod(tA(), tsu(), "-");
00816 tmp<errorEstimate<Type> > tC(tA.ptr());
00817 tC().res() += tsu().internalField();
00818 tsu.clear();
00819 return tC;
00820 }
00821
00822
00823 template<class Type>
00824 tmp<errorEstimate<Type> > operator-
00825 (
00826 const GeometricField<Type, fvPatchField, volMesh>& su,
00827 const errorEstimate<Type>& A
00828 )
00829 {
00830 checkMethod(A, su, "-");
00831 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00832 tC().negate();
00833 tC().res() -= su.internalField();
00834 return tC;
00835 }
00836
00837
00838 template<class Type>
00839 tmp<errorEstimate<Type> > operator-
00840 (
00841 const GeometricField<Type, fvPatchField, volMesh>& su,
00842 const tmp<errorEstimate<Type> >& tA
00843 )
00844 {
00845 checkMethod(tA(), su, "-");
00846 tmp<errorEstimate<Type> > tC(tA.ptr());
00847 tC().negate();
00848 tC().res() -= su.internalField();
00849 return tC;
00850 }
00851
00852 template<class Type>
00853 tmp<errorEstimate<Type> > operator-
00854 (
00855 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
00856 const errorEstimate<Type>& A
00857 )
00858 {
00859 checkMethod(A, tsu(), "-");
00860 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00861 tC().negate();
00862 tC().res() -= tsu().internalField();
00863 tsu.clear();
00864 return tC;
00865 }
00866
00867
00868 template<class Type>
00869 tmp<errorEstimate<Type> > operator-
00870 (
00871 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
00872 const tmp<errorEstimate<Type> >& tA
00873 )
00874 {
00875 checkMethod(tA(), tsu(), "-");
00876 tmp<errorEstimate<Type> > tC(tA.ptr());
00877 tC().negate();
00878 tC().res() -= tsu().internalField();
00879 tsu.clear();
00880 return tC;
00881 }
00882
00883
00884 template<class Type>
00885 tmp<errorEstimate<Type> > operator+
00886 (
00887 const errorEstimate<Type>& A,
00888 const dimensioned<Type>& su
00889 )
00890 {
00891 checkMethod(A, su, "+");
00892 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00893 tC().res() -= su.value();
00894 return tC;
00895 }
00896
00897
00898 template<class Type>
00899 tmp<errorEstimate<Type> > operator+
00900 (
00901 const tmp<errorEstimate<Type> >& tA,
00902 const dimensioned<Type>& su
00903 )
00904 {
00905 checkMethod(tA(), su, "+");
00906 tmp<errorEstimate<Type> > tC(tA.ptr());
00907 tC().res() -= su.value();
00908 return tC;
00909 }
00910
00911
00912 template<class Type>
00913 tmp<errorEstimate<Type> > operator+
00914 (
00915 const dimensioned<Type>& su,
00916 const errorEstimate<Type>& A
00917 )
00918 {
00919 checkMethod(A, su, "+");
00920 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00921 tC().res() -= su.value();
00922 return tC;
00923 }
00924
00925
00926 template<class Type>
00927 tmp<errorEstimate<Type> > operator+
00928 (
00929 const dimensioned<Type>& su,
00930 const tmp<errorEstimate<Type> >& tA
00931 )
00932 {
00933 checkMethod(tA(), su, "+");
00934 tmp<errorEstimate<Type> > tC(tA.ptr());
00935 tC().res() -= su.value();
00936 return tC;
00937 }
00938
00939
00940 template<class Type>
00941 tmp<errorEstimate<Type> > operator-
00942 (
00943 const errorEstimate<Type>& A,
00944 const dimensioned<Type>& su
00945 )
00946 {
00947 checkMethod(A, su, "-");
00948 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00949 tC().res() += su.value();
00950 return tC;
00951 }
00952
00953
00954 template<class Type>
00955 tmp<errorEstimate<Type> > operator-
00956 (
00957 const tmp<errorEstimate<Type> >& tA,
00958 const dimensioned<Type>& su
00959 )
00960 {
00961 checkMethod(tA(), su, "-");
00962 tmp<errorEstimate<Type> > tC(tA.ptr());
00963 tC().res() += su.value();
00964 return tC;
00965 }
00966
00967
00968 template<class Type>
00969 tmp<errorEstimate<Type> > operator-
00970 (
00971 const dimensioned<Type>& su,
00972 const errorEstimate<Type>& A
00973 )
00974 {
00975 checkMethod(A, su, "-");
00976 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
00977 tC().negate();
00978 tC().res() -= su.value();
00979 return tC;
00980 }
00981
00982
00983 template<class Type>
00984 tmp<errorEstimate<Type> > operator-
00985 (
00986 const dimensioned<Type>& su,
00987 const tmp<errorEstimate<Type> >& tA
00988 )
00989 {
00990 checkMethod(tA(), su, "-");
00991 tmp<errorEstimate<Type> > tC(tA.ptr());
00992 tC().negate();
00993 tC().res() -= su.value();
00994 return tC;
00995 }
00996
00997
00998 template<class Type>
00999 tmp<errorEstimate<Type> > operator==
01000 (
01001 const errorEstimate<Type>& A,
01002 const GeometricField<Type, fvPatchField, volMesh>& su
01003 )
01004 {
01005 checkMethod(A, su, "==");
01006 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
01007 tC().res() += su.internalField();
01008 return tC;
01009 }
01010
01011 template<class Type>
01012 tmp<errorEstimate<Type> > operator==
01013 (
01014 const tmp<errorEstimate<Type> >& tA,
01015 const GeometricField<Type, fvPatchField, volMesh>& su
01016 )
01017 {
01018 checkMethod(tA(), su, "==");
01019 tmp<errorEstimate<Type> > tC(tA.ptr());
01020 tC().res() += su.internalField();
01021 return tC;
01022 }
01023
01024 template<class Type>
01025 tmp<errorEstimate<Type> > operator==
01026 (
01027 const errorEstimate<Type>& A,
01028 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01029 )
01030 {
01031 checkMethod(A, tsu(), "==");
01032 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
01033 tC().res() += tsu().internalField();
01034 tsu.clear();
01035 return tC;
01036 }
01037
01038 template<class Type>
01039 tmp<errorEstimate<Type> > operator==
01040 (
01041 const tmp<errorEstimate<Type> >& tA,
01042 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01043 )
01044 {
01045 checkMethod(tA(), tsu(), "==");
01046 tmp<errorEstimate<Type> > tC(tA.ptr());
01047 tC().res() += tsu().internalField();
01048 tsu.clear();
01049 return tC;
01050 }
01051
01052
01053 template<class Type>
01054 tmp<errorEstimate<Type> > operator==
01055 (
01056 const errorEstimate<Type>& A,
01057 const dimensioned<Type>& su
01058 )
01059 {
01060 checkMethod(A, su, "==");
01061 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
01062 tC().res() += su.value();
01063 return tC;
01064 }
01065
01066
01067 template<class Type>
01068 tmp<errorEstimate<Type> > operator==
01069 (
01070 const tmp<errorEstimate<Type> >& tA,
01071 const dimensioned<Type>& su
01072 )
01073 {
01074 checkMethod(tA(), su, "==");
01075 tmp<errorEstimate<Type> > tC(tA.ptr());
01076 tC().res() += su.value();
01077 return tC;
01078 }
01079
01080
01081 template<class Type>
01082 tmp<errorEstimate<Type> > operator*
01083 (
01084 const volScalarField& vsf,
01085 const errorEstimate<Type>& A
01086 )
01087 {
01088 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
01089 tC() *= vsf;
01090 return tC;
01091 }
01092
01093 template<class Type>
01094 tmp<errorEstimate<Type> > operator*
01095 (
01096 const tmp<volScalarField>& tvsf,
01097 const errorEstimate<Type>& A
01098 )
01099 {
01100 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
01101 tC() *= tvsf;
01102 return tC;
01103 }
01104
01105 template<class Type>
01106 tmp<errorEstimate<Type> > operator*
01107 (
01108 const volScalarField& vsf,
01109 const tmp<errorEstimate<Type> >& tA
01110 )
01111 {
01112 tmp<errorEstimate<Type> > tC(tA.ptr());
01113 tC() *= vsf;
01114 return tC;
01115 }
01116
01117 template<class Type>
01118 tmp<errorEstimate<Type> > operator*
01119 (
01120 const tmp<volScalarField>& tvsf,
01121 const tmp<errorEstimate<Type> >& tA
01122 )
01123 {
01124 tmp<errorEstimate<Type> > tC(tA.ptr());
01125 tC() *= tvsf;
01126 return tC;
01127 }
01128
01129
01130 template<class Type>
01131 tmp<errorEstimate<Type> > operator*
01132 (
01133 const dimensioned<scalar>& ds,
01134 const errorEstimate<Type>& A
01135 )
01136 {
01137 tmp<errorEstimate<Type> > tC(new errorEstimate<Type>(A));
01138 tC() *= ds;
01139 return tC;
01140 }
01141
01142
01143 template<class Type>
01144 tmp<errorEstimate<Type> > operator*
01145 (
01146 const dimensioned<scalar>& ds,
01147 const tmp<errorEstimate<Type> >& tA
01148 )
01149 {
01150 tmp<errorEstimate<Type> > tC(tA.ptr());
01151 tC() *= ds;
01152 return tC;
01153 }
01154
01155
01156
01157
01158 }
01159
01160