Loading...
Searching...
No Matches
springRenumber.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2024 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
29#include "springRenumber.H"
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
38
40 (
44 );
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
50Foam::springRenumber::springRenumber(const dictionary& dict)
51:
53 coeffsDict_(dict.optionalSubDict(typeName + "Coeffs")),
54 maxIter_(coeffsDict_.get<label>("maxIter")),
55 maxCo_(coeffsDict_.get<scalar>("maxCo")),
56 freezeFraction_(coeffsDict_.get<scalar>("freezeFraction")),
57 verbose_(coeffsDict_.getOrDefault("verbose", true))
58{}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
63template<class ConnectionListListType>
64Foam::labelList Foam::springRenumber::renumberImpl
65(
66 const ConnectionListListType& cellCells
67) const
68{
69 const label nOldCells(cellCells.size());
70
71 // Look at cell index as a 1D position parameter.
72 // Move cells to the average 'position' of their neighbour.
73
74 scalarField position(nOldCells);
75 std::iota(position.begin(), position.end(), 0);
76
77 // Sum force per cell. Also reused for the displacement
78 scalarField sumForce(nOldCells);
79
80 labelList oldToNew(Foam::identity(nOldCells));
81
82 scalar maxCo = (maxCo_ * nOldCells);
83
84 for (label iter = 0; iter < maxIter_; ++iter)
85 {
86 //Pout<< "Iteration : " << iter << nl
87 // << "------------"
88 // << endl;
89
90 //Pout<< "Position :" << nl
91 // << " min : " << Foam::min(position) << nl
92 // << " max : " << Foam::max(position) << nl
93 // << " avg : " << Foam::average(position) << nl
94 // << endl;
95
96 // Sum force per cell.
97 sumForce = Zero;
98 for (label oldCelli = 0; oldCelli < nOldCells; ++oldCelli)
99 {
100 const label celli = oldToNew[oldCelli];
101
102 for (const label nbr : cellCells[oldCelli])
103 {
104 const label nbrCelli = oldToNew[nbr];
105
106 sumForce[celli] += (position[nbrCelli]-position[celli]);
107 }
108 }
109
110 //Pout<< "Force :" << nl
111 // << " min : " << Foam::min(sumForce) << nl
112 // << " max : " << Foam::max(sumForce) << nl
113 // << " avgMag : " << Foam::average(mag(sumForce)) << nl
114 // << "DeltaT : " << deltaT << nl
115 // << endl;
116
117 // Limit displacement
118 scalar deltaT = maxCo / max(mag(sumForce));
119
120 if (verbose_)
121 {
122 Info<< "Iter:" << iter
123 << " maxCo:" << maxCo
124 << " deltaT:" << deltaT
125 << " average force:" << average(mag(sumForce)) << endl;
126 }
127
128 // Determine displacement
129 scalarField& displacement = sumForce;
130 displacement *= deltaT;
131
132 //Pout<< "Displacement :" << nl
133 // << " min : " << Foam::min(displacement) << nl
134 // << " max : " << Foam::max(displacement) << nl
135 // << " avgMag : " << Foam::average(mag(displacement)) << nl
136 // << endl;
137
138 // Calculate new position and scale to be within original range
139 // (0..nCells-1) for ease of postprocessing.
140 position += displacement;
141 position -= Foam::min(position);
142 position *= (position.size()-1)/Foam::max(position);
143
144 // Slowly freeze.
145 maxCo *= freezeFraction_;
146 }
147
148 //writeOBJ("endPosition.obj", cellCells, position);
149
150 // Move cells to new position
152
153 // Reorder oldToNew
154 inplaceReorder(shuffle, oldToNew);
155
156 return invert(nOldCells, oldToNew);
157}
158
159
161(
162 const polyMesh& mesh
163) const
164{
165 // Local mesh connectivity
168
169 return renumberImpl(cellCells);
170}
171
172
174(
175 const CompactListList<label>& cellCells
176) const
177{
178 return renumberImpl(cellCells);
179}
180
181
183(
184 const labelListList& cellCells
185) const
186{
187 return renumberImpl(cellCells);
188}
189
190
191// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A packed storage of objects of type <T> using an offset table for access.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition UListI.H:454
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static void calcCellCells(const polyMesh &mesh, const labelUList &agglom, const label nLocalCoarse, const bool parallel, CompactListList< label > &cellCells)
Determine (local or global) cellCells from mesh agglomeration.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Abstract base class for renumbering.
renumberMethod()
Default construct.
Use spring analogy - attract neighbouring cells according to the distance of their cell indices.
virtual labelList renumber(const polyMesh &mesh) const
Return the cell visit order (from ordered back to original cell id) using the mesh to determine the c...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
scalar maxCo
Namespace for OpenFOAM.
void shuffle(UList< T > &list)
Randomise the list order.
Definition UList.C:315
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dictionary dict