tractionDisplacementFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const fvPatch& p,
36  const DimensionedField<vector, volMesh>& iF
37 )
38 :
39  fixedGradientFvPatchVectorField(p, iF),
40  traction_(p.size(), Zero),
41  pressure_(p.size(), 0.0)
42 {
43  fvPatchVectorField::operator=(patchInternalField());
44  gradient() = Zero;
45 }
46 
47 
50 (
52  const fvPatch& p,
53  const DimensionedField<vector, volMesh>& iF,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  fixedGradientFvPatchVectorField(tdpvf, p, iF, mapper),
58  traction_(mapper(tdpvf.traction_)),
59  pressure_(mapper(tdpvf.pressure_))
60 {}
61 
62 
65 (
66  const fvPatch& p,
67  const DimensionedField<vector, volMesh>& iF,
68  const dictionary& dict
69 )
70 :
71  fixedGradientFvPatchVectorField(p, iF),
72  traction_("traction", dict, p.size()),
73  pressure_("pressure", dict, p.size())
74 {
75  fvPatchVectorField::operator=(patchInternalField());
76  gradient() = Zero;
77 }
78 
79 
82 (
84 )
85 :
86  fixedGradientFvPatchVectorField(tdpvf),
87  traction_(tdpvf.traction_),
88  pressure_(tdpvf.pressure_)
89 {}
90 
91 
94 (
96  const DimensionedField<vector, volMesh>& iF
97 )
98 :
99  fixedGradientFvPatchVectorField(tdpvf, iF),
100  traction_(tdpvf.traction_),
101  pressure_(tdpvf.pressure_)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 (
109  const fvPatchFieldMapper& m
110 )
111 {
112  fixedGradientFvPatchVectorField::autoMap(m);
113  m(traction_, traction_);
114  m(pressure_, pressure_);
115 }
116 
117 
119 (
120  const fvPatchVectorField& ptf,
121  const labelList& addr
122 )
123 {
124  fixedGradientFvPatchVectorField::rmap(ptf, addr);
125 
127  refCast<const tractionDisplacementFvPatchVectorField>(ptf);
128 
129  traction_.rmap(dmptf.traction_, addr);
130  pressure_.rmap(dmptf.pressure_, addr);
131 }
132 
133 
135 {
136  if (updated())
137  {
138  return;
139  }
140 
141  const dictionary& mechanicalProperties =
142  db().lookupObject<IOdictionary>("mechanicalProperties");
143 
144  const dictionary& thermalProperties =
145  db().lookupObject<IOdictionary>("thermalProperties");
146 
147 
148  const fvPatchField<scalar>& rho =
149  patch().lookupPatchField<volScalarField, scalar>("rho");
150 
151  const fvPatchField<scalar>& rhoE =
152  patch().lookupPatchField<volScalarField, scalar>("E");
153 
154  const fvPatchField<scalar>& nu =
155  patch().lookupPatchField<volScalarField, scalar>("nu");
156 
157  scalarField E(rhoE/rho);
158  scalarField mu(E/(2.0*(1.0 + nu)));
159  scalarField lambda(nu*E/((1.0 + nu)*(1.0 - 2.0*nu)));
160  scalarField threeK(E/(1.0 - 2.0*nu));
161 
162  Switch planeStress(mechanicalProperties.lookup("planeStress"));
163 
164  if (planeStress)
165  {
166  lambda = nu*E/((1.0 + nu)*(1.0 - nu));
167  threeK = E/(1.0 - nu);
168  }
169 
170  scalarField twoMuLambda(2*mu + lambda);
171 
172  vectorField n(patch().nf());
173 
174  const fvPatchField<symmTensor>& sigmaD =
175  patch().lookupPatchField<volSymmTensorField, symmTensor>("sigmaD");
176 
177  gradient() =
178  (
179  (traction_ - pressure_*n)/rho
180  + twoMuLambda*fvPatchField<vector>::snGrad() - (n & sigmaD)
181  )/twoMuLambda;
182 
183  Switch thermalStress(thermalProperties.lookup("thermalStress"));
184 
185  if (thermalStress)
186  {
187  const fvPatchField<scalar>& threeKalpha=
188  patch().lookupPatchField<volScalarField, scalar>("threeKalpha");
189 
190  const fvPatchField<scalar>& T =
191  patch().lookupPatchField<volScalarField, scalar>("T");
192 
193  gradient() += n*threeKalpha*T/twoMuLambda;
194  }
195 
196  fixedGradientFvPatchVectorField::updateCoeffs();
197 }
198 
199 
201 {
203  writeEntry(os, "traction", traction_);
204  writeEntry(os, "pressure", pressure_);
205  writeEntry(os, "value", *this);
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 namespace Foam
212 {
214  (
217  );
218 }
219 
220 
221 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
fvPatchField< vector > fvPatchVectorField
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:186
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
tractionDisplacementFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionedScalar mu
Atomic mass unit.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:295
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
volScalarField & nu
Namespace for OpenFOAM.