uniformFixedValueFvPatchField.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  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
41  uniformValue_(nullptr)
42 {}
43 
44 
45 template<class Type>
47 (
48  const fvPatch& p,
50  const Field<Type>& fld
51 )
52 :
54  uniformValue_(nullptr)
55 {}
56 
57 
58 template<class Type>
60 (
61  const fvPatch& p,
63  const dictionary& dict
64 )
65 :
66  fixedValueFvPatchField<Type>(p, iF, dict, false),
67  uniformValue_(PatchFunction1<Type>::New(p.patch(), "uniformValue", dict))
68 {
69  if (dict.found("value"))
70  {
72  (
73  Field<Type>("value", dict, p.size())
74  );
75  }
76  else
77  {
78  this->evaluate();
79  }
80 }
81 
82 
83 template<class Type>
85 (
87  const fvPatch& p,
89  const fvPatchFieldMapper& mapper
90 )
91 :
92  fixedValueFvPatchField<Type>(p, iF), // Don't map
93  uniformValue_(ptf.uniformValue_.clone(p.patch()))
94 {
95  if (mapper.direct() && !mapper.hasUnmapped())
96  {
97  // Use mapping instead of re-evaluation
98  this->map(ptf, mapper);
99  }
100  else
101  {
102  // Evaluate since value not mapped
103  this->evaluate();
104  }
105 }
106 
107 
108 template<class Type>
110 (
112 )
113 :
115  uniformValue_(ptf.uniformValue_.clone(this->patch().patch()))
116 {}
117 
118 
119 template<class Type>
121 (
124 )
125 :
127  uniformValue_(ptf.uniformValue_.clone(this->patch().patch()))
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 template<class Type>
135 (
136  const fvPatchFieldMapper& mapper
137 )
138 {
140  uniformValue_().autoMap(mapper);
141 
142  if (uniformValue_().constant())
143  {
144  // If mapper is not dependent on time we're ok to evaluate
145  this->evaluate();
146  }
147 }
148 
149 
150 template<class Type>
152 (
153  const fvPatchField<Type>& ptf,
154  const labelList& addr
155 )
156 {
158 
159  const uniformFixedValueFvPatchField& tiptf =
160  refCast<const uniformFixedValueFvPatchField>(ptf);
161 
162  uniformValue_().rmap(tiptf.uniformValue_(), addr);
163 }
164 
165 
166 template<class Type>
168 {
169  if (this->updated())
170  {
171  return;
172  }
173 
174  const scalar t = this->db().time().timeOutputValue();
175  fvPatchField<Type>::operator==(uniformValue_->value(t));
177 }
178 
179 
180 template<class Type>
182 {
184  uniformValue_->writeData(os);
185  this->writeEntry("value", os);
186 }
187 
188 
189 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::uniformFixedValueFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: uniformFixedValueFvPatchField.C:135
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::uniformFixedValueFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: uniformFixedValueFvPatchField.C:152
Foam::uniformFixedValueFvPatchField::uniformFixedValueFvPatchField
uniformFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: uniformFixedValueFvPatchField.C:35
uniformFixedValueFvPatchField.H
Foam::FieldMapper::direct
virtual bool direct() const =0
Foam::FieldMapper::hasUnmapped
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:80
Foam::uniformFixedValueFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: uniformFixedValueFvPatchField.C:167
Foam::uniformFixedValueFvPatchField::write
virtual void write(Ostream &os) const
Write.
Definition: uniformFixedValueFvPatchField.C:181
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:60
Foam::uniformFixedValueFvPatchField
This boundary condition provides a uniform fixed value condition.
Definition: uniformFixedValueFvPatchField.H:90
Foam::List< label >
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
constant
constant condensation/saturation model.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::stringOps::evaluate
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Definition: stringOpsEvaluate.C:37