patchMeanVelocityForce.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) 2015-2016 OpenFOAM Foundation
9  Copyright (C) 2022 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 "patchMeanVelocityForce.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39  defineTypeNameAndDebug(patchMeanVelocityForce, 0);
40  addToRunTimeSelectionTable(option, patchMeanVelocityForce, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const word& sourceName,
50  const word& modelType,
51  const dictionary& dict,
52  const fvMesh& mesh
53 )
54 :
55  meanVelocityForce(sourceName, modelType, dict, mesh),
56  patch_(coeffs_.get<word>("patch")),
57  patchi_(mesh.boundaryMesh().findPatchID(patch_))
58 {
59  if (patchi_ < 0)
60  {
62  << "Cannot find patch " << patch_
63  << exit(FatalIOError);
64  }
65 }
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 (
72  const volVectorField& U
73 ) const
74 {
75  FixedList<scalar, 2> sumAmagUsumA(Zero);
76 
77  sumAmagUsumA[0] +=
78  sum
79  (
80  (flowDir_ & U.boundaryField()[patchi_])
81  * mesh_.boundary()[patchi_].magSf()
82  );
83  sumAmagUsumA[1] += sum(mesh_.boundary()[patchi_].magSf());
84 
85 
86  // If the mean velocity force is applied to a cyclic patch
87  // for parallel runs include contributions from processorCyclic patches
88  // generated from the decomposition of the cyclic patch
89  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
90 
91  if (Pstream::parRun() && isA<cyclicPolyPatch>(patches[patchi_]))
92  {
93  for
94  (
95  const label patchi
97  )
98  {
99  sumAmagUsumA[0] +=
100  sum
101  (
102  (flowDir_ & U.boundaryField()[patchi])
103  * mesh_.boundary()[patchi].magSf()
104  );
105  sumAmagUsumA[1] += sum(mesh_.boundary()[patchi].magSf());
106  }
107  }
108 
109  mesh_.reduce(sumAmagUsumA, sumOp<scalar>());
110 
111  return
112  (
113  sumAmagUsumA[0]
114  / stabilise(sumAmagUsumA[1], VSMALL)
115  );
116 }
117 
118 
119 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
patchMeanVelocityForce(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
dictionary dict
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
word patch_
Name of operand patch.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity averaged over the specified patch...
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
U
Definition: pEqn.H:72
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const polyBoundaryMesh & patches
static labelList patchIDs(const word &cyclicPolyPatchName, const polyBoundaryMesh &bm)
Return the indices of a processorCyclicPolyPatchs.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127