edgeVertex.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-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "edgeVertex.H"
30 #include "meshTools.H"
31 #include "refineCell.H"
32 
33 
34 // * * * * * * * * * * * * * * * Static Functions * * * * * * * * * * * * * //
35 
36 // Update stored refine list using map
38 (
39  const labelList& map,
40  List<refineCell>& refCells
41 )
42 {
43  label newRefI = 0;
44 
45  forAll(refCells, refI)
46  {
47  const refineCell& refCell = refCells[refI];
48 
49  label oldCelli = refCell.cellNo();
50 
51  label newCelli = map[oldCelli];
52 
53  if (newCelli != -1)
54  {
55  refCells[newRefI++] = refineCell(newCelli, refCell.direction());
56  }
57  }
58  refCells.setSize(newRefI);
59 }
60 
61 
62 // Update stored cell numbers using map.
63 // Do in two passes to prevent allocation if nothing changed.
65 (
66  const labelList& map,
67  Map<label>& cellPairs
68 )
69 {
70  // Iterate over map to see if anything changed
71  bool changed = false;
72 
73  forAllConstIters(cellPairs, iter)
74  {
75  label newMaster = map[iter.key()];
76 
77  label newSlave = -1;
78 
79  if (iter.val() != -1)
80  {
81  newSlave = map[iter.val()];
82  }
83 
84  if ((newMaster != iter.key()) || (newSlave != iter.val()))
85  {
86  changed = true;
87 
88  break;
89  }
90  }
91 
92  // Relabel (use second Map to prevent overlapping)
93  if (changed)
94  {
95  Map<label> newCellPairs(2*cellPairs.size());
96 
97  forAllConstIters(cellPairs, iter)
98  {
99  label newMaster = map[iter.key()];
100 
101  label newSlave = -1;
102 
103  if (iter.val() != -1)
104  {
105  newSlave = map[iter.val()];
106  }
107 
108  if (newMaster == -1)
109  {
111  << "master cell:" << iter.key()
112  << " has disappeared" << endl;
113  }
114  else
115  {
116  newCellPairs.insert(newMaster, newSlave);
117  }
118  }
119 
120  cellPairs = newCellPairs;
121  }
122 }
123 
124 
125 // Update stored cell numbers using map.
126 // Do in two passes to prevent allocation if nothing changed.
128 (
129  const labelList& map,
131 )
132 {
133  // Iterate over map to see if anything changed
134  bool changed = false;
135 
136  for (const label celli : cells)
137  {
138  const label newCelli = map[celli];
139 
140  if (newCelli != celli)
141  {
142  changed = true;
143 
144  break;
145  }
146  }
147 
148  // Relabel (use second Map to prevent overlapping)
149  if (changed)
150  {
151  labelHashSet newCells(2*cells.size());
152 
153  for (const label celli : cells)
154  {
155  const label newCelli = map[celli];
156 
157  if (newCelli != -1)
158  {
159  newCells.insert(newCelli);
160  }
161  }
162 
163  cells = newCells;
164  }
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
169 
171 (
172  const primitiveMesh& mesh,
173  const label cut,
174  const scalar weight
175 )
176 {
177  const pointField& pts = mesh.points();
178 
179  if (isEdge(mesh, cut))
180  {
181  const edge& e = mesh.edges()[getEdge(mesh, cut)];
182 
183  return weight*pts[e.end()] + (1-weight)*pts[e.start()];
184  }
185  else
186  {
187  return pts[getVertex(mesh, cut)];
188  }
189 }
190 
191 
193 (
194  const primitiveMesh& mesh,
195  const label cut0,
196  const label cut1
197 )
198 {
199  if (!isEdge(mesh, cut0) && !isEdge(mesh, cut1))
200  {
201  return meshTools::findEdge
202  (
203  mesh,
204  getVertex(mesh, cut0),
205  getVertex(mesh, cut1)
206  );
207  }
208  else
209  {
210  return -1;
211  }
212 }
213 
214 
216 (
217  Ostream& os,
218  const label cut,
219  const scalar weight
220 ) const
221 {
222  if (isEdge(cut))
223  {
224  label edgeI = getEdge(cut);
225 
226  const edge& e = mesh().edges()[edgeI];
227 
228  os << "edge:" << edgeI << e << ' ' << coord(cut, weight);
229  }
230  else
231  {
232  label vertI = getVertex(cut);
233 
234  os << "vertex:" << vertI << ' ' << coord(cut, weight);
235  }
236  return os;
237 }
238 
239 
241 (
242  Ostream& os,
243  const labelList& cuts,
244  const scalarField& weights
245 ) const
246 {
247  forAll(cuts, cutI)
248  {
249  if (cutI > 0)
250  {
251  os << ' ';
252  }
253  writeCut(os, cuts[cutI], weights[cutI]);
254  }
255  return os;
256 }
257 
258 
259 // ************************************************************************* //
meshTools.H
Foam::edgeVertex::writeCut
Ostream & writeCut(Ostream &os, const label cut, const scalar) const
Write cut description to Ostream.
Definition: edgeVertex.C:216
Foam::edgeVertex::updateLabels
static void updateLabels(const labelList &map, List< refineCell > &)
Update refine list from map. Used to update cell/face labels.
Definition: edgeVertex.C:38
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::Map< label >
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::HashSet< label, Hash< label > >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::label
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:62
Foam::meshTools::findEdge
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:359
Foam::Field< vector >
refineCell.H
Foam::edgeVertex::coord
static point coord(const primitiveMesh &, const label cut, const scalar weight)
Return coordinate of cut (uses weight if edgeCut)
Definition: edgeVertex.C:171
cut
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
Foam::refineCell::direction
const vector & direction() const
Definition: refineCell.H:89
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::refineCell
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:58
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::Vector< scalar >
Foam::edgeVertex::cutPairToEdge
static label cutPairToEdge(const primitiveMesh &, const label cut0, const label cut1)
Find mesh edge (or -1) between two cuts.
Definition: edgeVertex.C:193
Foam::List< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::refineCell::cellNo
label cellNo() const
Definition: refineCell.H:84
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::edgeVertex::writeCuts
Ostream & writeCuts(Ostream &os, const labelList &, const scalarField &) const
Write cut descriptions to Ostream.
Definition: edgeVertex.C:241
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
edgeVertex.H
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78