Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
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
26\*---------------------------------------------------------------------------*/
27
30#include "fvcGrad.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace regionModels
38namespace surfaceFilmModels
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43defineTypeNameAndDebug(thermocapillaryForce, 0);
44addToRunTimeSelectionTable(force, thermocapillaryForce, dictionary);
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48thermocapillaryForce::thermocapillaryForce
49(
51 const dictionary& dict
52)
54 force(film)
55{}
56
57
58// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61{}
62
63
64// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65
67{
70
72
73 tfvm.ref() += alpha*fvc::grad(sigma);
74
75 return tfvm;
76}
77
78
79// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80
81} // End namespace surfaceFilmModels
82} // End namespace regionModels
83} // End namespace Foam
84
85// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Base class for film (stress-based) force models.
Definition force.H:56
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
virtual tmp< fvVectorMatrix > correct(volVectorField &U)
Correct.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
Calculate the gradient of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimForce
const dimensionSet dimVolume(pow3(dimLength))
volScalarField & alpha
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)