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 "fvPatchFields.H"
32 #include "subCycle.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(KirchhoffShell, 0);
44 
45 addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary);
46 
47 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 bool KirchhoffShell::init(const dictionary& dict)
50 {
51  this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
52  return true;
53 }
54 
55 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
56 
58 {
59  if (debug)
60  {
62  }
63 
64  const Time& time = primaryMesh().time();
65 
66  areaScalarField solidMass(rho()*h_);
67  areaScalarField solidD(D()/solidMass);
68 
69  // Save old times
72 
73  if (nSubCycles_ > 1)
74  {
75  // Restore the oldTime in sub-cycling
76  w_.oldTime() = w0_;
77  w_.oldTime().oldTime() = w00_;
80  }
81 
82  for
83  (
84  subCycleTime wSubCycle
85  (
86  const_cast<Time&>(time),
88  );
89  !(++wSubCycle).end();
90  )
91  {
92 
95 
96  faScalarMatrix wEqn
97  (
98  fam::d2dt2(w_)
99  + f1_*fam::ddt(w_)
100  - f0_*sqrt(solidD)*fac::ddt(laplaceW_)
101  + solidD*(laplace2W_ + f2_*fac::ddt(laplace2W_))
102  ==
103  ps_/solidMass
104  + faOptions()(solidMass, w_, dimLength/sqr(dimTime))
105  );
106 
107  faOptions().constrain(wEqn);
108 
109  wEqn.solve();
110 
111  if (wSubCycle.index() >= wSubCycle.nSubCycles())
112  {
113  // Cache oldTimes inside the sub-cycling
114  w0_ = w_.oldTime();
115  w00_ = w_.oldTime().oldTime();
118 
119  // Update shell acceleration
120  a_ = fac::d2dt2(w_);
121  }
122  }
123 
124  Info<< w_.name() << " min/max = " << gMinMax(w_) << endl;
125 
126  // Restore old time in main time
127  w_.oldTime() = w0;
128  w_.oldTime().oldTime() = w00;
129 
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135 
137 (
138  const word& modelType,
139  const fvMesh& mesh,
140  const dictionary& dict
141 )
142 :
143  vibrationShellModel(modelType, mesh, dict),
144  f0_("f0", dimless, dict),
145  f1_("f1", inv(dimTime), dict),
146  f2_("f2", dimTime, dict),
147  nNonOrthCorr_(1),
148  nSubCycles_(1),
149  ps_
150  (
151  IOobject
152  (
153  "ps_" + regionName_,
154  regionMesh().time().timeName(),
155  regionMesh().thisDb(),
156  IOobject::READ_IF_PRESENT,
157  IOobject::NO_WRITE
158  ),
159  regionMesh(),
161  ),
162  h_
163  (
164  IOobject
165  (
166  "h_" + regionName_,
167  regionMesh().time().timeName(),
168  regionMesh().thisDb(),
169  IOobject::MUST_READ,
170  IOobject::AUTO_WRITE
171  ),
172  regionMesh()
173  ),
174  laplaceW_
175  (
176  IOobject
177  (
178  "laplaceW_" + regionName_,
179  regionMesh().time().timeName(),
180  regionMesh().thisDb(),
181  IOobject::NO_READ,
182  IOobject::NO_WRITE
183  ),
184  regionMesh(),
186  ),
187  laplace2W_
188  (
189  IOobject
190  (
191  "laplace2W_" + regionName_,
192  regionMesh().time().timeName(),
193  regionMesh().thisDb(),
194  IOobject::NO_READ,
195  IOobject::NO_WRITE
196  ),
197  regionMesh(),
199  ),
200  w0_
201  (
202  IOobject
203  (
204  "w0_" + regionName_,
205  regionMesh().time().timeName(),
206  regionMesh().thisDb(),
207  IOobject::NO_READ,
208  IOobject::NO_WRITE
209  ),
210  regionMesh(),
212  ),
213  w00_
214  (
215  IOobject
216  (
217  "w00_" + regionName_,
218  regionMesh().time().timeName(),
219  regionMesh().thisDb(),
220  IOobject::NO_READ,
221  IOobject::NO_WRITE
222  ),
223  regionMesh(),
225  ),
226  laplaceW0_
227  (
228  IOobject
229  (
230  "laplaceW0_" + regionName_,
231  regionMesh().time().timeName(),
232  regionMesh().thisDb(),
233  IOobject::NO_READ,
234  IOobject::NO_WRITE
235  ),
236  regionMesh(),
238  ),
239  laplace2W0_
240  (
241  IOobject
242  (
243  "laplace2W0_" + regionName_,
244  regionMesh().time().timeName(),
245  regionMesh().thisDb(),
246  IOobject::NO_READ,
247  IOobject::NO_WRITE
248  ),
249  regionMesh(),
251  )
252 {
253  init(dict);
254 }
256 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
257 
259 {}
260 
261 
263 {
264  nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
265  nSubCycles_ = solution().get<label>("nSubCycles");
266 
267  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
268  {
270  }
271 }
272 
273 
275 {
276  const dimensionedScalar E("E", dimForce/dimArea , solid().E());
277  const dimensionedScalar nu("nu", dimless, solid().nu());
278 
279  return tmp<areaScalarField>(E*pow3(h_)/(12*(1 - sqr(nu))));
280 }
281 
282 
284 {
285  return tmp<areaScalarField>
286  (
287  new areaScalarField
288  (
289  IOobject
290  (
291  "rhos",
292  regionMesh().time().timeName(),
293  regionMesh().thisDb(),
297  ),
299  dimensionedScalar("rho", dimDensity, solid().rho()),
301  )
302  );
303 }
304 
306 {}
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 } // End namespace regionModels
311 } // End namespace Foam
312 
313 // ************************************************************************* //
dictionary dict
const dictionary & solution() const
Return the solution dictionary.
areaScalarField h_
Thickness [m].
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
areaScalarField w00_
Cache w.oldTime.oldTime() in sub-cycling.
const tmp< areaScalarField > D() const
Return stiffness.
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
areaScalarField laplace2W_
Laplace of the Laplace for the displacement.
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
areaScalarField laplace2W0_
Cache laplace2.oldTime() in sub-cycling.
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
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
void solveDisplacement()
Solve energy equation.
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
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.
areaScalarField ps_
External surface source [Pa].
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
#define w0
Definition: blockCreate.C:26
const iterator & end()
End of list for forward iterators.
Definition: UILList.H:489
const dimensionSet dimPressure
virtual void info()
Provide some feedback.
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:41
const dimensionSet dimForce
int debug
Static debugging option.
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.
areaScalarField laplaceW_
Laplace of the displacement.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
Template specialisation for scalar faMatrix.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
areaScalarField w0_
Cache w.oldTime() in sub-cycling.
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:172
areaScalarField laplaceW0_
Cache laplaceW.oldTime() in sub-cycling.
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.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:80
label nNonOrthCorr_
Number of non orthogonal correctors.
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127