perturbedTemperatureDependentContactAngleForce.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) 2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace surfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(perturbedTemperatureDependentContactAngleForce, 0);
44 (
45  force,
46  perturbedTemperatureDependentContactAngleForce,
47  dictionary
48 );
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 perturbedTemperatureDependentContactAngleForce::
54 perturbedTemperatureDependentContactAngleForce
55 (
57  const dictionary& dict
58 )
59 :
60  contactAngleForce(typeName, film, dict),
61  thetaPtr_(Function1<scalar>::New("theta", coeffDict_, &film.regionMesh())),
62  rndGen_(label(0)),
63  distribution_
64  (
66  (
67  coeffDict_.subDict("distribution"),
68  rndGen_
69  )
70  )
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 
85 {
86  tmp<volScalarField> ttheta
87  (
88  new volScalarField
89  (
90  IOobject
91  (
92  typeName + ":theta",
95  ),
98  )
99  );
100 
101  volScalarField& theta = ttheta.ref();
102  volScalarField::Internal& thetai = theta.ref();
103 
104  const volScalarField& T = filmModel_.T();
105 
106  // Initialize with the function of temperature
107  thetai.field() = thetaPtr_->value(T());
108 
109  // Add the stochastic perturbation
110  forAll(thetai, celli)
111  {
112  thetai[celli] += distribution_->sample();
113  }
114 
115  forAll(theta.boundaryField(), patchi)
116  {
117  if (!filmModel_.isCoupledPatch(patchi))
118  {
119  fvPatchField<scalar>& thetaf = theta.boundaryFieldRef()[patchi];
120 
121  // Initialize with the function of temperature
122  thetaf = thetaPtr_->value(T.boundaryField()[patchi]);
123 
124  // Add the stochastic perturbation
125  forAll(thetaf, facei)
126  {
127  thetaf[facei] += distribution_->sample();
128  }
129  }
130  }
131 
132  return ttheta;
133 }
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 } // End namespace surfaceFilmModels
138 } // End namespace regionModels
139 } // End namespace Foam
140 
141 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DimensionedField::field
const Field< Type > & field() const
Return field.
Definition: DimensionedFieldI.H:90
Foam::regionModels::surfaceFilmModels::contactAngleForce
Base-class for film contact angle force models.
Definition: contactAngleForce.H:57
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::regionModels::surfaceFilmModels::perturbedTemperatureDependentContactAngleForce::~perturbedTemperatureDependentContactAngleForce
virtual ~perturbedTemperatureDependentContactAngleForce()
Destructor.
Definition: perturbedTemperatureDependentContactAngleForce.C:77
Foam::regionModels::surfaceFilmModels::perturbedTemperatureDependentContactAngleForce::theta
virtual tmp< volScalarField > theta() const
Return the contact angle field.
Definition: perturbedTemperatureDependentContactAngleForce.C:84
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::regionModels::regionModel::isCoupledPatch
bool isCoupledPatch(const label regionPatchi) const
Return true if patchi on the local region is a coupled.
Definition: regionModelI.H:134
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::T
virtual const volScalarField & T() const =0
Return the film mean temperature [K].
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distributionModel::New
static autoPtr< distributionModel > New(const dictionary &dict, Random &rndGen)
Selector.
Definition: distributionModelNew.C:34
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::regionModels::surfaceFilmModels::filmSubModelBase::filmModel_
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Definition: filmSubModelBase.H:65
Foam::GeometricField< scalar, fvPatchField, volMesh >
perturbedTemperatureDependentContactAngleForce.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189