sixDoFAccelerationSourceTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "fvMesh.H"
28 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class RhoFieldType>
34 void Foam::fv::sixDoFAccelerationSource::addSup
35 (
36  const RhoFieldType& rho,
37  fvMatrix<vector>& eqn,
38  const word& fieldName
39 ) const
40 {
41  Vector<vector> accelerations(accelerations_->value(mesh().time().value()));
42 
43  // If gravitational force is present combine with the linear acceleration
45  {
48 
49  const uniformDimensionedScalarField& hRef =
51 
52  g = g_ - dimensionedVector("a", dimAcceleration, accelerations.x());
53 
54  dimensionedScalar ghRef(- mag(g)*hRef);
55 
56  mesh().lookupObjectRef<volScalarField>("gh") = (g & mesh().C()) - ghRef;
57 
59  (g & mesh().Cf()) - ghRef;
60  }
61  // ... otherwise include explicitly in the momentum equation
62  else
63  {
64  eqn -= rho*dimensionedVector("a", dimAcceleration, accelerations.x());
65  }
66 
67  dimensionedVector Omega
68  (
69  "Omega",
70  dimensionSet(0, 0, -1, 0, 0),
71  accelerations.y()
72  );
73 
74  dimensionedVector dOmegaDT
75  (
76  "dOmegaDT",
77  dimensionSet(0, 0, -2, 0, 0),
78  accelerations.z()
79  );
80 
81  eqn -=
82  (
83  rho*(2*Omega ^ eqn.psi()) // Coriolis force
84  + rho*(Omega ^ (Omega ^ mesh().C())) // Centrifugal force
85  + rho*(dOmegaDT ^ mesh().C()) // Angular acceleration force
86  );
87 }
88 
89 
90 // ************************************************************************* //
UniformDimensionedField< vector > uniformDimensionedVectorField
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
bool foundObject(const word &name) const
Is the named Type found?
UniformDimensionedField< scalar > uniformDimensionedScalarField
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimAcceleration
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
dimensioned< scalar > mag(const dimensioned< Type > &)
const volVectorField & C() const
Return cell centres as volVectorField.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField