wallCellWallFunctionFvPatchScalarField.H
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) 2024 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 Class
25  Foam::wallCellWallFunctionFvPatchScalarField
26 
27 Description
28  Base class for wall functions that modify cell values
29 
30 SourceFiles
31  wallCellWallFunctionFvPatchScalarField.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef wallCellWallFunctionFvPatchScalarField_H
36 #define wallCellWallFunctionFvPatchScalarField_H
37 
38 #include "fixedValueFvPatchField.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 /*---------------------------------------------------------------------------*\
46  Class wallCellWallFunctionFvPatchScalarField Declaration
47 \*---------------------------------------------------------------------------*/
48 
50 :
51  public fixedValueFvPatchField<scalar>
52 {
53  // Private Data
54 
55  //- Tolerance used in weighted calculations
56  static scalar tol_;
57 
58  //- Master patch index
59  label masterPatchi_;
60 
61  //- Wall cell indices. Master patch only.
62  autoPtr<labelList> wallCellsPtr_;
63 
64  //- Wall cell fractions; i.e., what proportion of the wall cell value
65  // gets overridden. Master patch only.
66  autoPtr<scalarField> wallCellFractionPtr_;
67 
68 
69 protected:
70 
71  // Protected Member Functions
72 
73  // Access
74 
75  //- Return the master patch index
76  inline label masterPatchIndex() const
77  {
78  return masterPatchi_;
79  }
80 
81  //- Return the wall cell indices
82  inline const labelList& wallCells() const
83  {
84  return wallCellsPtr_();
85  }
86 
87  //- Return the wall cell fractions
88  inline const scalarField& wallCellFraction() const
89  {
90  return wallCellFractionPtr_();
91  }
92 
93 
94  //- Initialise the master cell indices and fractions, and allocate the
95  // turbulence fields
96  void initMaster();
97 
98  //- Sum values from a patch field into a cell field
99  static void patchFieldAddToCellField
100  (
101  const fvPatch& fvp,
102  const scalarField& pf,
103  scalarField& vf
104  );
105 
106  //- Average a set of patch fields into a wall cell field
108  (
109  const PtrList<scalarField>& pfs
110  ) const;
111 
112 
113 public:
114 
115  // Constructors
116 
117  //- Construct from patch, internal field and dictionary
119  (
120  const fvPatch&,
122  const dictionary&
123  );
124 
125  //- Construct by mapping given wallCellWallFunctionFvPatchScalarField
126  // onto a new patch
128  (
130  const fvPatch&,
132  const fieldMapper&
133  );
134 
135  //- Disallow copy without setting internal field reference
137  (
139  ) = delete;
140 
141  //- Copy constructor setting internal field reference
143  (
146  );
147 
148  //- Construct and return a clone setting internal field reference
150  (
152  ) const
153  {
155  (
157  );
158  }
159 
160 
161  //- Destructor
163  {}
164 
165 
166  // Member Functions
167 
168  // Mapping functions
169 
170  //- Map the given fvPatchField onto this fvPatchField
171  virtual void map(const fvPatchScalarField&, const fieldMapper&);
172 
173  //- Reset the fvPatchField to the given fvPatchField
174  // Used for mesh to mesh mapping
175  virtual void reset(const fvPatchScalarField&);
176 
177 
178  // Member Operators
179 
180  //- Inherit assignment operators to prevent warnings
182 };
183 
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 } // End namespace Foam
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 #endif
192 
193 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
tmp< fvPatchField< Type > > clone() const
Disallow clone without setting internal field reference.
Definition: fvPatchField.H:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for managing temporary objects.
Definition: tmp.H:55
Base class for wall functions that modify cell values.
const labelList & wallCells() const
Return the wall cell indices.
tmp< scalarField > patchFieldsToWallCellField(const PtrList< scalarField > &pfs) const
Average a set of patch fields into a wall cell field.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static void patchFieldAddToCellField(const fvPatch &fvp, const scalarField &pf, scalarField &vf)
Sum values from a patch field into a cell field.
wallCellWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
label masterPatchIndex() const
Return the master patch index.
void initMaster()
Initialise the master cell indices and fractions, and allocate the.
const scalarField & wallCellFraction() const
Return the wall cell fractions.
Namespace for OpenFOAM.
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