constantAlphaContactAngleFvPatchScalarField.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) 2011-2013 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 #include "volMesh.H"
29 #include "fvPatchFieldMapper.H"
30 
31 // Symbol to force loading at runtime
32 extern "C"
34 {}
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
46  theta0_(0.0)
47 {}
48 
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  alphaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
60  theta0_(gcpsf.theta0_)
61 {}
62 
63 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
73  theta0_(readScalar(dict.lookup("theta0")))
74 {
75  evaluate();
76 }
77 
78 
81 (
83 )
84 :
86  theta0_(gcpsf.theta0_)
87 {}
88 
89 
92 (
95 )
96 :
98  theta0_(gcpsf.theta0_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
106 (
107  const fvPatchVectorField&,
108  const fvsPatchVectorField&
109 ) const
110 {
111  return tmp<scalarField>(new scalarField(size(), theta0_));
112 }
113 
114 
116 (
117  Ostream& os
118 ) const
119 {
121  os.writeKeyword("theta0") << theta0_ << token::END_STATEMENT << nl;
122  writeEntry("value", os);
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 namespace Foam
129 {
131  (
134  );
135 }
136 
137 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
twoPhasePropertiesLoad
void twoPhasePropertiesLoad()
Definition: constantAlphaContactAngleFvPatchScalarField.C:33
p
p
Definition: pEqn.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::constantAlphaContactAngleFvPatchScalarField::constantAlphaContactAngleFvPatchScalarField
constantAlphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: constantAlphaContactAngleFvPatchScalarField.C:40
volMesh.H
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
fvPatchFieldMapper.H
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::constantAlphaContactAngleFvPatchScalarField::theta
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const
Return the equilibrium contact-angle.
Definition: constantAlphaContactAngleFvPatchScalarField.C:106
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::constantAlphaContactAngleFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: constantAlphaContactAngleFvPatchScalarField.C:116
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::alphaContactAngleFvPatchScalarField
Abstract base class for alphaContactAngle boundary conditions.
Definition: alphaContactAngleFvPatchScalarField.H:81
readScalar
#define readScalar
Definition: doubleScalar.C:38
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::constantAlphaContactAngleFvPatchScalarField::theta0_
scalar theta0_
Equilibrium contact angle.
Definition: constantAlphaContactAngleFvPatchScalarField.H:56
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::constantAlphaContactAngleFvPatchScalarField
A constant alphaContactAngle scalar boundary condition (alphaContactAngleFvPatchScalarField)
Definition: constantAlphaContactAngleFvPatchScalarField.H:49
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
write
Tcoeff write()
constantAlphaContactAngleFvPatchScalarField.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51