Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "dynamicAlphaContactAngleFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volMesh.H>
00030
00031
00032
00033 Foam::dynamicAlphaContactAngleFvPatchScalarField::
00034 dynamicAlphaContactAngleFvPatchScalarField
00035 (
00036 const fvPatch& p,
00037 const DimensionedField<scalar, volMesh>& iF
00038 )
00039 :
00040 alphaContactAngleFvPatchScalarField(p, iF),
00041 theta0_(0.0),
00042 uTheta_(0.0),
00043 thetaA_(0.0),
00044 thetaR_(0.0)
00045 {}
00046
00047
00048 Foam::dynamicAlphaContactAngleFvPatchScalarField::
00049 dynamicAlphaContactAngleFvPatchScalarField
00050 (
00051 const dynamicAlphaContactAngleFvPatchScalarField& gcpsf,
00052 const fvPatch& p,
00053 const DimensionedField<scalar, volMesh>& iF,
00054 const fvPatchFieldMapper& mapper
00055 )
00056 :
00057 alphaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
00058 theta0_(gcpsf.theta0_),
00059 uTheta_(gcpsf.uTheta_),
00060 thetaA_(gcpsf.thetaA_),
00061 thetaR_(gcpsf.thetaR_)
00062 {}
00063
00064
00065 Foam::dynamicAlphaContactAngleFvPatchScalarField::
00066 dynamicAlphaContactAngleFvPatchScalarField
00067 (
00068 const fvPatch& p,
00069 const DimensionedField<scalar, volMesh>& iF,
00070 const dictionary& dict
00071 )
00072 :
00073 alphaContactAngleFvPatchScalarField(p, iF, dict),
00074 theta0_(readScalar(dict.lookup("theta0"))),
00075 uTheta_(readScalar(dict.lookup("uTheta"))),
00076 thetaA_(readScalar(dict.lookup("thetaA"))),
00077 thetaR_(readScalar(dict.lookup("thetaR")))
00078 {
00079 evaluate();
00080 }
00081
00082
00083 Foam::dynamicAlphaContactAngleFvPatchScalarField::
00084 dynamicAlphaContactAngleFvPatchScalarField
00085 (
00086 const dynamicAlphaContactAngleFvPatchScalarField& gcpsf
00087 )
00088 :
00089 alphaContactAngleFvPatchScalarField(gcpsf),
00090 theta0_(gcpsf.theta0_),
00091 uTheta_(gcpsf.uTheta_),
00092 thetaA_(gcpsf.thetaA_),
00093 thetaR_(gcpsf.thetaR_)
00094 {}
00095
00096
00097 Foam::dynamicAlphaContactAngleFvPatchScalarField::
00098 dynamicAlphaContactAngleFvPatchScalarField
00099 (
00100 const dynamicAlphaContactAngleFvPatchScalarField& gcpsf,
00101 const DimensionedField<scalar, volMesh>& iF
00102 )
00103 :
00104 alphaContactAngleFvPatchScalarField(gcpsf, iF),
00105 theta0_(gcpsf.theta0_),
00106 uTheta_(gcpsf.uTheta_),
00107 thetaA_(gcpsf.thetaA_),
00108 thetaR_(gcpsf.thetaR_)
00109 {}
00110
00111
00112
00113
00114 Foam::tmp<Foam::scalarField>
00115 Foam::dynamicAlphaContactAngleFvPatchScalarField::theta
00116 (
00117 const fvPatchVectorField& Up,
00118 const fvsPatchVectorField& nHat
00119 ) const
00120 {
00121 if (uTheta_ < SMALL)
00122 {
00123 return tmp<scalarField>(new scalarField(size(), theta0_));
00124 }
00125
00126 vectorField nf = patch().nf();
00127
00128
00129 vectorField Uwall = Up.patchInternalField() - Up;
00130 Uwall -= (nf & Uwall)*nf;
00131
00132
00133 vectorField nWall = nHat - (nf & nHat)*nf;
00134
00135
00136 nWall /= (mag(nWall) + SMALL);
00137
00138
00139
00140 scalarField uwall = nWall & Uwall;
00141
00142 return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
00143 }
00144
00145
00146 void Foam::dynamicAlphaContactAngleFvPatchScalarField::write(Ostream& os) const
00147 {
00148 alphaContactAngleFvPatchScalarField::write(os);
00149 os.writeKeyword("theta0") << theta0_ << token::END_STATEMENT << nl;
00150 os.writeKeyword("uTheta") << uTheta_ << token::END_STATEMENT << nl;
00151 os.writeKeyword("thetaA") << thetaA_ << token::END_STATEMENT << nl;
00152 os.writeKeyword("thetaR") << thetaR_ << token::END_STATEMENT << nl;
00153 writeEntry("value", os);
00154 }
00155
00156
00157
00158
00159 namespace Foam
00160 {
00161 makePatchTypeField
00162 (
00163 fvPatchScalarField,
00164 dynamicAlphaContactAngleFvPatchScalarField
00165 );
00166 }
00167
00168
00169