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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 namespace surfaceFilmModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 distributionContactAngleForce::distributionContactAngleForce
50 (
52  const dictionary& dict
53 )
54 :
55  contactAngleForce(typeName, film, dict),
56  rndGen_(),
57  distribution_
58  (
60  (
61  coeffDict_.subDict("distribution"),
62  rndGen_
63  )
64  )
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 {
78  auto ttheta = volScalarField::New
79  (
80  IOobject::scopedName(typeName, "theta"),
84  );
85  auto& thetaIf = ttheta.ref().internalFieldRef();
86  auto& thetaBf = ttheta.ref().boundaryFieldRef();
87 
88  forAll(thetaIf, celli)
89  {
90  thetaIf[celli] = distribution_->sample();
91  }
92 
93  forAll(thetaBf, patchi)
94  {
95  if (!filmModel_.isCoupledPatch(patchi))
96  {
97  fvPatchField<scalar>& thetaf = thetaBf[patchi];
98 
99  forAll(thetaf, facei)
100  {
101  thetaf[facei] = distribution_->sample();
102  }
103  }
104  }
105 
106  return ttheta;
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111 
112 } // End namespace surfaceFilmModels
113 } // End namespace regionModels
114 } // End namespace Foam
115 
116 // ************************************************************************* //
dictionary dict
virtual tmp< volScalarField > theta() const
Return the contact angle field.
Base class for film (stress-based) force models.
Definition: force.H:53
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimless
Dimensionless.
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
bool isCoupledPatch(const label regionPatchi) const
True if patchi on the local region is a coupled patch to the primary region.
Definition: regionModelI.H:104
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Base-class for film contact angle force models.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127