tractionDisplacementFvPatchVectorField.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-2016 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 "volFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
41  const DimensionedField<vector, volMesh>& iF
42 )
43 :
44  fixedGradientFvPatchVectorField(p, iF),
45  traction_(p.size(), Zero),
46  pressure_(p.size(), 0.0)
47 {
48  fvPatchVectorField::operator=(patchInternalField());
49  gradient() = Zero;
50 }
51 
52 
53 tractionDisplacementFvPatchVectorField::
54 tractionDisplacementFvPatchVectorField
55 (
56  const tractionDisplacementFvPatchVectorField& tdpvf,
57  const fvPatch& p,
58  const DimensionedField<vector, volMesh>& iF,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  fixedGradientFvPatchVectorField(tdpvf, p, iF, mapper),
63  traction_(tdpvf.traction_, mapper),
64  pressure_(tdpvf.pressure_, mapper)
65 {}
66 
67 
68 tractionDisplacementFvPatchVectorField::
69 tractionDisplacementFvPatchVectorField
70 (
71  const fvPatch& p,
72  const DimensionedField<vector, volMesh>& iF,
73  const dictionary& dict
74 )
75 :
76  fixedGradientFvPatchVectorField(p, iF),
77  traction_("traction", dict, p.size()),
78  pressure_("pressure", dict, p.size())
79 {
80  fvPatchVectorField::operator=(patchInternalField());
81  gradient() = Zero;
82 }
83 
84 
85 tractionDisplacementFvPatchVectorField::
86 tractionDisplacementFvPatchVectorField
87 (
88  const tractionDisplacementFvPatchVectorField& tdpvf
89 )
90 :
91  fixedGradientFvPatchVectorField(tdpvf),
92  traction_(tdpvf.traction_),
93  pressure_(tdpvf.pressure_)
94 {}
95 
96 
97 tractionDisplacementFvPatchVectorField::
98 tractionDisplacementFvPatchVectorField
99 (
100  const tractionDisplacementFvPatchVectorField& tdpvf,
101  const DimensionedField<vector, volMesh>& iF
102 )
103 :
104  fixedGradientFvPatchVectorField(tdpvf, iF),
105  traction_(tdpvf.traction_),
106  pressure_(tdpvf.pressure_)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 void tractionDisplacementFvPatchVectorField::autoMap
113 (
114  const fvPatchFieldMapper& m
115 )
116 {
117  fixedGradientFvPatchVectorField::autoMap(m);
118  traction_.autoMap(m);
119  pressure_.autoMap(m);
120 }
121 
122 
123 void tractionDisplacementFvPatchVectorField::rmap
124 (
125  const fvPatchVectorField& ptf,
126  const labelList& addr
127 )
128 {
129  fixedGradientFvPatchVectorField::rmap(ptf, addr);
130 
131  const tractionDisplacementFvPatchVectorField& dmptf =
132  refCast<const tractionDisplacementFvPatchVectorField>(ptf);
133 
134  traction_.rmap(dmptf.traction_, addr);
135  pressure_.rmap(dmptf.pressure_, addr);
136 }
137 
138 
139 void tractionDisplacementFvPatchVectorField::updateCoeffs()
140 {
141  if (updated())
142  {
143  return;
144  }
145 
146  const dictionary& mechanicalProperties =
147  db().lookupObject<IOdictionary>("mechanicalProperties");
148 
149  const dictionary& thermalProperties =
150  db().lookupObject<IOdictionary>("thermalProperties");
151 
152 
153  const fvPatchField<scalar>& rho =
154  patch().lookupPatchField<volScalarField, scalar>("rho");
155 
156  const fvPatchField<scalar>& rhoE =
157  patch().lookupPatchField<volScalarField, scalar>("E");
158 
159  const fvPatchField<scalar>& nu =
160  patch().lookupPatchField<volScalarField, scalar>("nu");
161 
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));
166 
167  Switch planeStress(mechanicalProperties.lookup("planeStress"));
168 
169  if (planeStress)
170  {
171  lambda = nu*E/((1.0 + nu)*(1.0 - nu));
172  threeK = E/(1.0 - nu);
173  }
174 
175  scalarField twoMuLambda(2*mu + lambda);
176 
177  vectorField n(patch().nf());
178 
179  const fvPatchField<symmTensor>& sigmaD =
180  patch().lookupPatchField<volSymmTensorField, symmTensor>("sigmaD");
181 
182  gradient() =
183  (
184  (traction_ - pressure_*n)/rho
185  + twoMuLambda*fvPatchField<vector>::snGrad() - (n & sigmaD)
186  )/twoMuLambda;
187 
188  Switch thermalStress(thermalProperties.lookup("thermalStress"));
189 
190  if (thermalStress)
191  {
192  const fvPatchField<scalar>& threeKalpha=
193  patch().lookupPatchField<volScalarField, scalar>("threeKalpha");
194 
195  const fvPatchField<scalar>& T =
196  patch().lookupPatchField<volScalarField, scalar>("T");
197 
198  gradient() += n*threeKalpha*T/twoMuLambda;
199  }
200 
201  fixedGradientFvPatchVectorField::updateCoeffs();
202 }
203 
204 
205 void tractionDisplacementFvPatchVectorField::write(Ostream& os) const
206 {
208  traction_.writeEntry("traction", os);
209  pressure_.writeEntry("pressure", os);
210  writeEntry("value", os);
211 }
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
217 (
219  tractionDisplacementFvPatchVectorField
220 );
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace Foam
225 
226 // ************************************************************************* //
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:658
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
fvPatchField< vector > fvPatchVectorField
tractionDisplacementFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
volVectorField vectorField(fieldObject, mesh)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
volScalarField scalarField(fieldObject, mesh)
const volScalarField & T
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
label n
volScalarField & nu
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.