Merge branch 'upstream/OpenFOAM' into pu
[freefoam.git] / src / topoChangerFvMesh / linearValveLayersFvMesh / linearValveLayersFvMesh.C
blob3082aca3ee86996f7645083067ddced94bb64500
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "linearValveLayersFvMesh.H"
27 #include <OpenFOAM/Time.H>
28 #include <dynamicMesh/slidingInterface.H>
29 #include <dynamicMesh/layerAdditionRemoval.H>
30 #include <OpenFOAM/pointField.H>
31 #include <OpenFOAM/mapPolyMesh.H>
32 #include <dynamicMesh/polyTopoChange.H>
33 #include <OpenFOAM/addToRunTimeSelectionTable.H>
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 namespace Foam
39     defineTypeNameAndDebug(linearValveLayersFvMesh, 0);
41     addToRunTimeSelectionTable(topoChangerFvMesh, linearValveLayersFvMesh, IOobject);
45 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
47 void Foam::linearValveLayersFvMesh::addZonesAndModifiers()
49     // Add zones and modifiers for motion action
51     if
52     (
53         pointZones().size()
54      || faceZones().size()
55      || cellZones().size()
56      || topoChanger_.size()
57     )
58     {
59         Info<< "void linearValveLayersFvMesh::addZonesAndModifiers() : "
60             << "Zones and modifiers already present.  Skipping."
61             << endl;
63         return;
64     }
66     Info<< "Time = " << time().timeName() << endl
67         << "Adding zones and modifiers to the mesh" << endl;
69     // Add zones
70     List<pointZone*> pz(1);
71     List<faceZone*> fz(4);
72     List<cellZone*> cz(0);
75     // Add an empty zone for cut points
77     pz[0] = new pointZone
78     (
79         "cutPointZone",
80         labelList(0),
81         0,
82         pointZones()
83     );
86     // Do face zones for slider
88     // Inner slider
89     const word innerSliderName(motionDict_.subDict("slider").lookup("inside"));
90     const polyPatch& innerSlider =
91         boundaryMesh()[boundaryMesh().findPatchID(innerSliderName)];
93     labelList isf(innerSlider.size());
95     forAll (isf, i)
96     {
97         isf[i] = innerSlider.start() + i;
98     }
100     fz[0] = new faceZone
101     (
102         "insideSliderZone",
103         isf,
104         boolList(innerSlider.size(), false),
105         0,
106         faceZones()
107     );
109     // Outer slider
110     const word outerSliderName(motionDict_.subDict("slider").lookup("outside"));
111     const polyPatch& outerSlider =
112         boundaryMesh()[boundaryMesh().findPatchID(outerSliderName)];
114     labelList osf(outerSlider.size());
116     forAll (osf, i)
117     {
118         osf[i] = outerSlider.start() + i;
119     }
121     fz[1] = new faceZone
122     (
123         "outsideSliderZone",
124         osf,
125         boolList(outerSlider.size(), false),
126         1,
127         faceZones()
128     );
130     // Add empty zone for cut faces
131     fz[2] = new faceZone
132     (
133         "cutFaceZone",
134         labelList(0),
135         boolList(0, false),
136         2,
137         faceZones()
138     );
140     // Add face zone for layer addition
141     const word layerPatchName
142     (
143         motionDict_.subDict("layer").lookup("patch")
144     );
146     const polyPatch& layerPatch =
147         boundaryMesh()[boundaryMesh().findPatchID(layerPatchName)];
149     labelList lpf(layerPatch.size());
151     forAll (lpf, i)
152     {
153         lpf[i] = layerPatch.start() + i;
154     }
156     fz[3] = new faceZone
157     (
158         "valveLayerZone",
159         lpf,
160         boolList(layerPatch.size(), true),
161         0,
162         faceZones()
163     );
166     Info << "Adding point and face zones" << endl;
167     addZones(pz, fz, cz);
169     // Add a topology modifier
171     List<polyMeshModifier*> tm(2);
173     tm[0] = new slidingInterface
174     (
175         "valveSlider",
176         0,
177         topoChanger_,
178         outerSliderName + "Zone",
179         innerSliderName + "Zone",
180         "cutPointZone",
181         "cutFaceZone",
182         outerSliderName,
183         innerSliderName,
184         slidingInterface::INTEGRAL,
185         true                          // Attach-detach action
186     );
188     tm[1] =
189         new layerAdditionRemoval
190         (
191             "valveLayer",
192             1,
193             topoChanger_,
194             "valveLayerZone",
195             readScalar
196             (
197                 motionDict_.subDict("layer").lookup("minThickness")
198             ),
199             readScalar
200             (
201                 motionDict_.subDict("layer").lookup("maxThickness")
202             )
203         );
206     Info << "Adding topology modifiers" << endl;
207     addTopologyModifiers(tm);
209     // Write mesh
210     write();
214 void Foam::linearValveLayersFvMesh::makeLayersLive()
216     const polyTopoChanger& topoChanges = topoChanger_;
218     // Enable layering
219     forAll (topoChanges, modI)
220     {
221         if (isA<layerAdditionRemoval>(topoChanges[modI]))
222         {
223             topoChanges[modI].enable();
224         }
225         else if (isA<slidingInterface>(topoChanges[modI]))
226         {
227             topoChanges[modI].disable();
228         }
229         else
230         {
231             FatalErrorIn("void linearValveLayersFvMesh::makeLayersLive()")
232                 << "Don't know what to do with mesh modifier "
233                 << modI << " of type " << topoChanges[modI].type()
234                 << abort(FatalError);
235         }
236     }
240 void Foam::linearValveLayersFvMesh::makeSlidersLive()
242     const polyTopoChanger& topoChanges = topoChanger_;
244     // Enable sliding interface
245     forAll (topoChanges, modI)
246     {
247         if (isA<layerAdditionRemoval>(topoChanges[modI]))
248         {
249             topoChanges[modI].disable();
250         }
251         else if (isA<slidingInterface>(topoChanges[modI]))
252         {
253             topoChanges[modI].enable();
254         }
255         else
256         {
257             FatalErrorIn("void linearValveLayersFvMesh::makeLayersLive()")
258                 << "Don't know what to do with mesh modifier "
259                 << modI << " of type " << topoChanges[modI].type()
260                 << abort(FatalError);
261         }
262     }
266 bool Foam::linearValveLayersFvMesh::attached() const
268     const polyTopoChanger& topoChanges = topoChanger_;
270     bool result = false;
272     forAll (topoChanges, modI)
273     {
274         if (isA<slidingInterface>(topoChanges[modI]))
275         {
276             result =
277                 result
278              || refCast<const slidingInterface>(topoChanges[modI]).attached();
279         }
280     }
282     // Check thal all sliders are in sync (debug only)
283     forAll (topoChanges, modI)
284     {
285         if (isA<slidingInterface>(topoChanges[modI]))
286         {
287             if
288             (
289                 result
290              != refCast<const slidingInterface>(topoChanges[modI]).attached()
291             )
292             {
293                 FatalErrorIn("bool linearValveLayersFvMesh::attached() const")
294                     << "Slider " << modI << " named "
295                     << topoChanges[modI].name()
296                     << " out of sync: Should be" << result
297                     << abort(FatalError);
298             }
299         }
300     }
302     return result;
306 Foam::tmp<Foam::pointField> Foam::linearValveLayersFvMesh::newPoints() const
308     tmp<pointField> tnewPoints
309     (
310         new pointField(points())
311     );
313     pointField& np = tnewPoints();
315     const word layerPatchName
316     (
317         motionDict_.subDict("layer").lookup("patch")
318     );
320     const polyPatch& layerPatch =
321         boundaryMesh()[boundaryMesh().findPatchID(layerPatchName)];
323     const labelList& patchPoints = layerPatch.meshPoints();
325     const vector vel
326     (
327         motionDict_.lookup("pistonVelocity")
328     );
330     forAll (patchPoints, ppI)
331     {
332         np[patchPoints[ppI]] += vel*time().deltaT().value();
333     }
335     return tnewPoints;
340 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
342 // Construct from components
343 Foam::linearValveLayersFvMesh::linearValveLayersFvMesh(const IOobject& io)
345     topoChangerFvMesh(io),
346     motionDict_
347     (
348         IOdictionary
349         (
350             IOobject
351             (
352                 "dynamicMeshDict",
353                 time().constant(),
354                 *this,
355                 IOobject::MUST_READ,
356                 IOobject::NO_WRITE
357             )
358         ).subDict(typeName + "Coeffs")
359     )
361     addZonesAndModifiers();
365 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
367 Foam::linearValveLayersFvMesh::~linearValveLayersFvMesh()
370 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
372 void Foam::linearValveLayersFvMesh::update()
374     // Detaching the interface
375     if (attached())
376     {
377         Info << "Decoupling sliding interfaces" << endl;
378         makeSlidersLive();
380         // Changing topology
381         resetMorph();
382         setMorphTimeIndex(3*time().timeIndex());
383         updateMesh();
384     }
385     else
386     {
387         Info << "Sliding interfaces decoupled" << endl;
388     }
390     // Perform layer action and mesh motion
391     makeLayersLive();
393     // Changing topology
394     resetMorph();
395     setMorphTimeIndex(3*time().timeIndex() + 1);
396     updateMesh();
398     if (topoChangeMap.valid())
399     {
400         if (topoChangeMap().hasMotionPoints())
401         {
402             Info << "Topology change; executing pre-motion" << endl;
403             movePoints(topoChangeMap().preMotionPoints());
404         }
405     }
407     // Move points
408     movePoints(newPoints());
410     // Attach the interface
411     Info << "Coupling sliding interfaces" << endl;
412     makeSlidersLive();
414     // Changing topology
415     resetMorph();
416     setMorphTimeIndex(3*time().timeIndex() + 2);
417     updateMesh();
419     Info << "Moving points post slider attach" << endl;
420 //     const pointField p = allPoints();
421 //     movePoints(p);
423     Info << "Sliding interfaces coupled: " << attached() << endl;
427 // ************************ vim: set sw=4 sts=4 et: ************************ //