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-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(perturbedTemperatureDependentContactAngleForce, 0);
42 (
43  force,
44  perturbedTemperatureDependentContactAngleForce,
45  dictionary
46 );
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
53 (
55  const dictionary& dict
56 )
57 :
58  contactAngleForce(typeName, film, dict),
59  thetaPtr_(Function1<scalar>::New("theta", coeffDict_)),
60  rndGen_(label(0)),
61  distribution_
62  (
64  (
65  coeffDict_.subDict("distribution"),
66  rndGen_
67  )
68  )
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
76 {}
77 
78 
79 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
80 
83 {
84  tmp<volScalarField> ttheta
85  (
87  (
88  typeName + ":theta",
91  )
92  );
93 
94  volScalarField& theta = ttheta.ref();
95  volScalarField::Internal& thetai = theta.ref();
96 
97  const volScalarField& T = filmModel_.T();
98 
99  // Initialize with the function of temperature
100  thetai.field() = thetaPtr_->value(T());
101 
102  // Add the stochastic perturbation
103  forAll(thetai, celli)
104  {
105  thetai[celli] += distribution_->sample();
106  }
107 
108  forAll(theta.boundaryField(), patchi)
109  {
111  {
112  fvPatchField<scalar>& thetaf = theta.boundaryFieldRef()[patchi];
113 
114  // Initialize with the function of temperature
115  thetaf = thetaPtr_->value(T.boundaryField()[patchi]);
116 
117  // Add the stochastic perturbation
118  forAll(thetaf, facei)
119  {
120  thetaf[facei] += distribution_->sample();
121  }
122  }
123  }
124 
125  return ttheta;
126 }
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
130 } // End namespace surfaceFilmModels
131 } // End namespace regionModels
132 } // End namespace Foam
133 
134 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
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
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.
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
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.
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].