Loading...
Searching...
No Matches
uniformJumpFvPatchField.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
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
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const fvPatch& p,
38)
40 fixedJumpFvPatchField<Type>(p, iF),
42{}
43
44
45template<class Type>
47(
49 const fvPatch& p,
51 const fvPatchFieldMapper& mapper
52)
54 fixedJumpFvPatchField<Type>(ptf, p, iF, mapper),
56{}
57
58
59template<class Type>
61(
62 const fvPatch& p,
64 const dictionary& dict,
65 const bool needValue
66)
67:
68 fixedJumpFvPatchField<Type>(p, iF, dict, false), // needValue = false
69 jumpTable_()
70{
71 if (needValue)
72 {
73 if (this->cyclicPatch().owner())
74 {
75 jumpTable_ = Function1<Type>::New("jumpTable", dict, &this->db());
76 }
77
78 if (!this->readValueEntry(dict))
79 {
88(
90)
102)
103:
104 fixedJumpFvPatchField<Type>(ptf, iF),
106{}
107
108
109// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110
111template<class Type>
113{
114 if (this->updated())
115 {
116 return;
117 }
118
119 if (this->cyclicPatch().owner())
120 {
121 this->setJump(jumpTable_->value(this->db().time().value()));
123
125}
126
127
128template<class Type>
130{
132
133 if (this->cyclicPatch().owner())
134 {
135 jumpTable_->writeData(os);
136 }
137}
138
139
140// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
const cyclicFvPatch & cyclicPatch() const
Return local reference cast into the cyclic patch.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
virtual void setJump(const Field< Type > &jump)
Set the jump field.
fixedJumpFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const objectRegistry & db() const
The associated objectRegistry.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a jump condition, using the cyclic condition as a base....
uniformJumpFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients.
autoPtr< Function1< Type > > jumpTable_
The "jump" table.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
dictionary dict