Loading...
Searching...
No Matches
betaMax.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) 2020-2023 PCOpt/NTUA
9 Copyright (C) 2020-2023 FOSS GP
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 "betaMax.H"
30#include "EdgeMap.H"
31#include "syncTools.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
41
42// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43
44Foam::scalar Foam::betaMax::computeLength(const dictionary& dict) const
45{
46 scalar length = Zero;
47 // If length is not provided explicitly, loop over the provided patches
48 // and compute the hydraulic diamater
49 const dictionary& DarcyDict = dict.subDict(type() + "Coeffs");
50 if (!DarcyDict.readIfPresent("length", length))
51 {
52 const labelHashSet inletPatches =
54 (
55 DarcyDict.get<wordRes>("inletPatches")
56 );
57
58 // If 2D, use the inlet area divided by the depth in the empty direction
59 if (mesh_.nGeometricD() != label(3))
60 {
61 // Accumulate area
62 for (const label pI : inletPatches)
63 {
64 const fvPatch& patch = mesh_.boundary()[pI];
65 length += sum(patch.magSf());
66 }
67 reduce(length, sumOp<scalar>());
68
69 // Divide with the span in the empty direction
70 const Vector<label>& geometricD = mesh_.geometricD();
71 const boundBox& bounds = mesh_.bounds();
72 forAll(geometricD, iDir)
73 {
74 if (geometricD[iDir] == -1)
75 {
76 length /= bounds.span()[iDir];
77 }
78 }
79 }
80 // If 3D, use the inlet hydraulic diameter
81 else
82 {
83 scalar perimeter = Zero;
84 scalar area = Zero;
85 for (const label pI : inletPatches)
86 {
87 const polyPatch& patch = mesh_.boundaryMesh()[pI];
88 // Accumulate perimeter
89 const edgeList& edges = patch.edges();
90 const vectorField& points = patch.localPoints();
91 const label nInternalEdges = patch.nInternalEdges();
92 const label nEdges = patch.nEdges();
93 // Processor edges should not be accounted for
94 boolList isProcessorEdge = markProcessorEdges(patch);
95 label nActiveEdges(0);
96 forAll(isProcessorEdge, beI)
97 {
98 if (!isProcessorEdge[beI])
99 {
100 perimeter += edges[nInternalEdges + beI].mag(points);
101 nActiveEdges++;
102 }
103 }
104
105 if (debug > 1)
106 {
107 Info<< "patch:: " << patch.name() << nl
108 << "Number of boundary edges "
109 << returnReduce(nEdges - nInternalEdges, sumOp<label>())
110 << ", Number of non-processor edges "
111 << returnReduce(nActiveEdges, sumOp<label>()) << endl;
112 }
113
114 // Accumulate area
115 area += sum(patch.magFaceAreas());
116 }
117
118 reduce(perimeter, sumOp<scalar>());
119 reduce(area, sumOp<scalar>());
120
121 length = scalar(4)*area/perimeter;
123 }
124
125 return length;
126}
127
128
130(
131 const polyPatch& patch
132) const
133{
134 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
135 const label nNonProcessor = pbm.nNonProcessor();
136
137 // Processor edges will artificially increase the perimeter of the inlet.
138 // We need to purge them.
139
140 // Mark all edges connected to a processor patch
141 label nProcEdges(0);
142 for (label procI = nNonProcessor; procI < pbm.size() ; ++procI)
143 {
144 const polyPatch& procPatch = pbm[procI];
145 nProcEdges += procPatch.nEdges() - procPatch.nInternalEdges();
146 }
147
148 EdgeMap<bool> isInletEdge(nProcEdges);
149
150 for (label procI = nNonProcessor; procI < pbm.size() ; ++procI)
151 {
152 const polyPatch& procPatch = pbm[procI];
153 const labelList& procMp = procPatch.meshPoints();
154 const label procInternalEdges = procPatch.nInternalEdges();
155 const label procEdges = procPatch.nEdges();
156
157 for (label edgeI = procInternalEdges; edgeI < procEdges; ++edgeI)
158 {
159 const edge& e = procPatch.edges()[edgeI];
160 const edge meshE = edge(procMp[e[0]], procMp[e[1]]);
161 isInletEdge.insert(meshE, false);
162 }
163 }
164
165 // Loop over inlet edges to mark common edges with processor patches
166 const label nInternalEdges = patch.nInternalEdges();
167 const label nEdges = patch.nEdges();
168
169 const edgeList& edges = patch.edges();
170 const labelList& mp = patch.meshPoints();
171 for (label edgeI = nInternalEdges; edgeI < nEdges; ++edgeI)
172 {
173 const edge& e = edges[edgeI];
174 const edge meshE = edge(mp[e[0]], mp[e[1]]);
175 auto iter = isInletEdge.find(meshE);
176
177 if (iter.good())
178 {
179 iter.val() = true;
180 }
181 }
182
183 // A special case is a processor patch intersecting the inlet patches on a
184 // (true) boundary edge of the latter. Use syncEdgeMap to make sure these
185 // edges don't make it to the final list
187 (
188 pbm.mesh(),
189 isInletEdge,
191 );
192
193 // Now loop over the inlet faces once more and mark the common edges with
194 // processor boundaries
195 boolList isProcessorEdge(nEdges - nInternalEdges, false);
196 for (label edgeI = nInternalEdges; edgeI < nEdges; ++edgeI)
197 {
198 const edge& e = edges[edgeI];
199 const edge meshE = edge(mp[e[0]], mp[e[1]]);
200
201 if (isInletEdge.lookup(meshE, false))
202 {
203 isProcessorEdge[edgeI - nInternalEdges] = true;
204 }
205 }
207 return isProcessorEdge;
208}
209
210
211// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
212
213Foam::betaMax::betaMax
214(
215 const fvMesh& mesh,
216 const dictionary& dict
217)
218:
220 value_(dict.getOrDefault<scalar>("betaMax", Zero))
221{}
222
223
224// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
225
227(
228 const fvMesh& mesh,
229 const dictionary& dict
230)
231{
232 const word modelType(dict.getOrDefault<word>("betaMaxType", "value"));
233
234 auto* ctorPtr = dictionaryConstructorTable(modelType);
235
236 Info<< "betaMax type " << modelType << endl;
237
238 if (!ctorPtr)
239 {
241 (
242 dict,
243 "betaMaxType",
244 modelType,
245 *dictionaryConstructorTablePtr_
246 ) << exit(FatalIOError);
247 }
249 return autoPtr<betaMax>(ctorPtr(mesh, dict));
250}
251
252
253// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
254
255Foam::scalar Foam::betaMax::value() const
256{
257 return value_;
258}
259
260
261// ************************************************************************* //
const polyBoundaryMesh & pbm
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
const T & lookup(const Key &key, const T &deflt) const
Return hashed entry if it exists, or return the given default.
Definition HashTableI.H:222
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
label nEdges() const
Number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Base class for selecting the betaMax value, i.e. the value multiplying the Brinkman penalisation term...
Definition betaMax.H:50
const fvMesh & mesh_
Reference to mesh.
Definition betaMax.H:73
virtual scalar value() const
Get value.
Definition betaMax.C:248
scalar value_
betaMax value
Definition betaMax.H:78
static autoPtr< betaMax > New(const fvMesh &mesh, const dictionary &dict)
Construct and return the selected betaMax model.
Definition betaMax.C:220
boolList markProcessorEdges(const polyPatch &patch) const
Mark all common inlet - processor edges.
Definition betaMax.C:123
scalar computeLength(const dictionary &dict) const
Compute the characteristic length.
Definition betaMax.C:37
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label nNonProcessor() const
The number of patches before the first processor patch.
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
const polyMesh & mesh() const noexcept
Return the mesh reference.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition polyMesh.C:861
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition polyMesh.H:617
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition polyMesh.C:850
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static void syncEdgeMap(const polyMesh &mesh, EdgeMap< T > &edgeValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected edges.
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
const pointField & points
Namespace for bounding specifications. At the moment, mostly for tables.
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList area
Standard area field types (scalar, vector, tensor, etc).
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
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 POSIX.C:801
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
List< bool > boolList
A List of bools.
Definition List.H:60
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299