Go to the documentation of this file.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 "adjustPhi.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <finiteVolume/surfaceFields.H>
00029 #include <finiteVolume/processorFvsPatchFields.H>
00030 #include <finiteVolume/inletOutletFvPatchFields.H>
00031
00032
00033
00034 bool Foam::adjustPhi
00035 (
00036 surfaceScalarField& phi,
00037 const volVectorField& U,
00038 volScalarField& p
00039 )
00040 {
00041 if (p.needReference())
00042 {
00043
00044
00045
00046 scalar massIn = 0.0;
00047 scalar fixedMassOut = 0.0;
00048 scalar adjustableMassOut = 0.0;
00049
00050 forAll (phi.boundaryField(), patchi)
00051 {
00052 const fvPatchVectorField& Up = U.boundaryField()[patchi];
00053 const fvsPatchScalarField& phip = phi.boundaryField()[patchi];
00054
00055 if (!isA<processorFvsPatchScalarField>(phip))
00056 {
00057 if
00058 (
00059 Up.fixesValue()
00060 && !isA<inletOutletFvPatchVectorField>(Up)
00061 )
00062 {
00063 forAll(phip, i)
00064 {
00065 if (phip[i] < 0.0)
00066 {
00067 massIn -= phip[i];
00068 }
00069 else
00070 {
00071 fixedMassOut += phip[i];
00072 }
00073 }
00074 }
00075 else
00076 {
00077 forAll(phip, i)
00078 {
00079 if (phip[i] < 0.0)
00080 {
00081 massIn -= phip[i];
00082 }
00083 else
00084 {
00085 adjustableMassOut += phip[i];
00086 }
00087 }
00088 }
00089 }
00090 }
00091
00092
00093 scalar totalFlux = VSMALL + sum(mag(phi)).value();
00094
00095 reduce(massIn, sumOp<scalar>());
00096 reduce(fixedMassOut, sumOp<scalar>());
00097 reduce(adjustableMassOut, sumOp<scalar>());
00098
00099 scalar massCorr = 1.0;
00100 scalar magAdjustableMassOut = mag(adjustableMassOut);
00101
00102 if
00103 (
00104 magAdjustableMassOut > VSMALL
00105 && magAdjustableMassOut/totalFlux > SMALL
00106 )
00107 {
00108 massCorr = (massIn - fixedMassOut)/adjustableMassOut;
00109 }
00110 else if(mag(fixedMassOut - massIn)/totalFlux > 1e-10)
00111 {
00112 FatalErrorIn
00113 (
00114 "adjustPhi(surfaceScalarField& phi, const volVectorField& U,"
00115 "const volScalarField& p"
00116 ) << "Continuity error cannot be removed by adjusting the"
00117 " outflow.\nPlease check the velocity boundary conditions"
00118 " and/or run potentialFoam to initialise the outflow." << nl
00119 << "Total flux : " << totalFlux << nl
00120 << "Specified mass inflow : " << massIn << nl
00121 << "Specified mass outflow : " << fixedMassOut << nl
00122 << "Adjustable mass outflow : " << adjustableMassOut << nl
00123 << exit(FatalError);
00124 }
00125
00126 forAll (phi.boundaryField(), patchi)
00127 {
00128 const fvPatchVectorField& Up = U.boundaryField()[patchi];
00129 fvsPatchScalarField& phip = phi.boundaryField()[patchi];
00130
00131 if (!isA<processorFvsPatchScalarField>(phip))
00132 {
00133 if
00134 (
00135 !Up.fixesValue()
00136 || isA<inletOutletFvPatchVectorField>(Up)
00137 )
00138 {
00139 forAll(phip, i)
00140 {
00141 if (phip[i] > 0.0)
00142 {
00143 phip[i] *= massCorr;
00144 }
00145 }
00146 }
00147 }
00148 }
00149
00150 return mag(massIn)/totalFlux < SMALL
00151 && mag(fixedMassOut)/totalFlux < SMALL
00152 && mag(adjustableMassOut)/totalFlux < SMALL;
00153 }
00154 else
00155 {
00156 return false;
00157 }
00158 }
00159
00160
00161