1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "tractionDisplacementFvPatchVectorField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "volFields.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 tractionDisplacementFvPatchVectorField::
38 tractionDisplacementFvPatchVectorField
41 const DimensionedField<vector, volMesh>& iF
44 fixedGradientFvPatchVectorField(p, iF),
45 traction_(p.size(), vector::zero),
46 pressure_(p.size(), 0.0)
48 fvPatchVectorField::operator=(patchInternalField());
49 gradient() = vector::zero;
53 tractionDisplacementFvPatchVectorField::
54 tractionDisplacementFvPatchVectorField
56 const tractionDisplacementFvPatchVectorField& tdpvf,
58 const DimensionedField<vector, volMesh>& iF,
59 const fvPatchFieldMapper& mapper
62 fixedGradientFvPatchVectorField(tdpvf, p, iF, mapper),
63 traction_(tdpvf.traction_, mapper),
64 pressure_(tdpvf.pressure_, mapper)
68 tractionDisplacementFvPatchVectorField::
69 tractionDisplacementFvPatchVectorField
72 const DimensionedField<vector, volMesh>& iF,
73 const dictionary& dict
76 fixedGradientFvPatchVectorField(p, iF),
77 traction_("traction", dict, p.size()),
78 pressure_("pressure", dict, p.size())
80 fvPatchVectorField::operator=(patchInternalField());
81 gradient() = vector::zero;
85 tractionDisplacementFvPatchVectorField::
86 tractionDisplacementFvPatchVectorField
88 const tractionDisplacementFvPatchVectorField& tdpvf
91 fixedGradientFvPatchVectorField(tdpvf),
92 traction_(tdpvf.traction_),
93 pressure_(tdpvf.pressure_)
97 tractionDisplacementFvPatchVectorField::
98 tractionDisplacementFvPatchVectorField
100 const tractionDisplacementFvPatchVectorField& tdpvf,
101 const DimensionedField<vector, volMesh>& iF
104 fixedGradientFvPatchVectorField(tdpvf, iF),
105 traction_(tdpvf.traction_),
106 pressure_(tdpvf.pressure_)
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 void tractionDisplacementFvPatchVectorField::autoMap
114 const fvPatchFieldMapper& m
117 fixedGradientFvPatchVectorField::autoMap(m);
118 traction_.autoMap(m);
119 pressure_.autoMap(m);
123 void tractionDisplacementFvPatchVectorField::rmap
125 const fvPatchVectorField& ptf,
126 const labelList& addr
129 fixedGradientFvPatchVectorField::rmap(ptf, addr);
131 const tractionDisplacementFvPatchVectorField& dmptf =
132 refCast<const tractionDisplacementFvPatchVectorField>(ptf);
134 traction_.rmap(dmptf.traction_, addr);
135 pressure_.rmap(dmptf.pressure_, addr);
139 void tractionDisplacementFvPatchVectorField::updateCoeffs()
146 const dictionary& mechanicalProperties =
147 db().lookupObject<IOdictionary>("mechanicalProperties");
149 const dictionary& thermalProperties =
150 db().lookupObject<IOdictionary>("thermalProperties");
153 const fvPatchField<scalar>& rho =
154 patch().lookupPatchField<volScalarField, scalar>("rho");
156 const fvPatchField<scalar>& rhoE =
157 patch().lookupPatchField<volScalarField, scalar>("E");
159 const fvPatchField<scalar>& nu =
160 patch().lookupPatchField<volScalarField, scalar>("nu");
162 scalarField E = rhoE/rho;
163 scalarField mu = E/(2.0*(1.0 + nu));
164 scalarField lambda = nu*E/((1.0 + nu)*(1.0 - 2.0*nu));
165 scalarField threeK = E/(1.0 - 2.0*nu);
167 Switch planeStress(mechanicalProperties.lookup("planeStress"));
171 lambda = nu*E/((1.0 + nu)*(1.0 - nu));
172 threeK = E/(1.0 - nu);
175 scalarField twoMuLambda = (2*mu + lambda);
177 vectorField n(patch().nf());
179 const fvPatchField<symmTensor>& sigmaD =
180 patch().lookupPatchField<volSymmTensorField, symmTensor>("sigmaD");
184 (traction_ + pressure_*n)/rho
185 + twoMuLambda*fvPatchField<vector>::snGrad() - (n & sigmaD)
188 Switch thermalStress(thermalProperties.lookup("thermalStress"));
192 const fvPatchField<scalar>& threeKalpha=
193 patch().lookupPatchField<volScalarField, scalar>("threeKalpha");
195 const fvPatchField<scalar>& T =
196 patch().lookupPatchField<volScalarField, scalar>("T");
198 gradient() += n*threeKalpha*T/twoMuLambda;
201 fixedGradientFvPatchVectorField::updateCoeffs();
205 void tractionDisplacementFvPatchVectorField::write(Ostream& os) const
207 fvPatchVectorField::write(os);
208 traction_.writeEntry("traction", os);
209 pressure_.writeEntry("pressure", os);
210 writeEntry("value", os);
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 tractionDisplacementFvPatchVectorField
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 } // End namespace Foam
226 // ************************************************************************* //