uniformFixedValuePointPatchField.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) 2011-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 \*---------------------------------------------------------------------------*/
25 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
29 
30 template<class Type>
33 (
34  const pointPatch& p,
36  const dictionary& dict
37 )
38 :
39  fixedValuePointPatchField<Type>(p, iF, dict, false),
40  uniformValue_
41  (
42  Function1<Type>::New
43  (
44  "uniformValue",
45  this->db().time().userUnits(),
46  iF.dimensions(),
47  dict
48  )
49  )
50 {
51  if (dict.found("value"))
52  {
54  (
55  Field<Type>("value", iF.dimensions(), dict, p.size())
56  );
57  }
58  else
59  {
60  const scalar t = this->db().time().value();
61  fixedValuePointPatchField<Type>::operator=(uniformValue_->value(t));
62  }
63 }
64 
65 
66 template<class Type>
69 (
71  const pointPatch& p,
73  const fieldMapper& mapper
74 )
75 :
76  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
77  uniformValue_(ptf.uniformValue_, false)
78 {
79  // For safety re-evaluate
80  const scalar t = this->db().time().value();
81  fixedValuePointPatchField<Type>::operator=(uniformValue_->value(t));
82 }
83 
84 
85 template<class Type>
88 (
91 )
92 :
93  fixedValuePointPatchField<Type>(ptf, iF),
94  uniformValue_(ptf.uniformValue_, false)
95 {
96  // For safety re-evaluate
97  const scalar t = this->db().time().value();
98  fixedValuePointPatchField<Type>::operator==(uniformValue_->value(t));
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class Type>
106 {
107  if (this->updated())
108  {
109  return;
110  }
111 
112  const scalar t = this->db().time().value();
113  fixedValuePointPatchField<Type>::operator==(uniformValue_->value(t));
114 
116 }
117 
118 
119 template<class Type>
121 write(Ostream& os) const
122 {
123  // Note: write value
125  writeEntry
126  (
127  os,
128  this->db().time().userUnits(),
129  this->internalField().dimensions(),
130  uniformValue_())
131  ;
132 }
133 
134 
135 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Run-time selectable general function of one variable.
Definition: Function1.H:125
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
A FixedValue boundary condition for pointField.
virtual void operator=(const Field< Type > &)
const Time & time() const
Return time.
const objectRegistry & db() const
Return local objectRegistry.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
Enables the specification of a uniform fixed value boundary condition.
uniformFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void operator==(const valuePointPatchField< Type > &)
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dictionary dict
volScalarField & p