perturbedTemperatureDependentContactAngleForce.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-2019 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(perturbedTemperatureDependentContactAngleForce, 0);
43 (
44  force,
45  perturbedTemperatureDependentContactAngleForce,
46  dictionary
47 );
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
56  const dictionary& dict
57 )
58 :
59  contactAngleForce(typeName, film, dict),
60  thetaPtr_(Function1<scalar>::New("theta", coeffDict_)),
61  rndGen_(label(0)),
62  distribution_
63  (
65  (
66  coeffDict_.subDict("distribution"),
67  rndGen_
68  )
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
77 {}
78 
79 
80 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
81 
84 {
85  tmp<volScalarField> ttheta
86  (
88  (
89  IOobject::modelName("theta", typeName),
92  )
93  );
94 
95  volScalarField& theta = ttheta.ref();
96  volScalarField::Internal& thetai = theta.ref();
97 
98  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
99 
100  const volScalarField& T = film.T();
101 
102  // Initialize with the function of temperature
103  thetai.field() = thetaPtr_->value(T());
104 
105  // Add the stochastic perturbation
106  forAll(thetai, celli)
107  {
108  thetai[celli] += distribution_->sample();
109  }
110 
111  forAll(theta.boundaryField(), patchi)
112  {
114  {
115  fvPatchField<scalar>& thetaf = theta.boundaryFieldRef()[patchi];
116 
117  // Initialize with the function of temperature
118  thetaf = thetaPtr_->value(T.boundaryField()[patchi]);
119 
120  // Add the stochastic perturbation
121  forAll(thetaf, facei)
122  {
123  thetaf[facei] += distribution_->sample();
124  }
125  }
126  }
127 
128  return ttheta;
129 }
130 
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 
133 } // End namespace surfaceFilmModels
134 } // End namespace regionModels
135 } // End namespace Foam
136 
137 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:52
perturbedTemperatureDependentContactAngleForce(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
virtual Type value(const scalar x) const =0
Return value as a function of (scalar) independent variable.
bool isCoupledPatch(const label regionPatchi) const
Return true if patchi on the local region is a coupled.
Definition: regionModelI.H:138
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
static autoPtr< distributionModel > New(const dictionary &dict, Random &rndGen)
Selector.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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)
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
virtual const volScalarField & T() const
Return the film mean temperature [K].