displacementLayeredMotionSolver.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) 2024 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 "PointEdgeWave.H"
29 #include "pointConstraints.H"
30 #include "polyTopoChangeMap.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 
40  (
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::displacementLayeredMotionSolver::walkLayers
51 (
52  const polyPatch& startPatch,
53  scalarField& distance,
54  vectorField& displacement
55 ) const
56 {
57  // Get the start mesh points from the start patch
58  const labelList& startMeshPoints = startPatch.meshPoints();
59 
60  // Set the displacement of the start points
61  const vectorField startDisplacement
62  (
63  pointDisplacement_.boundaryField()[startPatch.index()]
64  .patchInternalField()
65  );
66 
67  // Set the wave info for the start patch
68  List<pointEdgeStructuredWalk> startInfo(startMeshPoints.size());
69  forAll(startMeshPoints, i)
70  {
71  startInfo[i] = pointEdgeStructuredWalk
72  (
73  points0()[startMeshPoints[i]], // Location of start point
74  points0()[startMeshPoints[i]], // Previous location
75  0,
76  startDisplacement[i] // Displacement of the start point
77  );
78  }
79 
80  // Initialise the wave info for all the points
81  List<pointEdgeStructuredWalk> allPointInfo(mesh().nPoints());
82  forAll(allPointInfo, pointi)
83  {
84  allPointInfo[pointi] = pointEdgeStructuredWalk
85  (
86  points0()[pointi], // Mesh point location
87  vector::max, // No valid previous location
88  0,
89  Zero // Initial displacement = 0
90  );
91  }
92 
93  // Initialise the wave info for all edges
94  List<pointEdgeStructuredWalk> allEdgeInfo(mesh().nEdges());
95  forAll(allEdgeInfo, edgei)
96  {
97  allEdgeInfo[edgei] = pointEdgeStructuredWalk
98  (
99  mesh().edges()[edgei].centre(points0()), // Edge centre location
100  vector::max, // No valid previous location
101  0,
102  Zero // Initial displacement = 0
103  );
104  }
105 
106  // Walk the distance and displacement from the startPatch
107  PointEdgeWave<pointEdgeStructuredWalk> walk
108  (
109  mesh(),
110  startMeshPoints,
111  startInfo,
112  allPointInfo,
113  allEdgeInfo,
114  mesh().globalData().nTotalPoints() // Max number of iterations
115  );
116 
117  // Extract distance and displacement from the wave info
118  forAll(allPointInfo, pointi)
119  {
120  distance[pointi] = allPointInfo[pointi].dist();
121  displacement[pointi] = allPointInfo[pointi].data();
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
129 (
130  const word& name,
131  const polyMesh& mesh,
132  const dictionary& dict
133 )
134 :
135  displacementMotionSolver(name, mesh, dict, typeName),
136  oppositePatchNames_(dict.lookup("oppositePatches")),
137  oppositePatches_
138  (
139  mesh.boundaryMesh().findIndex(oppositePatchNames_.first()),
140  mesh.boundaryMesh().findIndex(oppositePatchNames_.second())
141  )
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146 
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
155 {
156  tmp<pointField> tcurPoints
157  (
158  points0() + pointDisplacement_.primitiveField()
159  );
160 
161  return tcurPoints;
162 }
163 
164 
166 {
167  // The points have moved so before interpolation update the motionSolver
168  movePoints(mesh().points());
169 
170  // Update the displacement boundary conditions
171  pointDisplacement_.boundaryFieldRef().updateCoeffs();
172 
173  // Walk the layers from patch0 to patch1
174  const polyPatch& patch0 = mesh().boundaryMesh()[oppositePatches_.first()];
175  scalarField patchDist0(mesh().nPoints());
176  vectorField patchDisp0(pointDisplacement_);
177  walkLayers(patch0, patchDist0, patchDisp0);
178 
179  // Walk the layers from patch1 to patch0
180  const polyPatch& patch1 = mesh().boundaryMesh()[oppositePatches_.second()];
181  scalarField patchDist1(mesh().nPoints());
182  vectorField patchDisp1(pointDisplacement_);
183  walkLayers(patch1, patchDist1, patchDisp1);
184 
185  // Calculate the interpolation factor from the distance to the
186  // opposite patches
187  const scalarField w(patchDist0/(patchDist0 + patchDist1 + small));
188 
189  // Linearly interpolate the displacements of the opposite patches
190  pointDisplacement_.primitiveFieldRef() = (1 - w)*patchDisp0 + w*patchDisp1;
191 
192  // Constrain the pointDisplacement field
193  pointConstraints::New(pointDisplacement_.mesh())
194  .constrainDisplacement(pointDisplacement_, false);
195 }
196 
197 
199 (
200  const polyTopoChangeMap& map
201 )
202 {
203  // Pending implementation of the inverse transformation of points0
205 
206  // Alternatively using the current points instead of points0
207  // would support topology change but the motion might be less robust
208  // and potentially accumulate error
209 }
210 
211 
212 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
static pointConstraints & New(const word &name, const pointMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const Form max
Definition: VectorSpace.H:120
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Interpolating motion solver for extruded/layered meshes.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void topoChange(const polyTopoChangeMap &)
Update topology (not implemented)
displacementLayeredMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
Virtual base class for displacement motion solver.
pointVectorField pointDisplacement_
Point motion field.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:133
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
pointField & points0()
Return reference to the reference field.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
const pointField & points
label nPoints
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
defineTypeNameAndDebug(combustionModel, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
dictionary dict