patchMeanVelocityForce.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "patchMeanVelocityForce.H"
29 
30 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36  defineTypeNameAndDebug(patchMeanVelocityForce, 0);
37 
39  (
40  option,
41  patchMeanVelocityForce,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::fv::patchMeanVelocityForce::patchMeanVelocityForce
51 (
52  const word& sourceName,
53  const word& modelType,
54  const dictionary& dict,
55  const fvMesh& mesh
56 )
57 :
58  meanVelocityForce(sourceName, modelType, dict, mesh),
59  patch_(coeffs_.lookup("patch")),
60  patchi_(mesh.boundaryMesh().findPatchID(patch_))
61 {
62  if (patchi_ < 0)
63  {
65  << "Cannot find patch " << patch_
66  << exit(FatalError);
67  }
68 }
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
74 (
75  const volVectorField& U
76 ) const
77 {
78  vector2D sumAmagUsumA
79  (
80  sum
81  (
82  (flowDir_ & U.boundaryField()[patchi_])
83  *mesh_.boundary()[patchi_].magSf()
84  ),
85  sum(mesh_.boundary()[patchi_].magSf())
86  );
87 
88 
89  // If the mean velocity force is applied to a cyclic patch
90  // for parallel runs include contributions from processorCyclic patches
91  // generated from the decomposition of the cyclic patch
92  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
93 
94  if (Pstream::parRun() && isA<cyclicPolyPatch>(patches[patchi_]))
95  {
96  labelList processorCyclicPatches
97  (
98  processorCyclicPolyPatch::patchIDs(patch_, patches)
99  );
100 
101  forAll(processorCyclicPatches, pcpi)
102  {
103  const label patchi = processorCyclicPatches[pcpi];
104 
105  sumAmagUsumA.x() +=
106  sum
107  (
108  (flowDir_ & U.boundaryField()[patchi])
109  *mesh_.boundary()[patchi].magSf()
110  );
111 
112  sumAmagUsumA.y() += sum(mesh_.boundary()[patchi].magSf());
113  }
114  }
115 
116  mesh_.reduce(sumAmagUsumA, sumOp<vector2D>());
117 
118  return sumAmagUsumA.x()/sumAmagUsumA.y();
119 }
120 
121 
122 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Boundary & boundaryField() const
Return const-reference to the boundary field.
label findPatchID(const word &patchName) const
Find patch index given a name.
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
patches[0]
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity.
const Cmpt & x() const
Definition: Vector2DI.H:68
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const Cmpt & y() const
Definition: Vector2DI.H:74
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Foam::polyBoundaryMesh.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:397
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static labelList patchIDs(const word &cyclicPolyPatchName, const polyBoundaryMesh &bm)
Return the indices of a processorCyclicPolyPatchs.
Calculates and applies the force necessary to maintain the specified mean velocity.
Namespace for OpenFOAM.