Loading...
Searching...
No Matches
KirchhoffShell.H
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
26Class
27 Foam::regionModels::KirchhoffShell (Foam::regionFaModels)
28
29Description
30 Vibration-shell finite-area model.
31
32Usage
33 Example of the boundary condition specification:
34 \verbatim
35 <patchName>
36 {
37 // Mandatory entries
38 vibrationShellModel KirchhoffShell;
39 f0 <scalar>;
40 f1 <scalar>;
41 f2 <scalar>;
42
43 // Inherited entries
44 ...
45 nNonOrthCorr <int>; // read from another dict
46 nSubCycles <int>;
47 }
48 \endverbatim
49
50 where the entries mean:
51 \table
52 Property | Description | Type | Reqd | Deflt
53 vibrationShellModel | Type name: KirchhoffShell | word | yes | -
54 f0 | Damping coefficient [1/s] | scalar | yes | -
55 f1 | Damping coefficient [1/s] | scalar | yes | -
56 f2 | Damping coefficient [1/s] | scalar | yes | -
57 h | Name of thickness field | word | no | h (suffix)
58 ps | Name of pressure field | word | no | ps (suffix)
59 \endtable
60
61 Fields/variables used:
62 \table
63 Property | Description | Type | Deflt
64 h | Thickness | shell | h (suffix)
65 ps | Pressure on shell | shell | ps (suffix)
66 \endtable
67
68 Naming changes from 2056 and earlier:
69 \table
70 Keyword | Description | Keyword (old) | Deflt (old)
71 h | Thickness (shell) | - | "h_" + regionName
72 ps | Pressure (shell) | - | "ps_" + regionName
73 \endtable
74
75 The inherited entries are elaborated in:
76 - \link vibrationShellModel.H \endlink
77
78SourceFiles
79 KirchhoffShell.C
80
81\*---------------------------------------------------------------------------*/
82
83#ifndef Foam_regionModels_KirchhoffShell_H
84#define Foam_regionModels_KirchhoffShell_H
85
86#include "volFieldsFwd.H"
87#include "vibrationShellModel.H"
88#include "faCFD.H"
89
90// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91
92namespace Foam
93{
94namespace regionModels
95{
96
97/*---------------------------------------------------------------------------*\
98 Class KirchhoffShell Declaration
99\*---------------------------------------------------------------------------*/
100
101class KirchhoffShell
102:
104{
105 // Private Data
106
107 // Source term fields
108
109 //- Shell thickness [m]
111
112 //- External surface source [Pa]
113 const areaScalarField ps_;
114
115 //- Laplace of the displacement
116 areaScalarField laplaceW_;
117
118 //- Laplace of the Laplace for the displacement
119 areaScalarField laplace2W_;
120
121 //- Cache w.oldTime() in sub-cycling
122 areaScalarField w0_;
123
124 //- Cache w.oldTime.oldTime() in sub-cycling
125 areaScalarField w00_;
126
127 //- Cache laplaceW.oldTime() in sub-cycling
128 areaScalarField laplaceW0_;
129
130 //- Cache laplace2.oldTime() in sub-cycling
131 areaScalarField laplace2W0_;
132
133
134 // Solution parameters
135
136 //- Damping coefficients [1/s]
137 const dimensionedScalar f0_;
138 const dimensionedScalar f1_;
139 const dimensionedScalar f2_;
140
141 //- Number of non orthogonal correctors
142 label nNonOrthCorr_;
143
144 //- Sub cycles
145 label nSubCycles_;
146
147
148 // Private Member Functions
149
150 //- Initialise Kirchhoff shell model
151 bool init(const dictionary& dict);
152
153 //- Solve energy equation
154 void solveDisplacement();
155
156
157public:
158
159 //- Runtime type information
160 TypeName("KirchhoffShell");
161
162
163 // Constructors
164
165 //- Construct from components and dict
167 (
168 const word& modelType,
169 const fvMesh& mesh,
170 const dictionary& dict
171 );
172
173 //- No copy construct
174 KirchhoffShell(const KirchhoffShell&) = delete;
175
176 //- No copy assignment
177 void operator=(const KirchhoffShell&) = delete;
178
179
180 //- Destructor
181 virtual ~KirchhoffShell() = default;
182
183
184 // Member Functions
185
186 // Fields
187
188 //- Return stiffness
189 const tmp<areaScalarField> D() const;
190
191 //- Return density [Kg/m3]
192 const tmp<areaScalarField> rho() const;
193
194
195 // Evolution
196
197 //- Pre-evolve thermal baffle
198 virtual void preEvolveRegion();
199
200 //- Evolve the thermal baffle
201 virtual void evolveRegion();
202
203
204 // IO
205
206 //- Provide some feedback
207 virtual void info();
208};
209
210
211// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212
213} // End namespace regionModels
214} // End namespace Foam
215
216// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217
218
219#endif
220
221// ************************************************************************* //
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Vibration-shell finite-area model.
KirchhoffShell(const KirchhoffShell &)=delete
No copy construct.
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
virtual ~KirchhoffShell()=default
Destructor.
void operator=(const KirchhoffShell &)=delete
No copy assignment.
KirchhoffShell(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components and dict.
TypeName("KirchhoffShell")
Runtime type information.
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
virtual void info()
Provide some feedback.
const tmp< areaScalarField > D() const
Return stiffness.
virtual void evolveRegion()
Evolve the thermal baffle.
Intermediate class for vibration-shell finite-area models.
vibrationShellModel(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from type name and mesh and dict.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68
Forwards and collection of common volume field types.