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 "fixedPressureCompressibleDensityFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/surfaceFields.H>
00030 #include <finiteVolume/volFields.H>
00031
00032
00033
00034 namespace Foam
00035 {
00036
00037
00038
00039 fixedPressureCompressibleDensityFvPatchScalarField::
00040 fixedPressureCompressibleDensityFvPatchScalarField
00041 (
00042 const fvPatch& p,
00043 const DimensionedField<scalar, volMesh>& iF
00044 )
00045 :
00046 fixedValueFvPatchField<scalar>(p, iF),
00047 pName_("pNameIsUndefined")
00048 {}
00049
00050
00051 fixedPressureCompressibleDensityFvPatchScalarField::
00052 fixedPressureCompressibleDensityFvPatchScalarField
00053 (
00054 const fixedPressureCompressibleDensityFvPatchScalarField& ptf,
00055 const fvPatch& p,
00056 const DimensionedField<scalar, volMesh>& iF,
00057 const fvPatchFieldMapper& mapper
00058 )
00059 :
00060 fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
00061 pName_(ptf.pName_)
00062 {}
00063
00064
00065 fixedPressureCompressibleDensityFvPatchScalarField::
00066 fixedPressureCompressibleDensityFvPatchScalarField
00067 (
00068 const fvPatch& p,
00069 const DimensionedField<scalar, volMesh>& iF,
00070 const dictionary& dict
00071 )
00072 :
00073 fixedValueFvPatchField<scalar>(p, iF, dict),
00074 pName_(dict.lookup("p"))
00075 {}
00076
00077
00078 fixedPressureCompressibleDensityFvPatchScalarField::
00079 fixedPressureCompressibleDensityFvPatchScalarField
00080 (
00081 const fixedPressureCompressibleDensityFvPatchScalarField& ptf
00082 )
00083 :
00084 fixedValueFvPatchField<scalar>(ptf),
00085 pName_(ptf.pName_)
00086 {}
00087
00088
00089 fixedPressureCompressibleDensityFvPatchScalarField::
00090 fixedPressureCompressibleDensityFvPatchScalarField
00091 (
00092 const fixedPressureCompressibleDensityFvPatchScalarField& ptf,
00093 const DimensionedField<scalar, volMesh>& iF
00094 )
00095 :
00096 fixedValueFvPatchField<scalar>(ptf, iF),
00097 pName_(ptf.pName_)
00098 {}
00099
00100
00101
00102
00103 void fixedPressureCompressibleDensityFvPatchScalarField::updateCoeffs()
00104 {
00105 if (updated())
00106 {
00107 return;
00108 }
00109
00110 const fvPatchField<scalar>& pp =
00111 patch().lookupPatchField<volScalarField, scalar>(pName_);
00112
00113 const dictionary& thermoProps =
00114 db().lookupObject<IOdictionary>("thermodynamicProperties");
00115
00116 const scalar rholSat =
00117 dimensionedScalar(thermoProps.lookup("rholSat")).value();
00118
00119 const scalar pSat =
00120 dimensionedScalar(thermoProps.lookup("pSat")).value();
00121
00122 const scalar psil = dimensionedScalar(thermoProps.lookup("psil")).value();
00123
00124 operator==(rholSat + psil*(pp - pSat));
00125
00126 fixedValueFvPatchField<scalar>::updateCoeffs();
00127 }
00128
00129
00130 void fixedPressureCompressibleDensityFvPatchScalarField::write
00131 (
00132 Ostream& os
00133 ) const
00134 {
00135 fvPatchField<scalar>::write(os);
00136 os.writeKeyword("p") << pName_ << token::END_STATEMENT << nl;
00137 writeEntry("value", os);
00138 }
00139
00140
00141
00142
00143 makePatchTypeField
00144 (
00145 fvPatchScalarField,
00146 fixedPressureCompressibleDensityFvPatchScalarField
00147 );
00148
00149
00150
00151
00152 }
00153
00154