multiSolidBodyMotionSolver.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-2022 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 
28 #include "transformField.H"
29 #include "meshCellZones.H"
30 #include "boolList.H"
31 #include "syncTools.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
39  (
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const polyMesh& mesh,
53  const dictionary& dict
54 )
55 :
56  points0MotionSolver(name, mesh, dict, typeName)
57 {
58  zoneIDs_.setSize(coeffDict().size());
59  SBMFs_.setSize(coeffDict().size());
60  pointIDs_.setSize(coeffDict().size());
61  label zonei = 0;
62 
64  {
65  if (iter().isDict())
66  {
67  zoneIDs_[zonei] = mesh.cellZones().findZoneID(iter().keyword());
68 
69  if (zoneIDs_[zonei] == -1)
70  {
72  (
73  coeffDict()
74  ) << "Cannot find cellZone named " << iter().keyword()
75  << ". Valid zones are " << mesh.cellZones().names()
76  << exit(FatalIOError);
77  }
78 
79  const dictionary& subDict = iter().dict();
80 
81  SBMFs_.set
82  (
83  zonei,
85  );
86 
87  // Collect points of cell zone.
88  const cellZone& cz = mesh.cellZones()[zoneIDs_[zonei]];
89 
90  boolList movePts(mesh.nPoints(), false);
91 
92  forAll(cz, i)
93  {
94  label celli = cz[i];
95  const cell& c = mesh.cells()[celli];
96  forAll(c, j)
97  {
98  const face& f = mesh.faces()[c[j]];
99  forAll(f, k)
100  {
101  label pointi = f[k];
102  movePts[pointi] = true;
103  }
104  }
105  }
106 
107  syncTools::syncPointList(mesh, movePts, orEqOp<bool>(), false);
108 
110  forAll(movePts, i)
111  {
112  if (movePts[i])
113  {
114  ptIDs.append(i);
115  }
116  }
117 
118  pointIDs_[zonei].transfer(ptIDs);
119 
120  Info<< "Applying solid body motion " << SBMFs_[zonei].type()
121  << " to " << pointIDs_[zonei].size() << " points of cellZone "
122  << iter().keyword() << endl;
123 
124  zonei++;
125  }
126  }
127  zoneIDs_.setSize(zonei);
128  SBMFs_.setSize(zonei);
129  pointIDs_.setSize(zonei);
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
134 
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
142 {
143  tmp<pointField> ttransformedPts(new pointField(mesh().points()));
144  pointField& transformedPts = ttransformedPts.ref();
145 
146  forAll(zoneIDs_, i)
147  {
148  const labelList& zonePoints = pointIDs_[i];
149 
150  UIndirectList<point>(transformedPts, zonePoints) = transformPoints
151  (
152  SBMFs_[i].transformation(),
153  pointField(points0_, zonePoints)
154  );
155  }
156 
157  return ttransformedPts;
158 }
159 
160 
161 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
wordList names() const
Return a list of zone names.
Definition: MeshZones.C:256
A List with indirect addressing.
Definition: UIndirectList.H:60
A subset of mesh cells.
Definition: cellZone.H:64
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
const word & keyword() const
Return keyword.
Definition: motionSolver.H:131
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:143
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:137
Solid-body motion of the mesh specified by a run-time selectable motion function.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
multiSolidBodyMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from mesh and dictionary.
const Time & time() const
Return time.
Virtual base class for displacement motion solvers.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:451
label nPoints() const
const cellList & cells() const
static autoPtr< solidBodyMotionFunction > New(const dictionary &SBMFCoeffs, const Time &runTime)
Select constructed from the SBMFCoeffs dictionary and Time.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
const pointField & points
Foam::meshCellZones.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void transformPoints(vectorField &, const spatialTransform &, const vectorField &)
Transform given vectorField of coordinates with the given spatialTransform.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
defineTypeNameAndDebug(combustionModel, 0)
IOerror FatalIOError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
labelList f(nPoints)
dictionary dict
Spatial transformation functions for primitive fields.