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-2021 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 #include "dictionaryContent.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const fvPatch& p,
43 )
44 :
46  baffle_(),
47  dict_()
48 {
49  refValue() = 0;
50  refGrad() = 0;
51  valueFraction() = 1;
52 }
53 
54 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
64  (
65  ptf,
66  p,
67  iF,
68  mapper
69  ),
70  baffle_(),
71  dict_(ptf.dict_)
72 {}
73 
74 
76 (
77  const fvPatch& p,
79  const dictionary& dict
80 )
81 :
83  baffle_(),
84  dict_
85  (
86  // Copy dictionary, but without "heavy" data chunks
88  (
89  dict,
90  wordRes(), // allow
91  wordRes // deny
92  ({
93  "type", // redundant
94  "value", "refValue", "refGradient", "valueFraction"
95  })
96  )
97  )
98 {
100 
101  if (dict.found("refValue"))
102  {
103  // Full restart
104  refValue() = scalarField("refValue", dict, p.size());
105  refGrad() = scalarField("refGradient", dict, p.size());
106  valueFraction() = scalarField("valueFraction", dict, p.size());
107  }
108  else
109  {
110  // Start from user entered data. Assume fixedValue.
111  refValue() = *this;
112  refGrad() = 0;
113  valueFraction() = 1;
114  }
115 
116  if (!baffle_)
117  {
118  baffle_.reset(regionModels::vibrationShellModel::New(p, dict).ptr());
119  }
120 }
121 
122 
124 (
127 )
128 :
129  mixedFvPatchField<scalar>(ptf, iF),
130  baffle_(),
131  dict_(ptf.dict_)
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  if (this->updated())
140  {
141  return;
142  }
143 
144  baffle_->evolve();
145 
146  const auto& transportProperties =
147  db().lookupObject<IOdictionary>("transportProperties");
148 
150 
151  const areaScalarField aRho(rho*baffle_->a());
152 
153  baffle_->vsm().mapToField(aRho, this->refGrad());
154 
155  this->refValue() = Zero;
156  this->valueFraction() = Zero;
157 
159 }
160 
161 
163 {
165  dict_.write(os, false);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
172 (
175 );
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace Foam
180 
181 
182 // ************************************************************************* //
Foam::fvPatchField< scalar >
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:250
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dictionaryContent::copyDict
static dictionary copyDict(const dictionary &input, const wordList &allow=wordList(), const wordList &deny=wordList())
Definition: dictionaryContent.C:78
Foam::mixedFvPatchField< scalar >::refValue
virtual Field< scalar > & refValue()
Definition: mixedFvPatchField.H:230
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: dictionaryI.H:87
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:240
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:387
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:162
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mixedFvPatchField.C:235
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)
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:321
Foam::vibrationShellFvPatchScalarField::vibrationShellFvPatchScalarField
vibrationShellFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: vibrationShellFvPatchScalarField.C:40
dictionaryContent.H
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:206
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::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
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:404
Foam::vibrationShellFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: vibrationShellFvPatchScalarField.C:137
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54