linearValveLayersFvMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "Time.H"
28 #include "slidingInterface.H"
29 #include "layerAdditionRemoval.H"
30 #include "pointField.H"
31 #include "mapPolyMesh.H"
32 #include "polyTopoChange.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(linearValveLayersFvMesh, 0);
41  (
42  topoChangerFvMesh,
43  linearValveLayersFvMesh,
44  IOobject
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::linearValveLayersFvMesh::addZonesAndModifiers()
52 {
53  // Add zones and modifiers for motion action
54 
55  if
56  (
57  pointZones().size()
58  || faceZones().size()
59  || cellZones().size()
60  || topoChanger_.size()
61  )
62  {
64  << "Zones and modifiers already present. Skipping."
65  << endl;
66 
67  return;
68  }
69 
70  Info<< "Time = " << time().timeName() << endl
71  << "Adding zones and modifiers to the mesh" << endl;
72 
73  // Add zones
74  List<pointZone*> pz(1);
75  List<faceZone*> fz(4);
76  List<cellZone*> cz(0);
77 
78 
79  // Add an empty zone for cut points
80 
81  pz[0] = new pointZone
82  (
83  "cutPointZone",
84  labelList(0),
85  0,
86  pointZones()
87  );
88 
89 
90  // Do face zones for slider
91 
92  // Inner slider
93  const word innerSliderName(motionDict_.subDict("slider").lookup("inside"));
94  const polyPatch& innerSlider = boundaryMesh()[innerSliderName];
95 
96  labelList isf(innerSlider.size());
97 
98  forAll(isf, i)
99  {
100  isf[i] = innerSlider.start() + i;
101  }
102 
103  fz[0] = new faceZone
104  (
105  "insideSliderZone",
106  isf,
107  boolList(innerSlider.size(), false),
108  0,
109  faceZones()
110  );
111 
112  // Outer slider
113  const word outerSliderName(motionDict_.subDict("slider").lookup("outside"));
114  const polyPatch& outerSlider = boundaryMesh()[outerSliderName];
115 
116  labelList osf(outerSlider.size());
117 
118  forAll(osf, i)
119  {
120  osf[i] = outerSlider.start() + i;
121  }
122 
123  fz[1] = new faceZone
124  (
125  "outsideSliderZone",
126  osf,
127  boolList(outerSlider.size(), false),
128  1,
129  faceZones()
130  );
131 
132  // Add empty zone for cut faces
133  fz[2] = new faceZone
134  (
135  "cutFaceZone",
136  labelList(0),
137  boolList(0, false),
138  2,
139  faceZones()
140  );
141 
142  // Add face zone for layer addition
143  const word layerPatchName
144  (
145  motionDict_.subDict("layer").lookup("patch")
146  );
147 
148  const polyPatch& layerPatch = boundaryMesh()[layerPatchName];
149 
150  labelList lpf(layerPatch.size());
151 
152  forAll(lpf, i)
153  {
154  lpf[i] = layerPatch.start() + i;
155  }
156 
157  fz[3] = new faceZone
158  (
159  "valveLayerZone",
160  lpf,
161  boolList(layerPatch.size(), true),
162  0,
163  faceZones()
164  );
165 
166 
167  Info<< "Adding point and face zones" << endl;
168  addZones(pz, fz, cz);
169 
170  // Add a topology modifier
171 
172  List<polyMeshModifier*> tm(2);
173 
174  tm[0] = new slidingInterface
175  (
176  "valveSlider",
177  0,
178  topoChanger_,
179  outerSliderName + "Zone",
180  innerSliderName + "Zone",
181  "cutPointZone",
182  "cutFaceZone",
183  outerSliderName,
184  innerSliderName,
186  true // Attach-detach action
187  );
188 
189  tm[1] =
190  new layerAdditionRemoval
191  (
192  "valveLayer",
193  1,
194  topoChanger_,
195  "valveLayerZone",
196  motionDict_.subDict("layer").lookup<scalar>("minThickness"),
197  motionDict_.subDict("layer").lookup<scalar>("maxThickness")
198  );
199 
200 
201  Info<< "Adding topology modifiers" << endl;
202  addTopologyModifiers(tm);
203 
204  // Write mesh
205  write();
206 }
207 
208 
209 void Foam::linearValveLayersFvMesh::makeLayersLive()
210 {
211  const polyTopoChanger& topoChanges = topoChanger_;
212 
213  // Enable layering
214  forAll(topoChanges, modI)
215  {
216  if (isA<layerAdditionRemoval>(topoChanges[modI]))
217  {
218  topoChanges[modI].enable();
219  }
220  else if (isA<slidingInterface>(topoChanges[modI]))
221  {
222  topoChanges[modI].disable();
223  }
224  else
225  {
227  << "Don't know what to do with mesh modifier "
228  << modI << " of type " << topoChanges[modI].type()
229  << abort(FatalError);
230  }
231  }
232 }
233 
234 
235 void Foam::linearValveLayersFvMesh::makeSlidersLive()
236 {
237  const polyTopoChanger& topoChanges = topoChanger_;
238 
239  // Enable sliding interface
240  forAll(topoChanges, modI)
241  {
242  if (isA<layerAdditionRemoval>(topoChanges[modI]))
243  {
244  topoChanges[modI].disable();
245  }
246  else if (isA<slidingInterface>(topoChanges[modI]))
247  {
248  topoChanges[modI].enable();
249  }
250  else
251  {
253  << "Don't know what to do with mesh modifier "
254  << modI << " of type " << topoChanges[modI].type()
255  << abort(FatalError);
256  }
257  }
258 }
259 
260 
261 bool Foam::linearValveLayersFvMesh::attached() const
262 {
263  const polyTopoChanger& topoChanges = topoChanger_;
264 
265  bool result = false;
266 
267  forAll(topoChanges, modI)
268  {
269  if (isA<slidingInterface>(topoChanges[modI]))
270  {
271  result =
272  result
273  || refCast<const slidingInterface>(topoChanges[modI]).attached();
274  }
275  }
276 
277  // Check thal all sliders are in sync (debug only)
278  forAll(topoChanges, modI)
279  {
280  if (isA<slidingInterface>(topoChanges[modI]))
281  {
282  if
283  (
284  result
285  != refCast<const slidingInterface>(topoChanges[modI]).attached()
286  )
287  {
289  << "Slider " << modI << " named "
290  << topoChanges[modI].name()
291  << " out of sync: Should be" << result
292  << abort(FatalError);
293  }
294  }
295  }
296 
297  return result;
298 }
299 
300 
301 Foam::tmp<Foam::pointField> Foam::linearValveLayersFvMesh::newPoints() const
302 {
303  tmp<pointField> tnewPoints
304  (
305  new pointField(points())
306  );
307 
308  pointField& np = tnewPoints();
309 
310  const word layerPatchName
311  (
312  motionDict_.subDict("layer").lookup("patch")
313  );
314 
315  const polyPatch& layerPatch = boundaryMesh()[layerPatchName];
316 
317  const labelList& patchPoints = layerPatch.meshPoints();
318 
319  const vector vel
320  (
321  motionDict_.lookup("pistonVelocity")
322  );
323 
324  forAll(patchPoints, ppI)
325  {
326  np[patchPoints[ppI]] += vel*time().deltaTValue();
327  }
328 
329  return tnewPoints;
330 }
331 
332 
333 
334 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
335 
337 :
338  topoChangerFvMesh(io),
339  motionDict_(dynamicMeshDict().optionalSubDict(typeName + "Coeffs"))
340 {
341  addZonesAndModifiers();
342 }
343 
344 
345 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
346 
348 {}
349 
350 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
351 
353 {
354  // Detaching the interface
355  if (attached())
356  {
357  Info<< "Decoupling sliding interfaces" << endl;
358  makeSlidersLive();
359 
360  // Changing topology
361  resetMorph();
362  setMorphTimeIndex(3*time().timeIndex());
363  updateMesh();
364  }
365  else
366  {
367  Info<< "Sliding interfaces decoupled" << endl;
368  }
369 
370  // Perform layer action and mesh motion
371  makeLayersLive();
372 
373  // Changing topology
374  resetMorph();
375  setMorphTimeIndex(3*time().timeIndex() + 1);
376  updateMesh();
377 
378  if (topoChangeMap.valid())
379  {
380  if (topoChangeMap().hasMotionPoints())
381  {
382  Info<< "Topology change; executing pre-motion" << endl;
383  movePoints(topoChangeMap().preMotionPoints());
384  }
385  }
386 
387  // Move points
388  movePoints(newPoints());
389 
390  // Attach the interface
391  Info<< "Coupling sliding interfaces" << endl;
392  makeSlidersLive();
393 
394  // Changing topology
395  resetMorph();
396  setMorphTimeIndex(3*time().timeIndex() + 2);
397  updateMesh();
398 
399  // Info<< "Moving points post slider attach" << endl;
400  // const pointField p = allPoints();
401  // movePoints(p);
402 
403  Info<< "Sliding interfaces coupled: " << attached() << endl;
404 }
405 
406 
407 // ************************************************************************* //
virtual bool update()
Update the mesh for both mesh motion and topology change.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool movePoints()
Do what is necessary if the mesh has moved.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
Macros for easy insertion into run-time selection tables.
virtual ~linearValveLayersFvMesh()
Destructor.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:795
const pointField & points
Abstract base class for a topology changing fvMesh.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
linearValveLayersFvMesh(const IOobject &io)
Construct from database.
label timeIndex
Definition: getTimeIndex.H:4
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.