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-2020 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 
44 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 {
51  this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
52  return true;
53 }
54 
55 
57 {
58  if (debug)
59  {
61  }
62 
63  const Time& time = primaryMesh().time();
64 
65  areaScalarField solidMass(rho()*h_);
66  areaScalarField solidD(D()/solidMass);
67 
68  // Save old times
71 
72  if (nSubCycles_ > 1)
73  {
74  // Restore the oldTime in sub-cycling
75  w_.oldTime() = w0_;
76  w_.oldTime().oldTime() = w00_;
79  }
80 
81  for
82  (
83  subCycleTime wSubCycle
84  (
85  const_cast<Time&>(time),
87  );
88  !(++wSubCycle).end();
89  )
90  {
91 
94 
95  faScalarMatrix wEqn
96  (
97  fam::d2dt2(w_)
98  + f1_*fam::ddt(w_)
99  - f0_*sqrt(solidD)*fac::ddt(laplaceW_)
100  + solidD*(laplace2W_ + f2_*fac::ddt(laplace2W_))
101  ==
102  ps_/solidMass
103  + faOptions()(solidMass, w_, dimLength/sqr(dimTime))
104  );
105 
106  faOptions().constrain(wEqn);
107 
108  wEqn.solve();
109 
110  Info<< "w min/max = " << min(w_) << ", " << max(w_) << endl;
111 
112  if (wSubCycle.index() >= wSubCycle.nSubCycles())
113  {
114  // Cache oldTimes inside the sub-cycling
115  w0_ = w_.oldTime();
116  w00_ = w_.oldTime().oldTime();
119 
120  // Update shell acceleration
121  a_ = fac::d2dt2(w_);
122  }
123  }
124 
125  // Restore old time in main time
126  w_.oldTime() = w0;
127  w_.oldTime().oldTime() = w00;
128 
129  faOptions().correct(w_);
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
136 (
137  const word& modelType,
138  const fvPatch& patch,
139  const dictionary& dict
140 )
141 :
142  vibrationShellModel(modelType, patch, dict),
143  f0_("f0", dimless, dict),
144  f1_("f1", inv(dimTime), dict),
145  f2_("f2", dimTime, dict),
146  nNonOrthCorr_(1),
147  nSubCycles_(1),
148  ps_
149  (
150  IOobject
151  (
152  "ps_" + regionName_,
153  primaryMesh().time().timeName(),
154  primaryMesh(),
157  ),
158  regionMesh(),
160  ),
161  h_
162  (
163  IOobject
164  (
165  "h_" + regionName_,
166  primaryMesh().time().timeName(),
167  primaryMesh(),
170  ),
171  regionMesh()
172  ),
173  laplaceW_
174  (
175  IOobject
176  (
177  "laplaceW_" + regionName_,
178  primaryMesh().time().timeName(),
179  primaryMesh(),
182  ),
183  regionMesh(),
185  ),
186  laplace2W_
187  (
188  IOobject
189  (
190  "laplace2W_" + regionName_,
191  primaryMesh().time().timeName(),
192  primaryMesh(),
195  ),
196  regionMesh(),
198  ),
199  w0_
200  (
201  IOobject
202  (
203  "w0_" + regionName_,
204  primaryMesh().time().timeName(),
205  primaryMesh(),
208  ),
209  regionMesh(),
211  ),
212  w00_
213  (
214  IOobject
215  (
216  "w00_" + regionName_,
217  primaryMesh().time().timeName(),
218  primaryMesh(),
221  ),
222  regionMesh(),
224  ),
225  laplaceW0_
226  (
227  IOobject
228  (
229  "laplaceW0_" + regionName_,
230  primaryMesh().time().timeName(),
231  primaryMesh(),
234  ),
235  regionMesh(),
237  ),
238  laplace2W0_
239  (
240  IOobject
241  (
242  "laplace2W0_" + regionName_,
243  primaryMesh().time().timeName(),
244  primaryMesh(),
247  ),
248  regionMesh(),
250  )
251 {
252  init();
253 }
254 
255 
256 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
257 
258 void KirchhoffShell::init()
259 {}
260 
261 
263 {}
264 
265 
267 {
268  nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
269  nSubCycles_ = solution().get<label>("nSubCycles");
270 
271  for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
272  {
274  }
275 }
276 
277 
279 {
280  const dimensionedScalar E("E", dimForce/dimArea , solid().E());
281  const dimensionedScalar nu("nu", dimless, solid().nu());
282 
283  return tmp<areaScalarField>(E*pow3(h_)/(12*(1 - sqr(nu))));
284 }
285 
286 
288 {
289  return tmp<areaScalarField>
290  (
291  new areaScalarField
292  (
293  IOobject
294  (
295  "rhos",
296  primaryMesh().time().timeName(),
297  primaryMesh(),
300  false
301  ),
302  regionMesh(),
303  dimensionedScalar("rho", dimDensity, solid().rho()),
304  zeroGradientFaPatchScalarField::typeName
305  )
306  );
307 }
308 
309 
311 {}
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 } // End namespace regionModels
317 } // End namespace Foam
318 
319 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::regionModels::regionFaModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionFaModelI.H:106
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::regionModels::KirchhoffShell::rho
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
Definition: KirchhoffShell.C:287
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:325
Foam::regionModels::KirchhoffShell::KirchhoffShell
KirchhoffShell(const word &modelType, const fvPatch &patch, const dictionary &dict)
Construct from components and dict.
Definition: KirchhoffShell.C:136
vibrationShellModel
Intermediate class for vibration-shell finite-area models.
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::faMatrix
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:59
Foam::fac::d2dt2
tmp< GeometricField< Type, faPatchField, areaMesh > > d2dt2(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facD2dt2.C:46
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::regionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(KirchhoffShell, 0)
KirchhoffShell.H
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::regionModels::regionFaModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionFaModelI.H:33
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimDensity
const dimensionSet dimDensity
Foam::fam::d2dt2
tmp< faMatrix< Type > > d2dt2(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famD2dt2.C:47
subCycle.H
Foam::regionModels::vibrationShellModel::w_
areaScalarField w_
Shell displacement.
Definition: vibrationShellModel.H:142
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:850
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dimForce
const dimensionSet dimForce
Foam::subCycleTime
A class for managing sub-cycling times.
Definition: subCycleTime.H:50
Foam::regionModels::vibrationShellModel::solid
const solidProperties & solid() const
Return solid properties.
Definition: vibrationShellModel.C:137
Foam::regionModels::KirchhoffShell::info
virtual void info()
Provide some feedback.
Definition: KirchhoffShell.C:310
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::faMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: faMatrixSolve.C:55
Foam::fam::ddt
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:48
Foam::regionModels::KirchhoffShell::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
Definition: KirchhoffShell.C:262
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::regionModels::vibrationShellModel
Definition: vibrationShellModel.H:118
Foam::regionModels::KirchhoffShell::read
virtual bool read(const dictionary &dict)
Read control parameters from dictionary.
Definition: KirchhoffShell.C:49
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::regionModels::KirchhoffShell::laplace2W0_
areaScalarField laplace2W0_
Cache laplace2.oldTime() in sub-cycling.
Definition: KirchhoffShell.H:167
Foam::regionModels::KirchhoffShell::D
const tmp< areaScalarField > D() const
Return stiffness.
Definition: KirchhoffShell.C:278
Foam::regionModels::KirchhoffShell::nSubCycles_
label nSubCycles_
Sub cycles.
Definition: KirchhoffShell.H:140
Foam::regionModels::KirchhoffShell::laplace2W_
areaScalarField laplace2W_
Laplace of the Laplace for the displacement.
Definition: KirchhoffShell.H:155
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::regionModels::KirchhoffShell::solveDisplacement
void solveDisplacement()
Solve energy equation.
Definition: KirchhoffShell.C:56
Foam::regionModels::vibrationShellModel::a_
areaScalarField a_
Shell acceleration.
Definition: vibrationShellModel.H:145
KirchhoffShell
Foam::regionModels::KirchhoffShell::w00_
areaScalarField w00_
Cache w.oldTime.oldTime() in sub-cycling.
Definition: KirchhoffShell.H:161
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
init
mesh init(true)
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::fa::optionList::constrain
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
Definition: faOptionListTemplates.C:224
Foam::regionModels::KirchhoffShell::nNonOrthCorr_
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: KirchhoffShell.H:137
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::regionModels::KirchhoffShell::laplaceW0_
areaScalarField laplaceW0_
Cache laplaceW.oldTime() in sub-cycling.
Definition: KirchhoffShell.H:164
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::regionFaModel::regionMesh
const faMesh & regionMesh() const
Return the region mesh database.
Definition: regionFaModelI.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::regionFaModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionFaModelI.H:39
Foam::regionModels::KirchhoffShell::laplaceW_
areaScalarField laplaceW_
Laplace of the displacement.
Definition: KirchhoffShell.H:152
Foam::regionModels::KirchhoffShell::w0_
areaScalarField w0_
Cache w.oldTime() in sub-cycling.
Definition: KirchhoffShell.H:158
zeroGradientFaPatchFields.H
Foam::regionModels::KirchhoffShell::ps_
areaScalarField ps_
External surface source [Pa].
Definition: KirchhoffShell.H:146
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
fvPatchFields.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::regionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary)
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
Foam::regionModels::KirchhoffShell::evolveRegion
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: KirchhoffShell.C:266
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
w0
#define w0
Definition: blockCreate.C:33
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::regionModels::vibrationShellModel::faOptions
Foam::fa::options & faOptions()
Return faOptions.
Definition: vibrationShellModel.C:131
Foam::GeometricField< scalar, faPatchField, areaMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::regionModels::KirchhoffShell::h_
areaScalarField h_
Thickness [m].
Definition: KirchhoffShell.H:149
Foam::fa::optionList::correct
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
Definition: faOptionListTemplates.C:257
Foam::IOobject::MUST_READ
Definition: IOobject.H:120