Loading...
Searching...
No Matches
linearValveFvMesh.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "linearValveFvMesh.H"
30#include "Time.H"
31#include "slidingInterface.H"
32#include "mapPolyMesh.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48void Foam::linearValveFvMesh::addZonesAndModifiers()
49{
50 // Add zones and modifiers for motion action
51
52 if
53 (
55 || faceZones().size()
56 || cellZones().size()
58 )
59 {
61 << "Zones and modifiers already present. Skipping."
62 << endl;
63
64 return;
65 }
66
67 Info<< "Time = " << time().timeName() << endl
68 << "Adding zones and modifiers to the mesh" << endl;
69
70 // Add zones
71 List<pointZone*> pz(1);
72
73 // An empty zone for cut points
74 pz[0] = new pointZone("cutPointZone", 0, pointZones());
75
76
77 // Do face zones for slider
78
79 List<faceZone*> fz(3);
80
81 // Inner slider
82 const word innerSliderName
83 (
84 motionDict_.subDict("slider").get<word>("inside")
85 );
86 const polyPatch& innerSlider = boundaryMesh()[innerSliderName];
87
88 fz[0] = new faceZone
89 (
90 "insideSliderZone",
91 identity(innerSlider.range()),
92 false, // none are flipped
93 0,
94 faceZones()
95 );
96
97 // Outer slider
98 const word outerSliderName
99 (
100 motionDict_.subDict("slider").get<word>("outside")
101 );
102 const polyPatch& outerSlider = boundaryMesh()[outerSliderName];
103
104 fz[1] = new faceZone
105 (
106 "outsideSliderZone",
107 identity(outerSlider.range()),
108 false, // none are flipped
109 1,
110 faceZones()
111 );
112
113 // An empty zone for cut faces
114 fz[2] = new faceZone("cutFaceZone", 2, faceZones());
115
116 List<cellZone*> cz(0);
117
118 Info<< "Adding point, face and cell zones" << endl;
119 addZones(pz, fz, cz);
120
121 // Add a topology modifier
122 Info<< "Adding topology modifiers" << endl;
123 topoChanger_.setSize(1);
124 topoChanger_.set
125 (
126 0,
127 new slidingInterface
128 (
129 "mixerSlider",
130 0,
132 outerSliderName + "Zone",
133 innerSliderName + "Zone",
134 "cutPointZone",
135 "cutFaceZone",
136 outerSliderName,
137 innerSliderName,
139 true // Attach-detach action
140 )
141 );
143
144 // Write mesh
145 write();
146}
147
148
149void Foam::linearValveFvMesh::makeSlidersDead()
150{
151 const polyTopoChanger& topoChanges = topoChanger_;
152
153 // Enable layering
154 forAll(topoChanges, modI)
155 {
156 if (isA<slidingInterface>(topoChanges[modI]))
157 {
158 topoChanges[modI].disable();
159 }
160 else
161 {
163 << "Don't know what to do with mesh modifier "
164 << modI << " of type " << topoChanges[modI].type()
165 << abort(FatalError);
166 }
167 }
168}
169
170
171void Foam::linearValveFvMesh::makeSlidersLive()
172{
173 const polyTopoChanger& topoChanges = topoChanger_;
174
175 // Enable sliding interface
176 forAll(topoChanges, modI)
177 {
178 if (isA<slidingInterface>(topoChanges[modI]))
179 {
180 topoChanges[modI].enable();
181 }
182 else
183 {
185 << "Don't know what to do with mesh modifier "
186 << modI << " of type " << topoChanges[modI].type()
187 << abort(FatalError);
188 }
189 }
190}
191
192
193bool Foam::linearValveFvMesh::attached() const
194{
195 const polyTopoChanger& topoChanges = topoChanger_;
196
197 bool result = false;
198
199 forAll(topoChanges, modI)
200 {
201 if (isA<slidingInterface>(topoChanges[modI]))
202 {
203 result =
204 result
205 || refCast<const slidingInterface>(topoChanges[modI]).attached();
206 }
207 }
208
209 // Check thal all sliders are in sync (debug only)
210 forAll(topoChanges, modI)
211 {
212 if (isA<slidingInterface>(topoChanges[modI]))
213 {
214 if
215 (
216 result
217 != refCast<const slidingInterface>(topoChanges[modI]).attached()
218 )
219 {
221 << "Slider " << modI
222 << " named " << topoChanges[modI].name()
223 << " out of sync: Should be" << result
224 << abort(FatalError);
225 }
226 }
227 }
228
229 if (result)
230 {
231 Info<< "linearValveFvMesh: attached!" << endl;
232 }
233 else
234 {
235 Info<< "linearValveFvMesh: detached!" << endl;
236 }
237
238 return result;
239}
240
241
242// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
243
244// Construct from components
245Foam::linearValveFvMesh::linearValveFvMesh(const IOobject& io)
246:
248 motionDict_
249 (
250 IOdictionary::readContents
251 (
253 (
254 "dynamicMeshDict",
255 time().constant(),
256 *this,
257 IOobject::MUST_READ
258 )
259 ).optionalSubDict(typeName + "Coeffs")
260 ),
261 msPtr_(motionSolver::New(*this))
263 addZonesAndModifiers();
264}
265
266
267// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
270{}
271
272
273// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
274
276{
277 // Detaching the interface
278 if (attached())
279 {
280 Info<< "Decoupling sliding interfaces" << endl;
281 makeSlidersLive();
282
283 // Changing topology by hand
284 resetMorph();
285 setMorphTimeIndex(3*time().timeIndex());
286 updateMesh();
287
288 msPtr_->updateMesh();
289 }
290 else
291 {
292 Info<< "Sliding interfaces decoupled" << endl;
293 }
294
295 // Perform mesh motion
296 makeSlidersDead();
297
298 // Changing topology by hand
299 resetMorph();
300 setMorphTimeIndex(3*time().timeIndex() + 1);
301 updateMesh();
302
303 msPtr_->updateMesh();
304
305 if (topoChangeMap)
306 {
307 if (topoChangeMap().hasMotionPoints())
308 {
309 Info<< "Topology change; executing pre-motion" << endl;
310 movePoints(topoChangeMap().preMotionPoints());
311 }
312 }
313
314 // Solve for motion
315 msPtr_->solve();
316
317 movePoints(msPtr_->curPoints());
318
319 // Attach the interface
320 Info<< "Coupling sliding interfaces" << endl;
321 makeSlidersLive();
322 resetMorph();
323 setMorphTimeIndex(3*time().timeIndex() + 2);
324 updateMesh();
325
326 Info<< "Moving points post slider attach" << endl;
327
328 msPtr_->updateMesh();
329
330 Info<< "Sliding interfaces coupled: " << attached() << endl;
331}
332
333
334// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
Definition HashTable.H:358
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
static dictionary readContents(const IOobject &io)
Read and return contents, testing for "dictionary" type. The IOobject will not be registered.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition dictionary.C:560
static autoPtr< dynamicFvMesh > New(const IOobject &io)
Select, construct and return the dynamicFvMesh.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition fvMesh.C:960
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition fvMesh.C:884
A sliding linear valve.
virtual ~linearValveFvMesh()
Destructor.
virtual bool update()
Update the mesh for both mesh motion and topology change.
Virtual base class for mesh motion solver.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition polyMesh.C:994
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition polyMesh.H:663
List of mesh modifiers defining the mesh dynamics.
Abstract base class for a topology changing fvMesh.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
#define InfoInFunction
Report an information message using Foam::Info.
Different types of constants.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
runTime write()
label timeIndex
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299