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 \*---------------------------------------------------------------------------*/
29 #include "polyTopoChanger.H"
30 #include "layerAdditionRemoval.H"
32 #include "transformField.H"
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 void Foam::topoBody::addZones
38 DynamicList<pointZone*>& pz,
39 DynamicList<faceZone*>& fz,
40 DynamicList<cellZone*>& cz
43 // Check presence of moving cell zone
44 const label mcZoneID = mesh_.cellZones().findZoneID(movingCellsName_);
50 "void topoBody::addZones\n"
52 " DynamicList<pointZone*>& pz,\n"
53 " DynamicList<faceZone*>& fz,\n"
54 " DynamicList<cellZone*>& cz\n"
56 ) << "Cannot find moving cell zone " << movingCellsName_
57 << " for body " << name_
61 forAll (layerFacesNames_, lfI)
64 mesh_.faceZones().findZoneID(layerFacesNames_[lfI]);
70 "void topoBody::addZones\n"
72 " DynamicList<pointZone*>& pz,\n"
73 " DynamicList<faceZone*>& fz,\n"
74 " DynamicList<cellZone*>& cz\n"
76 ) << "Cannot find layering zone " << layerFacesNames_[lfI]
77 << " for body " << name_
84 void Foam::topoBody::addModifiers
90 // Add a topology modifier
91 forAll (layerFacesNames_, lfI)
93 Info<< "Creating layering topology modifier " << layerFacesNames_[lfI]
94 << " on object " << name_ << endl;
99 new layerAdditionRemoval
101 layerFacesNames_[lfI] + "Layer" + name_,
104 layerFacesNames_[lfI],
115 void Foam::topoBody::calcMovingMask() const
117 if (movingPointsMaskPtr_)
119 FatalErrorIn("void topoBody::calcMovingMask() const")
120 << "point mask already calculated"
121 << abort(FatalError);
124 // Set the point mask
125 movingPointsMaskPtr_ = new scalarField(mesh_.allPoints().size(), 0);
126 scalarField& movingPointsMask = *movingPointsMaskPtr_;
128 const cellList& c = mesh_.cells();
129 const faceList& f = mesh_.allFaces();
131 const labelList& cellAddr = mesh_.cellZones()
132 [mesh_.cellZones().findZoneID(movingCellsName_)];
134 forAll (cellAddr, cellI)
136 const cell& curCell = c[cellAddr[cellI]];
138 forAll (curCell, faceI)
140 // Mark all the points as moving
141 const face& curFace = f[curCell[faceI]];
143 forAll (curFace, pointI)
145 movingPointsMask[curFace[pointI]] = 1;
152 // Return moving points mask. Moving points marked with 1
153 const Foam::scalarField& Foam::topoBody::movingPointsMask() const
155 if (!movingPointsMaskPtr_)
160 return *movingPointsMaskPtr_;
164 void Foam::topoBody::clearPointMask()
166 deleteDemandDrivenData(movingPointsMaskPtr_);
170 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
172 Foam::topoBody::topoBody
175 const polyMesh& mesh,
176 const dictionary& dict
181 movingCellsName_(dict.lookup("movingCells")),
182 layerFacesNames_(dict.lookup("layerFaces")),
183 minThickness_(readScalar(dict.lookup("minThickness"))),
184 maxThickness_(readScalar(dict.lookup("maxThickness"))),
185 SBMFPtr_(solidBodyMotionFunction::New(dict, mesh.time())),
188 dict.lookupOrDefault<bool>("invertMotionMask", false)
190 movingPointsMaskPtr_(NULL)
192 Info<< "Moving body " << name << ":" << nl
193 << " moving cells: " << movingCellsName_ << nl
194 << " layer faces : " << layerFacesNames_ << nl
195 << " invert mask : " << invertMotionMask_ << nl
200 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
202 Foam::topoBody::~topoBody()
208 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210 Foam::tmp<Foam::vectorField> Foam::topoBody::pointMotion() const
212 // Rotational speed needs to be converted from rpm
213 scalarField mpm = movingPointsMask();
215 if (invertMotionMask_)
217 Info << "Inverting motion mask" << endl;
221 return mpm*transform(SBMFPtr_().velocity(), mesh_.allPoints())*
222 mesh_.time().deltaT().value();
226 void Foam::topoBody::updateTopology()
232 // ************************************************************************* //