distributionContactAngleForce.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 | Copyright (C) 2017 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(distributionContactAngleForce, 0);
41 addToRunTimeSelectionTable(force, distributionContactAngleForce, dictionary);
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 distributionContactAngleForce::distributionContactAngleForce
47 (
48  surfaceFilmModel& film,
49  const dictionary& dict
50 )
51 :
52  contactAngleForce(typeName, film, dict),
53  rndGen_(label(0), -1),
54  distribution_
55  (
57  (
58  coeffDict_.subDict("distribution"),
59  rndGen_
60  )
61  )
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 {
75  tmp<volScalarField> ttheta
76  (
77  new volScalarField
78  (
79  IOobject
80  (
81  typeName + ":theta",
84  ),
86  dimensionedScalar("0", dimless, 0)
87  )
88  );
89 
90  volScalarField& theta = ttheta.ref();
91  volScalarField::Internal& thetai = theta.ref();
92 
93  forAll(thetai, celli)
94  {
95  thetai[celli] = distribution_->sample();
96  }
97 
98  forAll(theta.boundaryField(), patchi)
99  {
101  {
102  fvPatchField<scalar>& thetaf = theta.boundaryFieldRef()[patchi];
103 
104  forAll(thetaf, facei)
105  {
106  thetaf[facei] = distribution_->sample();
107  }
108  }
109  }
110 
111  return ttheta;
112 }
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 } // End namespace surfaceFilmModels
117 } // End namespace regionModels
118 } // End namespace Foam
119 
120 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
virtual tmp< volScalarField > theta() const
Return the contact angle field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
surfaceFilmModel & filmModel_
Reference to the film surface film model.
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
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
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
Base-class for film contact angle force models.
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.