patchPatchDist.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 "patchPatchDist.H"
30 #include "PatchEdgeFaceWave.H"
31 #include "syncTools.H"
32 #include "polyMesh.H"
33 #include "patchEdgeFaceInfo.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const polyPatch& patch,
40  const labelHashSet& nbrPatchIDs
41 )
42 :
43  patch_(patch),
44  nbrPatchIDs_(nbrPatchIDs),
45  nUnset_(0)
46 {
48 }
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
54 {
55  // Mark all edge connected to a nbrPatch.
56  label nBnd = 0;
57  for (const label nbrPatchi : nbrPatchIDs_)
58  {
59  const polyPatch& nbrPatch = patch_.boundaryMesh()[nbrPatchi];
60  nBnd += nbrPatch.nEdges()-nbrPatch.nInternalEdges();
61  }
62 
63  // Mark all edges. Note: should use HashSet but have no syncTools
64  // functionality for these.
65  EdgeMap<label> nbrEdges(2*nBnd);
66 
67  for (const label nbrPatchi : nbrPatchIDs_)
68  {
69  const polyPatch& nbrPatch = patch_.boundaryMesh()[nbrPatchi];
70  const labelList& nbrMp = nbrPatch.meshPoints();
71 
72  for
73  (
74  label edgeI = nbrPatch.nInternalEdges();
75  edgeI < nbrPatch.nEdges();
76  ++edgeI
77  )
78  {
79  const edge& e = nbrPatch.edges()[edgeI];
80  const edge meshE = edge(nbrMp[e[0]], nbrMp[e[1]]);
81  nbrEdges.insert(meshE, nbrPatchi);
82  }
83  }
84 
85 
86  // Make sure these boundary edges are marked everywhere.
88  (
89  patch_.boundaryMesh().mesh(),
90  nbrEdges,
91  maxEqOp<label>()
92  );
93 
94 
95  // Data on all edges and faces
96  List<patchEdgeFaceInfo> allEdgeInfo(patch_.nEdges());
97  List<patchEdgeFaceInfo> allFaceInfo(patch_.size());
98 
99  // Initial seed
100  label nBndEdges = patch_.nEdges() - patch_.nInternalEdges();
101  DynamicList<label> initialEdges(2*nBndEdges);
102  DynamicList<patchEdgeFaceInfo> initialEdgesInfo(2*nBndEdges);
103 
104 
105  // Seed all my edges that are also nbrEdges
106 
107  const labelList& mp = patch_.meshPoints();
108 
109  for
110  (
111  label edgeI = patch_.nInternalEdges();
112  edgeI < patch_.nEdges();
113  edgeI++
114  )
115  {
116  const edge& e = patch_.edges()[edgeI];
117  const edge meshE = edge(mp[e[0]], mp[e[1]]);
118 
119  if (nbrEdges.found(meshE))
120  {
121  initialEdges.append(edgeI);
122  initialEdgesInfo.append
123  (
124  patchEdgeFaceInfo
125  (
126  e.centre(patch_.localPoints()),
127  0.0
128  )
129  );
130  }
131  }
132 
133 
134  // Walk
135  PatchEdgeFaceWave
136  <
138  patchEdgeFaceInfo
139  > calc
140  (
141  patch_.boundaryMesh().mesh(),
142  patch_,
143  initialEdges,
144  initialEdgesInfo,
145  allEdgeInfo,
146  allFaceInfo,
147  returnReduce(patch_.nEdges(), sumOp<label>())
148  );
149 
150 
151  // Extract into *this
152  setSize(patch_.size());
153  nUnset_ = 0;
154  forAll(allFaceInfo, facei)
155  {
156  if (allFaceInfo[facei].valid(calc.data()))
157  {
158  operator[](facei) = Foam::sqrt(allFaceInfo[facei].distSqr());
159  }
160  else
161  {
162  nUnset_++;
163  }
164  }
165 }
166 
167 
168 // ************************************************************************* //
const Field< point_type > & localPoints() const
Return pointField of points in patch.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label nInternalEdges() const
Number of internal edges.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:361
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: polyPatch.C:309
const polyMesh & mesh() const noexcept
Return the mesh reference.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nEdges() const
Number of edges in patch.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
patchPatchDist(const polyPatch &pp, const labelHashSet &nbrPatchIDs)
Construct from patch and neighbour patches.
virtual void correct()
Correct for mesh geom/topo changes.
const std::string patch
OpenFOAM patch number as a std::string.
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
static void syncEdgeMap(const polyMesh &mesh, EdgeMap< T > &edgeValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected edges.
const dimensionedScalar mp
Proton mass.