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-2020 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 )
37 :
39  uniformValue_()
40 {}
41 
42 
43 template<class Type>
46 (
47  const pointPatch& p,
49  const dictionary& dict
50 )
51 :
53  uniformValue_(Function1<Type>::New("uniformValue", dict))
54 {
55  if (dict.found("value"))
56  {
58  (
59  Field<Type>("value", dict, p.size())
60  );
61  }
62  else
63  {
64  const scalar t = this->db().time().timeOutputValue();
65  fixedValuePointPatchField<Type>::operator=(uniformValue_->value(t));
66  }
67 }
68 
69 
70 template<class Type>
73 (
75  const pointPatch& p,
77  const pointPatchFieldMapper& mapper
78 )
79 :
80  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
81  uniformValue_(ptf.uniformValue_, false)
82 {
83  // For safety re-evaluate
84  const scalar t = this->db().time().timeOutputValue();
85  fixedValuePointPatchField<Type>::operator=(uniformValue_->value(t));
86 }
87 
88 
89 template<class Type>
92 (
95 )
96 :
98  uniformValue_(ptf.uniformValue_, false)
99 {
100  // For safety re-evaluate
101  const scalar t = this->db().time().timeOutputValue();
102  fixedValuePointPatchField<Type>::operator==(uniformValue_->value(t));
103 }
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
108 template<class Type>
110 {
111  if (this->updated())
112  {
113  return;
114  }
115 
116  const scalar t = this->db().time().timeOutputValue();
117  fixedValuePointPatchField<Type>::operator==(uniformValue_->value(t));
118 
120 }
121 
122 
123 template<class Type>
125 write(Ostream& os) const
126 {
127  // Note: write value
129  writeEntry(os, uniformValue_());
130 }
131 
132 
133 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A FixedValue boundary condition for pointField.
Foam::pointPatchFieldMapper.
label size() const
Return size.
Pre-declare SubField and related Field type.
Definition: Field.H:56
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
Enables the specification of a uniform fixed value boundary condition.
uniformFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.