Go to the documentation of this file.
38 namespace regionModels
111 if (wSubCycle.index() >= wSubCycle.nSubCycles())
124 Info<<
"ws_vibrationShell: "
125 <<
"min = " <<
min(
w_).value() <<
", "
140 const word& modelType,
180 "laplaceW_" + regionName_,
193 "laplace2W_" + regionName_,
219 "w00_" + regionName_,
232 "laplaceW0_" + regionName_,
245 "laplace2W0_" + regionName_,
302 zeroGradientFaPatchScalarField::typeName
int debug
Static debugging option.
const dictionary & solution() const
Return the solution dictionary.
const tmp< areaScalarField > rho() const
Return density [Kg/m3].
const dimensionSet dimPressure
Defines the attributes of an object for which implicit objectRegistry management is supported,...
#define InfoInFunction
Report an information message using Foam::Info.
KirchhoffShell(const word &modelType, const fvPatch &patch, const dictionary &dict)
Construct from components and dict.
Intermediate class for vibration-shell finite-area models.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A special matrix type and solver, designed for finite area solutions of scalar equations....
tmp< GeometricField< Type, faPatchField, areaMesh > > d2dt2(const dimensioned< Type > dt, const faMesh &mesh)
A class for handling words, derived from Foam::string.
const solidProperties & solid() const noexcept
Return solid properties.
defineTypeNameAndDebug(KirchhoffShell, 0)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimDensity
tmp< faMatrix< Type > > d2dt2(const GeometricField< Type, faPatchField, areaMesh > &vf)
areaScalarField w_
Shell displacement.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dimensionSet dimForce
A class for managing sub-cycling times.
virtual void info()
Provide some feedback.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
areaScalarField laplace2W0_
Cache laplace2.oldTime() in sub-cycling.
const tmp< areaScalarField > D() const
Return stiffness.
label nSubCycles_
Sub cycles.
areaScalarField laplace2W_
Laplace of the Laplace for the displacement.
const dimensionSet dimArea(sqr(dimLength))
Foam::fa::options & faOptions() noexcept
Return faOptions.
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
messageStream Info
Information stream (stdout output on master, null elsewhere)
void solveDisplacement()
Solve energy equation.
areaScalarField a_
Shell acceleration.
areaScalarField w00_
Cache w.oldTime.oldTime() in sub-cycling.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
label nNonOrthCorr_
Number of non orthogonal correctors.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
areaScalarField laplaceW0_
Cache laplaceW.oldTime() in sub-cycling.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const faMesh & regionMesh() const
Return the region mesh database.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
const Time & time() const
Return the reference to the time database.
areaScalarField laplaceW_
Laplace of the displacement.
areaScalarField w0_
Cache w.oldTime() in sub-cycling.
areaScalarField ps_
External surface source [Pa].
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedScalar sqrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(vibrationShellModel, KirchhoffShell, dictionary)
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
virtual void evolveRegion()
Evolve the thermal baffle.
const Time & time() const
Return the top-level database.
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
areaScalarField h_
Thickness [m].
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
const dimensionSet dimless
Dimensionless.