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-2022 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 
29 #include "dictionaryContent.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const fvPatch& p,
43 )
44 :
46  baffle_(nullptr),
47  dict_()
48 {
49  refValue() = Zero;
50  refGrad() = Zero;
51  valueFraction() = 1;
52 }
53 
54 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
63  mixedFvPatchField<scalar>
64  (
65  ptf,
66  p,
67  iF,
68  mapper
69  ),
70  baffle_(nullptr),
71  dict_(ptf.dict_)
72 {}
73 
74 
76 (
77  const fvPatch& p,
79  const dictionary& dict
80 )
81 :
82  mixedFvPatchField<scalar>(p, iF),
83  baffle_(nullptr),
84  dict_
85  (
86  // Copy dictionary, but without "heavy" data chunks
87  dictionaryContent::copyDict
88  (
89  dict,
90  wordList(), // allow
91  wordList // deny
92  ({
93  "type", // redundant
94  "value", "refValue", "refGradient", "valueFraction"
95  })
96  )
97  )
98 {
99  this->readValueEntry(dict, IOobjectOption::MUST_READ);
100 
101  if (this->readMixedEntries(dict))
102  {
103  // Full restart
104  }
105  else
106  {
107  // Start from user entered data. Assume fixedValue.
108  refValue() = *this;
109  refGrad() = Zero;
110  valueFraction() = 1;
111  }
112 
113  if (!baffle_)
114  {
115  baffle_.reset(baffleType::New(p.boundaryMesh().mesh(), dict_));
116  }
117 }
118 
119 
121 (
122  const vibrationShellFvPatchScalarField& ptf,
123  const DimensionedField<scalar, volMesh>& iF
124 )
125 :
126  mixedFvPatchField<scalar>(ptf, iF),
127  baffle_(nullptr),
128  dict_(ptf.dict_)
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 {
136  if (this->updated())
137  {
138  return;
139  }
140 
141  const auto& transportProperties =
142  db().lookupObject<IOdictionary>("transportProperties");
143 
145 
146  baffle_->evolve();
147 
148  // rho * acceleration
149 
150  refGrad() = Zero; // safety (for any unmapped values)
151 
152  baffle_->vsm().mapToVolumePatch(baffle_->a(), refGrad(), patch().index());
153 
154  refGrad() *= rho.value();
155 
157  valueFraction() = Zero;
158 
160 }
161 
162 
164 {
166  dict_.write(os, false);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
173 (
175  vibrationShellFvPatchScalarField
176 );
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 } // End namespace Foam
181 
182 
183 // ************************************************************************* //
dictionary dict
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const objectRegistry & db() const
The associated objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Macros for easy insertion into run-time selection tables.
virtual Field< scalar > & refValue()
vibrationShellFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
bool updated() const noexcept
True if the boundary condition has already been updated.
Definition: fvPatchField.H:296
fvPatchField< scalar > fvPatchScalarField
virtual scalarField & valueFraction()
A FieldMapper for finite-volume patch fields.
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A wrapper for dictionary content, without operators that could affect inheritance patterns...
const Mesh & mesh() const noexcept
Return mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimDensity
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:204
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
virtual void write(Ostream &) const
Write.
volScalarField & p
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual Field< scalar > & refGrad()
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127