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-------------------------------------------------------------------------------
10License
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
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
40(
41 const fvPatch& p,
43)
44:
45 mixedFvPatchField<scalar>(p, iF),
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:
63 mixedFvPatchField<scalar>
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:
82 mixedFvPatchField<scalar>(p, iF),
83 baffle_(),
84 dict_
85 (
86 // Copy dictionary, but without "heavy" data chunks
87 dictionaryContent::copyDict
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A wrapper for dictionary content, without operators that could affect inheritance patterns.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:206
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:408
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:392
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual Field< scalar > & refGrad()
virtual Field< scalar > & refValue()
virtual scalarField & valueFraction()
const Type & lookupObject(const word &name, const bool recursive=false) const
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimDensity
dictionary dict
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))