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 "mixerFvMesh.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <meshTools/regionSplit.H>
00029 #include <dynamicMesh/slidingInterface.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 #include <OpenFOAM/mapPolyMesh.H>
00032
00033
00034
00035 namespace Foam
00036 {
00037 defineTypeNameAndDebug(mixerFvMesh, 0);
00038
00039 addToRunTimeSelectionTable(topoChangerFvMesh, mixerFvMesh, IOobject);
00040 }
00041
00042
00043
00044
00045 void Foam::mixerFvMesh::addZonesAndModifiers()
00046 {
00047
00048
00049 if
00050 (
00051 pointZones().size()
00052 || faceZones().size()
00053 || cellZones().size()
00054 || topoChanger_.size()
00055 )
00056 {
00057 Info<< "void mixerFvMesh::addZonesAndModifiers() : "
00058 << "Zones and modifiers already present. Skipping."
00059 << endl;
00060
00061 return;
00062 }
00063
00064 Info<< "Time = " << time().timeName() << endl
00065 << "Adding zones and modifiers to the mesh" << endl;
00066
00067
00068 List<pointZone*> pz(1);
00069
00070
00071
00072 pz[0] = new pointZone
00073 (
00074 "cutPointZone",
00075 labelList(0),
00076 0,
00077 pointZones()
00078 );
00079
00080
00081
00082
00083 List<faceZone*> fz(3);
00084
00085
00086 const word innerSliderName(motionDict_.subDict("slider").lookup("inside"));
00087 const polyPatch& innerSlider =
00088 boundaryMesh()[boundaryMesh().findPatchID(innerSliderName)];
00089
00090 labelList isf(innerSlider.size());
00091
00092 forAll (isf, i)
00093 {
00094 isf[i] = innerSlider.start() + i;
00095 }
00096
00097 fz[0] = new faceZone
00098 (
00099 "insideSliderZone",
00100 isf,
00101 boolList(innerSlider.size(), false),
00102 0,
00103 faceZones()
00104 );
00105
00106
00107 const word outerSliderName(motionDict_.subDict("slider").lookup("outside"));
00108 const polyPatch& outerSlider =
00109 boundaryMesh()[boundaryMesh().findPatchID(outerSliderName)];
00110
00111 labelList osf(outerSlider.size());
00112
00113 forAll (osf, i)
00114 {
00115 osf[i] = outerSlider.start() + i;
00116 }
00117
00118 fz[1] = new faceZone
00119 (
00120 "outsideSliderZone",
00121 osf,
00122 boolList(outerSlider.size(), false),
00123 1,
00124 faceZones()
00125 );
00126
00127
00128 fz[2] = new faceZone
00129 (
00130 "cutFaceZone",
00131 labelList(0),
00132 boolList(0, false),
00133 2,
00134 faceZones()
00135 );
00136
00137 List<cellZone*> cz(1);
00138
00139
00140 regionSplit rs(*this);
00141
00142
00143 label originRegion = rs[findNearestCell(cs().origin())];
00144
00145 labelList movingCells(nCells());
00146 label nMovingCells = 0;
00147
00148 forAll(rs, cellI)
00149 {
00150 if (rs[cellI] == originRegion)
00151 {
00152 movingCells[nMovingCells] = cellI;
00153 nMovingCells++;
00154 }
00155 }
00156
00157 movingCells.setSize(nMovingCells);
00158 Info << "Number of cells in the moving region: " << nMovingCells << endl;
00159
00160 cz[0] = new cellZone
00161 (
00162 "movingCells",
00163 movingCells,
00164 0,
00165 cellZones()
00166 );
00167
00168 Info << "Adding point, face and cell zones" << endl;
00169 addZones(pz, fz, cz);
00170
00171
00172 Info << "Adding topology modifiers" << endl;
00173 topoChanger_.setSize(1);
00174 topoChanger_.set
00175 (
00176 0,
00177 new slidingInterface
00178 (
00179 "mixerSlider",
00180 0,
00181 topoChanger_,
00182 outerSliderName + "Zone",
00183 innerSliderName + "Zone",
00184 "cutPointZone",
00185 "cutFaceZone",
00186 outerSliderName,
00187 innerSliderName,
00188 slidingInterface::INTEGRAL
00189 )
00190 );
00191 topoChanger_.writeOpt() = IOobject::AUTO_WRITE;
00192
00193 write();
00194 }
00195
00196
00197 void Foam::mixerFvMesh::calcMovingMasks() const
00198 {
00199 if (debug)
00200 {
00201 Info<< "void mixerFvMesh::calcMovingMasks() const : "
00202 << "Calculating point and cell masks"
00203 << endl;
00204 }
00205
00206 if (movingPointsMaskPtr_)
00207 {
00208 FatalErrorIn("void mixerFvMesh::calcMovingMasks() const")
00209 << "point mask already calculated"
00210 << abort(FatalError);
00211 }
00212
00213
00214 movingPointsMaskPtr_ = new scalarField(points().size(), 0);
00215 scalarField& movingPointsMask = *movingPointsMaskPtr_;
00216
00217 const cellList& c = cells();
00218 const faceList& f = faces();
00219
00220 const labelList& cellAddr =
00221 cellZones()[cellZones().findZoneID("movingCells")];
00222
00223 forAll (cellAddr, cellI)
00224 {
00225 const cell& curCell = c[cellAddr[cellI]];
00226
00227 forAll (curCell, faceI)
00228 {
00229
00230 const face& curFace = f[curCell[faceI]];
00231
00232 forAll (curFace, pointI)
00233 {
00234 movingPointsMask[curFace[pointI]] = 1;
00235 }
00236 }
00237 }
00238
00239 const word innerSliderZoneName
00240 (
00241 word(motionDict_.subDict("slider").lookup("inside"))
00242 + "Zone"
00243 );
00244
00245 const labelList& innerSliderAddr =
00246 faceZones()[faceZones().findZoneID(innerSliderZoneName)];
00247
00248 forAll (innerSliderAddr, faceI)
00249 {
00250 const face& curFace = f[innerSliderAddr[faceI]];
00251
00252 forAll (curFace, pointI)
00253 {
00254 movingPointsMask[curFace[pointI]] = 1;
00255 }
00256 }
00257
00258 const word outerSliderZoneName
00259 (
00260 word(motionDict_.subDict("slider").lookup("outside"))
00261 + "Zone"
00262 );
00263
00264 const labelList& outerSliderAddr =
00265 faceZones()[faceZones().findZoneID(outerSliderZoneName)];
00266
00267 forAll (outerSliderAddr, faceI)
00268 {
00269 const face& curFace = f[outerSliderAddr[faceI]];
00270
00271 forAll (curFace, pointI)
00272 {
00273 movingPointsMask[curFace[pointI]] = 0;
00274 }
00275 }
00276 }
00277
00278
00279
00280
00281
00282 Foam::mixerFvMesh::mixerFvMesh
00283 (
00284 const IOobject& io
00285 )
00286 :
00287 topoChangerFvMesh(io),
00288 motionDict_
00289 (
00290 IOdictionary
00291 (
00292 IOobject
00293 (
00294 "dynamicMeshDict",
00295 time().constant(),
00296 *this,
00297 IOobject::MUST_READ,
00298 IOobject::NO_WRITE
00299 )
00300 ).subDict(typeName + "Coeffs")
00301 ),
00302 csPtr_
00303 (
00304 coordinateSystem::New
00305 (
00306 "coordinateSystem",
00307 motionDict_.subDict("coordinateSystem")
00308 )
00309 ),
00310 rpm_(readScalar(motionDict_.lookup("rpm"))),
00311 movingPointsMaskPtr_(NULL)
00312 {
00313 addZonesAndModifiers();
00314
00315 Info<< "Mixer mesh:" << nl
00316 << " origin: " << cs().origin() << nl
00317 << " axis: " << cs().axis() << nl
00318 << " rpm: " << rpm_ << endl;
00319 }
00320
00321
00322
00323
00324 Foam::mixerFvMesh::~mixerFvMesh()
00325 {
00326 deleteDemandDrivenData(movingPointsMaskPtr_);
00327 }
00328
00329
00330
00331
00332 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
00333 {
00334 if (!movingPointsMaskPtr_)
00335 {
00336 calcMovingMasks();
00337 }
00338
00339 return *movingPointsMaskPtr_;
00340 }
00341
00342
00343 bool Foam::mixerFvMesh::update()
00344 {
00345
00346 movePoints
00347 (
00348 csPtr_->globalPosition
00349 (
00350 csPtr_->localPosition(points())
00351 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
00352 *movingPointsMask()
00353 )
00354 );
00355
00356
00357 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
00358
00359 if (topoChangeMap.valid())
00360 {
00361 if (debug)
00362 {
00363 Info << "Mesh topology is changing" << endl;
00364 }
00365
00366 deleteDemandDrivenData(movingPointsMaskPtr_);
00367 }
00368
00369 movePoints
00370 (
00371 csPtr_->globalPosition
00372 (
00373 csPtr_->localPosition(oldPoints())
00374 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
00375 *movingPointsMask()
00376 )
00377 );
00378
00379 return true;
00380 }
00381
00382
00383