temperatureDependentContactAngleForce.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) 2017-2018 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 
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace regionModels
34 {
35 namespace surfaceFilmModels
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(temperatureDependentContactAngleForce, 0);
42 (
43  force,
44  temperatureDependentContactAngleForce,
45  dictionary
46 );
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
54  const dictionary& dict
55 )
56 :
57  contactAngleForce(typeName, film, dict),
58  thetaPtr_(Function1<scalar>::New("theta", coeffDict_))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 
71 {
72  tmp<volScalarField> ttheta
73  (
75  (
76  typeName + ":theta",
79  )
80  );
81 
82  volScalarField& theta = ttheta.ref();
83 
84  const volScalarField& T = filmModel_.T();
85 
86  theta.ref().field() = thetaPtr_->value(T());
87 
88  forAll(theta.boundaryField(), patchi)
89  {
91  {
92  theta.boundaryFieldRef()[patchi] =
93  thetaPtr_->value(T.boundaryField()[patchi]);
94  }
95  }
96 
97  return ttheta;
98 }
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 } // End namespace surfaceFilmModels
103 } // End namespace regionModels
104 } // End namespace Foam
105 
106 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual tmp< volScalarField > theta() const
Return the contact angle field.
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
bool isCoupledPatch(const label regionPatchi) const
Return true if patchi on the local region is a coupled.
Definition: regionModelI.H:138
temperatureDependentContactAngleForce(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
Base-class for film contact angle force models.
const Field< Type > & field() const
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
A class for managing temporary objects.
Definition: PtrList.H:53
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.
virtual const volScalarField & T() const =0
Return the film mean temperature [K].