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 "pressureGradientExplicitSource.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <OpenFOAM/IFstream.H>
00029
00030
00031
00032 void Foam::pressureGradientExplicitSource::writeGradP() const
00033 {
00034
00035 if (mesh_.time().outputTime())
00036 {
00037 IOdictionary propsDict
00038 (
00039 IOobject
00040 (
00041 sourceName_ + "Properties",
00042 mesh_.time().timeName(),
00043 "uniform",
00044 mesh_,
00045 IOobject::NO_READ,
00046 IOobject::NO_WRITE
00047 )
00048 );
00049 propsDict.add("gradient", gradP_);
00050 propsDict.write();
00051 }
00052 }
00053
00054
00055
00056
00057 Foam::pressureGradientExplicitSource::pressureGradientExplicitSource
00058 (
00059 const word& sourceName,
00060 volVectorField& U
00061 )
00062 :
00063 sourceName_(sourceName),
00064 mesh_(U.mesh()),
00065 U_(U),
00066 dict_
00067 (
00068 IOobject
00069 (
00070 sourceName + "Properties",
00071 mesh_.time().constant(),
00072 mesh_,
00073 IOobject::MUST_READ,
00074 IOobject::NO_WRITE
00075 )
00076 ),
00077 Ubar_(dict_.lookup("Ubar")),
00078 gradPini_(dict_.lookup("gradPini")),
00079 gradP_(gradPini_),
00080 flowDir_(Ubar_/mag(Ubar_)),
00081 cellSource_(dict_.lookup("cellSource")),
00082 cellSelector_
00083 (
00084 topoSetSource::New
00085 (
00086 cellSource_,
00087 mesh_,
00088 dict_.subDict(cellSource_ + "Coeffs")
00089 )
00090 ),
00091 selectedCellSet_
00092 (
00093 mesh_,
00094 sourceName_ + "CellSet",
00095 mesh_.nCells()/10 + 1
00096 )
00097 {
00098
00099 cellSelector_->applyToSet
00100 (
00101 topoSetSource::NEW,
00102 selectedCellSet_
00103 );
00104
00105
00106 Info<< " Selected "
00107 << returnReduce(selectedCellSet_.size(), sumOp<label>())
00108 << " cells" << endl;
00109
00110
00111 IFstream propsFile
00112 (
00113 mesh_.time().timeName()/"uniform"/(sourceName_ + "Properties")
00114 );
00115
00116 if (propsFile.good())
00117 {
00118 Info<< " Reading pressure gradient from file" << endl;
00119 dictionary propsDict(dictionary::null, propsFile);
00120 propsDict.lookup("gradient") >> gradP_;
00121 }
00122
00123 Info<< " Initial pressure gradient = " << gradP_ << nl << endl;
00124 }
00125
00126
00127
00128
00129 Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
00130 Foam::pressureGradientExplicitSource::Su() const
00131 {
00132 tmp<DimensionedField<vector, volMesh> > tSource
00133 (
00134 new DimensionedField<vector, volMesh>
00135 (
00136 IOobject
00137 (
00138 sourceName_,
00139 mesh_.time().timeName(),
00140 mesh_,
00141 IOobject::NO_READ,
00142 IOobject::NO_WRITE
00143 ),
00144 mesh_,
00145 dimensionedVector("zero", gradP_.dimensions(), vector::zero)
00146 )
00147 );
00148
00149 DimensionedField<vector, volMesh>& sourceField = tSource();
00150
00151 forAllConstIter(cellSet, selectedCellSet_, iter)
00152 {
00153 label cellI = iter.key();
00154
00155 sourceField[cellI] = flowDir_*gradP_.value();
00156 }
00157
00158 return tSource;
00159 }
00160
00161
00162 void Foam::pressureGradientExplicitSource::update()
00163 {
00164 const volScalarField& rUA =
00165 mesh_.lookupObject<volScalarField>("(1|A(" + U_.name() + "))");
00166
00167
00168 scalar volTot = 0.0;
00169 scalar magUbarAve = 0.0;
00170 scalar rUAave = 0.0;
00171 forAllConstIter(cellSet, selectedCellSet_, iter)
00172 {
00173 label cellI = iter.key();
00174
00175 scalar volCell = mesh_.V()[cellI];
00176 volTot += volCell;
00177
00178 magUbarAve += (flowDir_ & U_[cellI])*volCell;
00179 rUAave += rUA[cellI]*volCell;
00180 }
00181
00182
00183 reduce(volTot, sumOp<scalar>());
00184 reduce(magUbarAve, sumOp<scalar>());
00185 reduce(rUAave, sumOp<scalar>());
00186
00187
00188 magUbarAve /= volTot;
00189 rUAave /= volTot;
00190
00191
00192
00193 scalar gradPplus = (mag(Ubar_) - magUbarAve)/rUAave;
00194
00195
00196 forAllConstIter(cellSet, selectedCellSet_, iter)
00197 {
00198 label cellI = iter.key();
00199 U_[cellI] += flowDir_*rUA[cellI]*gradPplus;
00200 }
00201
00202
00203 gradP_.value() += gradPplus;
00204
00205 Info<< "Uncorrected Ubar = " << magUbarAve << tab
00206 << "Pressure gradient = " << gradP_.value() << endl;
00207
00208 writeGradP();
00209 }
00210
00211
00212