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 "tractionDisplacementCorrectionFvPatchVectorField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/volFields.H>
00029
00030
00031
00032 namespace Foam
00033 {
00034
00035
00036
00037 tractionDisplacementCorrectionFvPatchVectorField::
00038 tractionDisplacementCorrectionFvPatchVectorField
00039 (
00040 const fvPatch& p,
00041 const DimensionedField<vector, volMesh>& iF
00042 )
00043 :
00044 fixedGradientFvPatchVectorField(p, iF),
00045 traction_(p.size(), vector::zero),
00046 pressure_(p.size(), 0.0)
00047 {
00048 fvPatchVectorField::operator=(patchInternalField());
00049 gradient() = vector::zero;
00050 }
00051
00052
00053 tractionDisplacementCorrectionFvPatchVectorField::
00054 tractionDisplacementCorrectionFvPatchVectorField
00055 (
00056 const tractionDisplacementCorrectionFvPatchVectorField& tdpvf,
00057 const fvPatch& p,
00058 const DimensionedField<vector, volMesh>& iF,
00059 const fvPatchFieldMapper& mapper
00060 )
00061 :
00062 fixedGradientFvPatchVectorField(tdpvf, p, iF, mapper),
00063 traction_(tdpvf.traction_, mapper),
00064 pressure_(tdpvf.pressure_, mapper)
00065 {}
00066
00067
00068 tractionDisplacementCorrectionFvPatchVectorField::
00069 tractionDisplacementCorrectionFvPatchVectorField
00070 (
00071 const fvPatch& p,
00072 const DimensionedField<vector, volMesh>& iF,
00073 const dictionary& dict
00074 )
00075 :
00076 fixedGradientFvPatchVectorField(p, iF),
00077 traction_("traction", dict, p.size()),
00078 pressure_("pressure", dict, p.size())
00079 {
00080 fvPatchVectorField::operator=(patchInternalField());
00081 gradient() = vector::zero;
00082 }
00083
00084
00085 tractionDisplacementCorrectionFvPatchVectorField::
00086 tractionDisplacementCorrectionFvPatchVectorField
00087 (
00088 const tractionDisplacementCorrectionFvPatchVectorField& tdpvf
00089 )
00090 :
00091 fixedGradientFvPatchVectorField(tdpvf),
00092 traction_(tdpvf.traction_),
00093 pressure_(tdpvf.pressure_)
00094 {}
00095
00096
00097 tractionDisplacementCorrectionFvPatchVectorField::
00098 tractionDisplacementCorrectionFvPatchVectorField
00099 (
00100 const tractionDisplacementCorrectionFvPatchVectorField& tdpvf,
00101 const DimensionedField<vector, volMesh>& iF
00102 )
00103 :
00104 fixedGradientFvPatchVectorField(tdpvf, iF),
00105 traction_(tdpvf.traction_),
00106 pressure_(tdpvf.pressure_)
00107 {}
00108
00109
00110
00111
00112 void tractionDisplacementCorrectionFvPatchVectorField::autoMap
00113 (
00114 const fvPatchFieldMapper& m
00115 )
00116 {
00117 fixedGradientFvPatchVectorField::autoMap(m);
00118 traction_.autoMap(m);
00119 pressure_.autoMap(m);
00120 }
00121
00122
00123
00124 void tractionDisplacementCorrectionFvPatchVectorField::rmap
00125 (
00126 const fvPatchVectorField& ptf,
00127 const labelList& addr
00128 )
00129 {
00130 fixedGradientFvPatchVectorField::rmap(ptf, addr);
00131
00132 const tractionDisplacementCorrectionFvPatchVectorField& dmptf =
00133 refCast<const tractionDisplacementCorrectionFvPatchVectorField>(ptf);
00134
00135 traction_.rmap(dmptf.traction_, addr);
00136 pressure_.rmap(dmptf.pressure_, addr);
00137 }
00138
00139
00140
00141 void tractionDisplacementCorrectionFvPatchVectorField::updateCoeffs()
00142 {
00143 if (updated())
00144 {
00145 return;
00146 }
00147
00148 const dictionary& mechanicalProperties = db().lookupObject<IOdictionary>
00149 (
00150 "mechanicalProperties"
00151 );
00152
00153 dimensionedScalar rho(mechanicalProperties.lookup("rho"));
00154 dimensionedScalar rhoE(mechanicalProperties.lookup("E"));
00155 dimensionedScalar nu(mechanicalProperties.lookup("nu"));
00156
00157 dimensionedScalar E = rhoE/rho;
00158 dimensionedScalar mu = E/(2.0*(1.0 + nu));
00159 dimensionedScalar lambda = nu*E/((1.0 + nu)*(1.0 - 2.0*nu));
00160
00161 Switch planeStress(mechanicalProperties.lookup("planeStress"));
00162
00163 if (planeStress)
00164 {
00165 lambda = nu*E/((1.0 + nu)*(1.0 - nu));
00166 }
00167
00168 vectorField n = patch().nf();
00169
00170 const fvPatchField<symmTensor>& sigmaD =
00171 patch().lookupPatchField<volSymmTensorField, symmTensor>("sigmaD");
00172
00173 const fvPatchField<tensor>& sigmaExp =
00174 patch().lookupPatchField<volTensorField, tensor>("sigmaExp");
00175
00176 gradient() =
00177 (
00178 (traction_ + pressure_*n)/rho.value() - (n & (sigmaD + sigmaExp))
00179 )/(2.0*mu + lambda).value();
00180
00181 fixedGradientFvPatchVectorField::updateCoeffs();
00182 }
00183
00184
00185
00186 void tractionDisplacementCorrectionFvPatchVectorField::write(Ostream& os) const
00187 {
00188 fvPatchVectorField::write(os);
00189 traction_.writeEntry("traction", os);
00190 pressure_.writeEntry("pressure", os);
00191 writeEntry("value", os);
00192 }
00193
00194
00195
00196
00197 makePatchTypeField
00198 (
00199 fvPatchVectorField,
00200 tractionDisplacementCorrectionFvPatchVectorField
00201 );
00202
00203
00204
00205 }
00206
00207