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-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 "thermoSingleLayer.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(temperatureDependentContactAngleForce, 0);
43 (
44  force,
45  temperatureDependentContactAngleForce,
46  dictionary
47 );
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
55  const dictionary& dict
56 )
57 :
58  contactAngleForce(typeName, film, dict),
59  thetaPtr_(Function1<scalar>::New("theta", coeffDict_))
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
70 
72 {
73  tmp<volScalarField> ttheta
74  (
76  (
77  IOobject::modelName("theta", typeName),
80  )
81  );
82 
83  volScalarField& theta = ttheta.ref();
84 
85  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
86 
87  const volScalarField& T = film.thermo().T();
88 
89  theta.ref().field() = thetaPtr_->value(T());
90 
91  forAll(theta.boundaryField(), patchi)
92  {
94  {
95  theta.boundaryFieldRef()[patchi] =
96  thetaPtr_->value(T.boundaryField()[patchi]);
97  }
98  }
99 
100  return ttheta;
101 }
102 
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105 } // End namespace surfaceFilmModels
106 } // End namespace regionModels
107 } // End namespace Foam
108 
109 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
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:156
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:181
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
virtual tmp< volScalarField > theta() const
Return the contact angle field.
Macros for easy insertion into run-time selection tables.
virtual Type value(const scalar x) const =0
Return value as a function of scalar x.
bool isCoupledPatch(const label regionPatchi) const
Return true if patchi on the local region is a coupled.
Definition: regionModelI.H:132
temperatureDependentContactAngleForce(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:55
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
virtual const volScalarField & T() const =0
Temperature [K].
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.
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
A class for managing temporary objects.
Definition: PtrList.H:53
defineTypeNameAndDebug(kinematicSingleLayer, 0)
addToRunTimeSelectionTable(surfaceFilmModel, noFilm, mesh)
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.