Loading...
Searching...
No Matches
cellMotionFvPatchField.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "fvMesh.H"
30#include "volMesh.H"
31#include "pointFields.H"
32
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const fvPatch& p,
42:
43 fixedValueFvPatchField<Type>(p, iF)
44{}
45
46
47template<class Type>
49(
51 const fvPatch& p,
53 const fvPatchFieldMapper& mapper
55:
56 fixedValueFvPatchField<Type>(ptf, p, iF, mapper)
57{}
58
59
60template<class Type>
62(
63 const fvPatch& p,
65 const dictionary& dict
67:
68 fixedValueFvPatchField<Type>(p, iF, dict)
69{}
70
71
72template<class Type>
74(
77:
78 fixedValueFvPatchField<Type>(ptf)
79{}
80
81
82template<class Type>
84(
87)
88:
90{}
91
92
93// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94
95template<class Type>
97{
98 if (this->updated())
99 {
100 return;
101 }
102
103 const fvPatch& p = this->patch();
104 const polyPatch& pp = p.patch();
105 const fvMesh& mesh = this->internalField().mesh();
106 const pointField& points = mesh.points();
107
108 word pfName = this->internalField().name();
109 pfName.replace("cell", "point");
110
111 const GeometricField<Type, pointPatchField, pointMesh>& pointMotion =
112 this->db().objectRegistry::template
113 lookupObject<GeometricField<Type, pointPatchField, pointMesh>>
114 (pfName);
115
116 forAll(p, i)
117 {
118 this->operator[](i) = pp[i].average(points, pointMotion);
124
125template<class Type>
127{
130}
131
132// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Type & operator[](const label i)
Definition UListI.H:363
Foam::cellMotionFvPatchField.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
cellMotionFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
string & replace(const std::string &s1, const std::string &s2, size_type pos=0)
Replace first occurrence of sub-string s1 with s2, beginning at pos.
Definition string.C:101
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
vectorField pointField
pointField is a vectorField.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299