tabulatedAccelerationSourceTemplates.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 | Copyright (C) 2015 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::tabulatedAccelerationSource::addSup
35 (
36  const RhoFieldType& rho,
37  fvMatrix<vector>& eqn,
38  const label fieldi
39 )
40 {
41  Vector<vector> acceleration(motion_.acceleration());
42 
43  // If gravitational force is present combine with the linear acceleration
45  {
48  (
50  );
51 
52  const uniformDimensionedScalarField& hRef =
54 
55  g = g0_ - dimensionedVector("a", dimAcceleration, acceleration.x());
56 
57  dimensionedScalar ghRef
58  (
59  mag(g.value()) > SMALL
60  ? g & (cmptMag(g.value())/mag(g.value()))*hRef
61  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
62  );
63 
64  const_cast<volScalarField&>
65  (
67  ) = (g & mesh_.C()) - ghRef;
68 
69  const_cast<surfaceScalarField&>
70  (
72  ) = (g & mesh_.Cf()) - ghRef;
73  }
74  // ... otherwise include explicitly in the momentum equation
75  else
76  {
77  eqn -= rho*dimensionedVector("a", dimAcceleration, acceleration.x());
78  }
79 
80  dimensionedVector Omega
81  (
82  "Omega",
83  dimensionSet(0, 0, -1, 0, 0),
84  acceleration.y()
85  );
86 
87  dimensionedVector dOmegaDT
88  (
89  "dOmegaDT",
90  dimensionSet(0, 0, -2, 0, 0),
91  acceleration.z()
92  );
93 
94  eqn -=
95  (
96  rho*(2*Omega ^ eqn.psi()) // Coriolis force
97  + rho*(Omega ^ (Omega ^ mesh_.C())) // Centrifugal force
98  + rho*(dOmegaDT ^ mesh_.C()) // Angular tabulatedAcceleration force
99  );
100 }
101 
102 
103 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
UniformDimensionedField< vector > uniformDimensionedVectorField
bool foundObject(const word &name) const
Is the named Type found?
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:79
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
virtual Vector< vector > acceleration() const
Return the solid-body accelerations.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const volVectorField & C() const
Return cell centres as volVectorField.
UniformDimensionedField< scalar > uniformDimensionedScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const dimensionSet dimAcceleration
tabulated6DoFAcceleration motion_
Run-time selectable acceleration model.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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 > &)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField