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 -------------------------------------------------------------------------------
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(distributionContactAngleForce, 0);
43 addToRunTimeSelectionTable(force, distributionContactAngleForce, dictionary);
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 distributionContactAngleForce::distributionContactAngleForce
49 (
51  const dictionary& dict
52 )
53 :
54  contactAngleForce(typeName, film, dict),
55  rndGen_(),
56  distribution_
57  (
59  (
60  coeffDict_.subDict("distribution"),
61  rndGen_
62  )
63  )
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
70 {}
71 
72 
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
74 
76 {
77  tmp<volScalarField> ttheta
78  (
79  new volScalarField
80  (
81  IOobject
82  (
83  typeName + ":theta",
86  ),
89  )
90  );
91 
92  volScalarField& theta = ttheta.ref();
93  volScalarField::Internal& thetai = theta.ref();
94 
95  forAll(thetai, celli)
96  {
97  thetai[celli] = distribution_->sample();
98  }
99 
100  forAll(theta.boundaryField(), patchi)
101  {
102  if (!filmModel_.isCoupledPatch(patchi))
103  {
104  fvPatchField<scalar>& thetaf = theta.boundaryFieldRef()[patchi];
105 
106  forAll(thetaf, facei)
107  {
108  thetaf[facei] = distribution_->sample();
109  }
110  }
111  }
112 
113  return ttheta;
114 }
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 } // End namespace surfaceFilmModels
119 } // End namespace regionModels
120 } // End namespace Foam
121 
122 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
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::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)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
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::distributionContactAngleForce::~distributionContactAngleForce
virtual ~distributionContactAngleForce()
Destructor.
Definition: distributionContactAngleForce.C:69
distributionContactAngleForce.H
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
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:121
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::distributionContactAngleForce::theta
virtual tmp< volScalarField > theta() const
Return the contact angle field.
Definition: distributionContactAngleForce.C:75
Foam::regionModels::surfaceFilmModels::filmSubModelBase::filmModel_
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Definition: filmSubModelBase.H:65
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54