thermocapillaryForce.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) 2011-2014 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 
26 #include "thermocapillaryForce.H"
28 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  surfaceFilmModel& owner,
49  const dictionary& dict
50 )
51 :
52  force(owner)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
63 
65 {
66  const volScalarField& alpha = owner_.alpha();
67  const volScalarField& sigma = owner_.sigma();
68 
71 
72  tfvm() += alpha*fvc::grad(sigma);
73 
74  return tfvm;
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 } // End namespace surfaceFilmModels
81 } // End namespace regionModels
82 } // End namespace Foam
83 
84 // ************************************************************************* //
thermocapillaryForce
Thermocapillary force.
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::sigma
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
Foam::regionModels::surfaceFilmModels::thermocapillaryForce::~thermocapillaryForce
virtual ~thermocapillaryForce()
Destructor.
Definition: thermocapillaryForce.C:58
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::dimForce
const dimensionSet dimForce
U
U
Definition: pEqn.H:46
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::regionModels::surfaceFilmModels::thermocapillaryForce::thermocapillaryForce
thermocapillaryForce(const thermocapillaryForce &)
Disallow default bitwise copy construct.
thermocapillaryForce.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::regionModels::surfaceFilmModels::thermocapillaryForce::correct
virtual tmp< fvVectorMatrix > correct(volVectorField &U)
Correct.
Definition: thermocapillaryForce.C:64
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::regionModels::surfaceFilmModels::force
Base class for film (stress-based) force models.
Definition: force.H:55
fvcGrad.H
Calculate the gradient of the given field.
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner_
surfaceFilmModel & owner_
Reference to the owner surface film model.
Definition: filmSubModelBase.H:63
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::alpha
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].