vibrationShellFvPatchScalarField.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) 2019-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
45  baffle_(),
46  dict_(dictionary::null)
47 {
48  refValue() = 0;
49  refGrad() = 0;
50  valueFraction() = 1;
51 }
52 
53 
55 (
57  const fvPatch& p,
59  const fvPatchFieldMapper& mapper
60 )
61 :
63  (
64  ptf,
65  p,
66  iF,
67  mapper
68  ),
69  baffle_(),
70  dict_(ptf.dict_)
71 {}
72 
73 
75 (
76  const fvPatch& p,
78  const dictionary& dict
79 )
80 :
82  baffle_(),
83  dict_(dict)
84 {
86 
87  if (dict.found("refValue"))
88  {
89  // Full restart
90  refValue() = scalarField("refValue", dict, p.size());
91  refGrad() = scalarField("refGradient", dict, p.size());
92  valueFraction() = scalarField("valueFraction", dict, p.size());
93  }
94  else
95  {
96  // Start from user entered data. Assume fixedValue.
97  refValue() = *this;
98  refGrad() = 0;
99  valueFraction() = 1;
100  }
101 
102  if (!baffle_)
103  {
104  baffle_.reset(regionModels::vibrationShellModel::New(p, dict).ptr());
105  }
106 }
107 
108 
110 (
113 )
114 :
115  mixedFvPatchField<scalar>(ptf, iF),
116  baffle_(),
117  dict_(ptf.dict_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 {
125  if (this->updated())
126  {
127  return;
128  }
129 
130  baffle_->evolve();
131 
132  const auto& transportProperties =
133  db().lookupObject<IOdictionary>("transportProperties");
134 
136 
137  const areaScalarField aRho(rho*baffle_->a());
138 
139  baffle_->vsm().mapToField(aRho, this->refGrad());
140 
141  this->refValue() = Zero;
142  this->valueFraction() = Zero;
143 
145 }
146 
147 
149 {
151 
152  dict_.write(os, false);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
159 (
162 );
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace Foam
167 
168 
169 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::mixedFvPatchField< scalar >::valueFraction
virtual scalarField & valueFraction()
Definition: mixedFvPatchField.H:247
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::mixedFvPatchField< scalar >::refValue
virtual Field< scalar > & refValue()
Definition: mixedFvPatchField.H:227
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimDensity
const dimensionSet dimDensity
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:364
Foam::regionModels::vibrationShellModel::New
static autoPtr< vibrationShellModel > New(const fvPatch &patch, const dictionary &dict)
Return a reference to the selected model using dictionary.
Definition: vibrationShellModelNew.C:40
transportProperties
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::mixedFvPatchField< scalar >::refGrad
virtual Field< scalar > & refGrad()
Definition: mixedFvPatchField.H:237
rho
rho
Definition: readInitialConditions.H:88
Foam::fvPatchField< scalar >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:373
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:385
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::vibrationShellFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: vibrationShellFvPatchScalarField.C:148
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mixedFvPatchField.C:229
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::mixedFvPatchField< scalar >
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:313
Foam::vibrationShellFvPatchScalarField::vibrationShellFvPatchScalarField
vibrationShellFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: vibrationShellFvPatchScalarField.C:39
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
Foam::vibrationShellFvPatchScalarField
Definition: vibrationShellFvPatchScalarField.H:89
Foam::dictionary::write
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:206
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
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
vibrationShellFvPatchScalarField.H
Foam::GeometricField< scalar, faPatchField, areaMesh >
Foam::fvPatchField< scalar >::operator=
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:379
Foam::vibrationShellFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: vibrationShellFvPatchScalarField.C:123
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54