KirchhoffShell.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-2023 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 
28 #include "KirchhoffShell.H"
30 #include "subCycle.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(KirchhoffShell, 0);
42 addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary);
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 bool KirchhoffShell::init(const dictionary& dict)
47 {
48  this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
49  return true;
50 }
51 
52 
53 void KirchhoffShell::solveDisplacement()
54 {
56 
57  const Time& time = primaryMesh().time();
58 
59  // Create operand fields for solid physics
60  const areaScalarField solidMass(rho()*h_);
61  const areaScalarField solidD(D()/solidMass);
62 
63  // Deep copy old times of shell displacement field
64  auto tw0 = tmp<areaScalarField>::New(w_.oldTime());
66 
67  // Flag terms to avoid redundant time-derivative calculations
68  const bool f0_enabled = (f0_.value() != scalar(0));
69  const bool f1_enabled = (f1_.value() != scalar(0));
70  const bool f2_enabled = (f2_.value() != scalar(0));
71 
72  // Restore various old fields in sub-cycling, if need be
73  if (nSubCycles_ > 1)
74  {
75  w_.oldTime() = w0_;
76  w_.oldTime().oldTime() = w00_;
77 
78  if (f0_enabled) laplaceW_.oldTime() = laplaceW0_;
79  if (f2_enabled) laplace2W_.oldTime() = laplace2W0_;
80  }
81 
82  for
83  (
84  subCycleTime wSubCycle
85  (
86  const_cast<Time&>(time),
87  nSubCycles_
88  );
89  !(++wSubCycle).end();
90  )
91  {
92  laplaceW_ = fac::laplacian(w_);
93  laplace2W_ = fac::laplacian(laplaceW_);
94 
95  faScalarMatrix wEqn
96  (
97  fam::d2dt2(w_)
98  + solidD*laplace2W_
99  ==
100  ps_/solidMass
101  + faOptions()(solidMass, w_, dimLength/sqr(dimTime))
102  );
103 
104  // Avoid time-derivative calculations for f terms, if possible
105  if (f0_enabled) wEqn -= f0_*sqrt(solidD)*fac::ddt(laplaceW_);
106  if (f1_enabled) wEqn += f1_*fam::ddt(w_);
107  if (f2_enabled) wEqn += f2_*solidD*fac::ddt(laplace2W_);
108 
109  faOptions().constrain(wEqn);
110 
111  wEqn.solve();
112 
113  // Cache various old fields inside the sub-cycling
114  if (wSubCycle.index() >= wSubCycle.nSubCycles())
115  {
116  w0_ = w_.oldTime();
117  w00_ = w_.oldTime().oldTime();
118 
119  if (f0_enabled) laplaceW0_ = laplaceW_.oldTime();
120  if (f2_enabled) laplace2W0_ = laplace2W_.oldTime();
121 
122  // Update shell acceleration
123  a_ = fac::d2dt2(w_);
124  }
125  }
126 
127  Info<< w_.name() << " min/max = " << gMinMax(w_) << endl;
128 
129  // Steal the deep-copy of old times to restore the shell displacement
130  w_.oldTime() = tw0;
131  w_.oldTime().oldTime() = tw00;
132 
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
138 
140 (
141  const word& modelType,
142  const fvMesh& mesh,
143  const dictionary& dict
144 )
145 :
146  vibrationShellModel(modelType, mesh, dict),
147  ps_
148  (
149  IOobject
150  (
151  "ps_" + regionName_,
152  regionMesh().time().timeName(),
153  regionMesh().thisDb(),
154  IOobject::READ_IF_PRESENT,
155  IOobject::NO_WRITE
156  ),
157  regionMesh(),
159  ),
160  h_
161  (
162  IOobject
163  (
164  "h_" + regionName_,
165  regionMesh().time().timeName(),
166  regionMesh().thisDb(),
167  IOobject::MUST_READ,
168  IOobject::AUTO_WRITE
169  ),
170  regionMesh()
171  ),
172  laplaceW_
173  (
174  IOobject
175  (
176  "laplaceW_" + regionName_,
177  regionMesh().time().timeName(),
178  regionMesh().thisDb(),
179  IOobject::NO_READ,
180  IOobject::NO_WRITE
181  ),
182  regionMesh(),
184  ),
185  laplace2W_
186  (
187  IOobject
188  (
189  "laplace2W_" + regionName_,
190  regionMesh().time().timeName(),
191  regionMesh().thisDb(),
192  IOobject::NO_READ,
193  IOobject::NO_WRITE
194  ),
195  regionMesh(),
197  ),
198  w0_
199  (
200  IOobject
201  (
202  "w0_" + regionName_,
203  regionMesh().time().timeName(),
204  regionMesh().thisDb(),
205  IOobject::NO_READ,
206  IOobject::NO_WRITE
207  ),
208  regionMesh(),
210  ),
211  w00_
212  (
213  IOobject
214  (
215  "w00_" + regionName_,
216  regionMesh().time().timeName(),
217  regionMesh().thisDb(),
218  IOobject::NO_READ,
219  IOobject::NO_WRITE
220  ),
221  regionMesh(),
223  ),
224  laplaceW0_
225  (
226  IOobject
227  (
228  "laplaceW0_" + regionName_,
229  regionMesh().time().timeName(),
230  regionMesh().thisDb(),
231  IOobject::NO_READ,
232  IOobject::NO_WRITE
233  ),
234  regionMesh(),
236  ),
237  laplace2W0_
238  (
239  IOobject
240  (
241  "laplace2W0_" + regionName_,
242  regionMesh().time().timeName(),
243  regionMesh().thisDb(),
244  IOobject::NO_READ,
245  IOobject::NO_WRITE
246  ),
247  regionMesh(),
249  ),
250  f0_("f0", dimless, dict),
251  f1_("f1", inv(dimTime), dict),
252  f2_("f2", dimTime, dict),
253  nNonOrthCorr_(1),
254  nSubCycles_(1)
255 {
256  init(dict);
257 }
259 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
260 
262 {}
263 
264 
266 {
267  nNonOrthCorr_ = solution().getLabel("nNonOrthCorr");
268  nSubCycles_ = solution().getLabel("nSubCycles");
269 
270  for (int nonOrth=0; nonOrth <= nNonOrthCorr_; ++nonOrth)
271  {
272  solveDisplacement();
273  }
274 }
275 
276 
278 {
279  const dimensionedScalar E("E", dimForce/dimArea , solid().E());
280  const dimensionedScalar nu("nu", dimless, solid().nu());
281 
282  return tmp<areaScalarField>(E*pow3(h_)/(12*(1 - sqr(nu))));
283 }
284 
285 
287 {
288  return areaScalarField::New
289  (
290  "rhos",
292  regionMesh(),
293  dimensionedScalar("rho", dimDensity, solid().rho()),
295  );
296 }
297 
299 {}
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 } // End namespace regionModels
304 } // End namespace Foam
305 
306 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
const dictionary & solution() const
Return the solution dictionary.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const tmp< areaScalarField > D() const
Return stiffness.
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
Foam::fa::options & faOptions() noexcept
Return faOptions.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
areaScalarField a_
Shell acceleration.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
const solidProperties & solid() const noexcept
Return solid properties.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
const dimensionSet dimless
Dimensionless.
areaScalarField w_
Shell displacement.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
label getLabel(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get< label >(const word&, keyType::option)
Definition: dictionary.H:1677
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: faPatchField.H:191
Macros for easy insertion into run-time selection tables.
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
tmp< GeometricField< Type, faPatchField, areaMesh > > d2dt2(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facD2dt2.C:39
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
tmp< faMatrix< Type > > d2dt2(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famD2dt2.C:40
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
const faMesh & regionMesh() const
Return the region mesh database.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
const iterator & end()
End of list for forward iterators.
Definition: UILList.H:489
const dimensionSet dimPressure
virtual void info()
Provide some feedback.
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:41
const dimensionSet dimForce
const dimensionSet dimDensity
dimensionedScalar pow3(const dimensionedScalar &ds)
addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Template specialisation for scalar faMatrix.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
const Time & time() const noexcept
Return the reference to the time database.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
virtual void evolveRegion()
Evolve the thermal baffle.
A class for managing sub-cycling times.
Definition: subCycleTime.H:46
volScalarField & nu
defineTypeNameAndDebug(KirchhoffShell, 0)
Do not request registration (bool: false)
KirchhoffShell(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components and dict.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127