externalDisplacementMeshMover.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2015 OpenFOAM Foundation
9 Copyright (C) 2015-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "mapPolyMesh.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
41
42// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43
45(
47)
48{
49 DynamicList<label> adaptPatchIDs;
50
51 forAll(field.boundaryField(), patchI)
52 {
53 const pointPatchField<vector>& patchField =
54 field.boundaryField()[patchI];
55
56 if (isA<valuePointPatchField<vector>>(patchField))
57 {
59 {
60 // Special condition of fixed boundary condition. Does not
61 // get adapted
62 }
63 else
64 {
65 adaptPatchIDs.append(patchI);
66 }
67 }
68 }
69
70 return adaptPatchIDs;
71}
72
73
76(
77 const polyMesh& mesh,
78 const labelList& patchIDs
79)
80{
82
83 // Count faces.
84 label nFaces = 0;
85
86 forAll(patchIDs, i)
87 {
88 const polyPatch& pp = patches[patchIDs[i]];
89
90 nFaces += pp.size();
91 }
92
93 // Collect faces.
94 labelList addressing(nFaces);
95 nFaces = 0;
96
97 forAll(patchIDs, i)
98 {
99 const polyPatch& pp = patches[patchIDs[i]];
100
101 label meshFaceI = pp.start();
102
103 forAll(pp, i)
104 {
105 addressing[nFaces++] = meshFaceI++;
106 }
107 }
108
110 (
111 IndirectList<face>(mesh.faces(), addressing),
112 mesh.points()
113 );
114}
115
116
117// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118
120(
121 const dictionary& dict,
122 const List<labelPair>& baffles,
123 pointVectorField& pointDisplacement,
124 const bool dryRun
125)
126:
127 baffles_(baffles),
128 pointDisplacement_(pointDisplacement),
129 dryRun_(dryRun)
130{}
131
132
133// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
134
137(
138 const word& type,
139 const dictionary& dict,
140 const List<labelPair>& baffles,
141 pointVectorField& pointDisplacement,
142 const bool dryRun
143)
144{
145 Info<< "Selecting externalDisplacementMeshMover " << type << endl;
146
147 auto* ctorPtr = dictionaryConstructorTable(type);
148
149 if (!ctorPtr)
150 {
152 (
153 dict,
154 "externalDisplacementMeshMover",
155 type,
156 *dictionaryConstructorTablePtr_
157 ) << exit(FatalIOError);
158 }
159
161 (
162 ctorPtr(dict, baffles, pointDisplacement, dryRun)
163 );
164}
165
166
167// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
168
170{}
171
172
173// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
174
176{
177 // No local data to update
178}
179
180
182{
183 // Renumber baffles
184 DynamicList<labelPair> newBaffles(baffles_.size());
185 forAll(baffles_, i)
186 {
187 label f0 = mpm.reverseFaceMap()[baffles_[i].first()];
188 label f1 = mpm.reverseFaceMap()[baffles_[i].second()];
189
190 if (f0 >= 0 && f1 >= 0)
191 {
192 newBaffles.append(labelPair(f0, f1));
193 }
194 }
195 newBaffles.shrink();
196 baffles_.transfer(newBaffles);
197}
198
199
200// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A List with indirect addressing.
Definition: IndirectList.H:119
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Virtual base class for mesh movers with externally provided displacement field giving the boundary co...
static labelList getFixedValueBCs(const pointVectorField &)
Extract fixed-value patchfields.
static autoPtr< indirectPrimitivePatch > getPatch(const polyMesh &, const labelList &)
Construct patch on selected patches.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Abstract base class for point-mesh patch fields.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
Foam::valuePointPatchField.
A class for handling words, derived from Foam::string.
Definition: word.H:68
Enables the specification of a zero fixed value boundary condition.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
rDeltaTY field()
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:57
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333