uniformFixedGradientFvPatchField.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 | Copyright (C) 2013-2015 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>
32 (
33  const fvPatch& p,
35 )
36 :
38  uniformGradient_()
39 {}
40 
41 
42 template<class Type>
44 (
45  const fvPatch& p,
47  const Field<Type>& fld
48 )
49 :
51  uniformGradient_()
52 {}
53 
54 
55 template<class Type>
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedGradientFvPatchField<Type>(ptf, p, iF, mapper),
65  uniformGradient_(ptf.uniformGradient_().clone().ptr())
66 {
67  // For safety re-evaluate
68  const scalar t = this->db().time().timeOutputValue();
69  this->gradient() = uniformGradient_->value(t);
70 }
71 
72 
73 template<class Type>
75 (
76  const fvPatch& p,
78  const dictionary& dict
79 )
80 :
82  uniformGradient_(DataEntry<Type>::New("uniformGradient", dict))
83 {
84  if (dict.found("gradient"))
85  {
86  this->gradient() = Field<Type>("gradient", dict, p.size());
87  }
88  else
89  {
90  const scalar t = this->db().time().timeOutputValue();
91  this->gradient() = uniformGradient_->value(t);
92  }
93 
95 }
96 
97 
98 template<class Type>
100 (
102 )
103 :
105  uniformGradient_
106  (
107  ptf.uniformGradient_.valid()
108  ? ptf.uniformGradient_().clone().ptr()
109  : NULL
110  )
111 {}
112 
113 
114 template<class Type>
116 (
119 )
120 :
122  uniformGradient_
123  (
124  ptf.uniformGradient_.valid()
125  ? ptf.uniformGradient_().clone().ptr()
126  : NULL
127  )
128 {
129  // For safety re-evaluate
130  const scalar t = this->db().time().timeOutputValue();
131 
132  if (ptf.uniformGradient_.valid())
133  {
134  this->gradient() = uniformGradient_->value(t);
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
141 template<class Type>
143 {
144  if (this->updated())
145  {
146  return;
147  }
148 
149  const scalar t = this->db().time().timeOutputValue();
150  this->gradient() = uniformGradient_->value(t);
151 
153 }
154 
155 
156 template<class Type>
158 {
160  uniformGradient_->writeData(os);
161  this->writeEntry("value", os);
162 }
163 
164 
165 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::fvPatchFieldMapper.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
runTime write()
uniformFixedGradientFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
dictionary dict
This boundary condition supplies a fixed gradient condition, such that the patch values are calculate...
volScalarField & p
Definition: createFields.H:51
Pre-declare SubField and related Field type.
Definition: Field.H:57
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject( name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE ))
This boundary condition provides a uniform fixed gradient condition.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
virtual void updateCoeffs()
Update the coefficients associated with the patch field.