linearValveLayersFvMesh.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 | Copyright (C) 2011-2017 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  readScalar
197  (
198  motionDict_.subDict("layer").lookup("minThickness")
199  ),
200  readScalar
201  (
202  motionDict_.subDict("layer").lookup("maxThickness")
203  )
204  );
205 
206 
207  Info<< "Adding topology modifiers" << endl;
208  addTopologyModifiers(tm);
209 
210  // Write mesh
211  write();
212 }
213 
214 
215 void Foam::linearValveLayersFvMesh::makeLayersLive()
216 {
217  const polyTopoChanger& topoChanges = topoChanger_;
218 
219  // Enable layering
220  forAll(topoChanges, modI)
221  {
222  if (isA<layerAdditionRemoval>(topoChanges[modI]))
223  {
224  topoChanges[modI].enable();
225  }
226  else if (isA<slidingInterface>(topoChanges[modI]))
227  {
228  topoChanges[modI].disable();
229  }
230  else
231  {
233  << "Don't know what to do with mesh modifier "
234  << modI << " of type " << topoChanges[modI].type()
235  << abort(FatalError);
236  }
237  }
238 }
239 
240 
241 void Foam::linearValveLayersFvMesh::makeSlidersLive()
242 {
243  const polyTopoChanger& topoChanges = topoChanger_;
244 
245  // Enable sliding interface
246  forAll(topoChanges, modI)
247  {
248  if (isA<layerAdditionRemoval>(topoChanges[modI]))
249  {
250  topoChanges[modI].disable();
251  }
252  else if (isA<slidingInterface>(topoChanges[modI]))
253  {
254  topoChanges[modI].enable();
255  }
256  else
257  {
259  << "Don't know what to do with mesh modifier "
260  << modI << " of type " << topoChanges[modI].type()
261  << abort(FatalError);
262  }
263  }
264 }
265 
266 
267 bool Foam::linearValveLayersFvMesh::attached() const
268 {
269  const polyTopoChanger& topoChanges = topoChanger_;
270 
271  bool result = false;
272 
273  forAll(topoChanges, modI)
274  {
275  if (isA<slidingInterface>(topoChanges[modI]))
276  {
277  result =
278  result
279  || refCast<const slidingInterface>(topoChanges[modI]).attached();
280  }
281  }
282 
283  // Check thal all sliders are in sync (debug only)
284  forAll(topoChanges, modI)
285  {
286  if (isA<slidingInterface>(topoChanges[modI]))
287  {
288  if
289  (
290  result
291  != refCast<const slidingInterface>(topoChanges[modI]).attached()
292  )
293  {
295  << "Slider " << modI << " named "
296  << topoChanges[modI].name()
297  << " out of sync: Should be" << result
298  << abort(FatalError);
299  }
300  }
301  }
302 
303  return result;
304 }
305 
306 
307 Foam::tmp<Foam::pointField> Foam::linearValveLayersFvMesh::newPoints() const
308 {
309  tmp<pointField> tnewPoints
310  (
311  new pointField(points())
312  );
313 
314  pointField& np = tnewPoints();
315 
316  const word layerPatchName
317  (
318  motionDict_.subDict("layer").lookup("patch")
319  );
320 
321  const polyPatch& layerPatch = boundaryMesh()[layerPatchName];
322 
323  const labelList& patchPoints = layerPatch.meshPoints();
324 
325  const vector vel
326  (
327  motionDict_.lookup("pistonVelocity")
328  );
329 
330  forAll(patchPoints, ppI)
331  {
332  np[patchPoints[ppI]] += vel*time().deltaTValue();
333  }
334 
335  return tnewPoints;
336 }
337 
338 
339 
340 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341 
342 // Construct from components
343 Foam::linearValveLayersFvMesh::linearValveLayersFvMesh(const IOobject& io)
344 :
345  topoChangerFvMesh(io),
346  motionDict_
347  (
349  (
350  IOobject
351  (
352  "dynamicMeshDict",
353  time().constant(),
354  *this,
355  IOobject::MUST_READ_IF_MODIFIED,
356  IOobject::NO_WRITE,
357  false
358  )
359  ).optionalSubDict(typeName + "Coeffs")
360  )
361 {
362  addZonesAndModifiers();
363 }
364 
365 
366 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
367 
369 {}
370 
371 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
372 
374 {
375  // Detaching the interface
376  if (attached())
377  {
378  Info<< "Decoupling sliding interfaces" << endl;
379  makeSlidersLive();
380 
381  // Changing topology
382  resetMorph();
383  setMorphTimeIndex(3*time().timeIndex());
384  updateMesh();
385  }
386  else
387  {
388  Info<< "Sliding interfaces decoupled" << endl;
389  }
390 
391  // Perform layer action and mesh motion
392  makeLayersLive();
393 
394  // Changing topology
395  resetMorph();
396  setMorphTimeIndex(3*time().timeIndex() + 1);
397  updateMesh();
398 
399  if (topoChangeMap.valid())
400  {
401  if (topoChangeMap().hasMotionPoints())
402  {
403  Info<< "Topology change; executing pre-motion" << endl;
404  movePoints(topoChangeMap().preMotionPoints());
405  }
406  }
407 
408  // Move points
409  movePoints(newPoints());
410 
411  // Attach the interface
412  Info<< "Coupling sliding interfaces" << endl;
413  makeSlidersLive();
414 
415  // Changing topology
416  resetMorph();
417  setMorphTimeIndex(3*time().timeIndex() + 2);
418  updateMesh();
419 
420  //Info<< "Moving points post slider attach" << endl;
421  //const pointField p = allPoints();
422  //movePoints(p);
423 
424  Info<< "Sliding interfaces coupled: " << attached() << endl;
425 }
426 
427 
428 // ************************************************************************* //
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:428
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:253
bool movePoints()
Do what is neccessary 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:243
Macros for easy insertion into run-time selection tables.
virtual ~linearValveLayersFvMesh()
Destructor.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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:800
const pointField & points
Abstract base class for a topology changing fvMesh.
List< label > labelList
A List of labels.
Definition: labelList.H:56
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
errorManip< error > abort(error &err)
Definition: errorManip.H:131
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
Constant dispersed-phase particle diameter model.
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.