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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
43  const DimensionedField<vector, volMesh>& iF
44 )
45 :
46  fixedGradientFvPatchVectorField(p, iF),
47  traction_(p.size(), Zero),
48  pressure_(p.size(), Zero)
49 {
50  fvPatchVectorField::operator=(patchInternalField());
51  gradient() = Zero;
52 }
53 
54 
55 tractionDisplacementFvPatchVectorField::
56 tractionDisplacementFvPatchVectorField
57 (
58  const tractionDisplacementFvPatchVectorField& tdpvf,
59  const fvPatch& p,
60  const DimensionedField<vector, volMesh>& iF,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedGradientFvPatchVectorField(tdpvf, p, iF, mapper),
65  traction_(tdpvf.traction_, mapper),
66  pressure_(tdpvf.pressure_, mapper)
67 {}
68 
69 
70 tractionDisplacementFvPatchVectorField::
71 tractionDisplacementFvPatchVectorField
72 (
73  const fvPatch& p,
74  const DimensionedField<vector, volMesh>& iF,
75  const dictionary& dict
76 )
77 :
78  fixedGradientFvPatchVectorField(p, iF),
79  traction_("traction", dict, p.size()),
80  pressure_("pressure", dict, p.size())
81 {
82  fvPatchVectorField::operator=(patchInternalField());
83  gradient() = Zero;
84 }
85 
86 
87 tractionDisplacementFvPatchVectorField::
88 tractionDisplacementFvPatchVectorField
89 (
90  const tractionDisplacementFvPatchVectorField& tdpvf
91 )
92 :
93  fixedGradientFvPatchVectorField(tdpvf),
94  traction_(tdpvf.traction_),
95  pressure_(tdpvf.pressure_)
96 {}
97 
98 
99 tractionDisplacementFvPatchVectorField::
100 tractionDisplacementFvPatchVectorField
101 (
102  const tractionDisplacementFvPatchVectorField& tdpvf,
103  const DimensionedField<vector, volMesh>& iF
104 )
105 :
106  fixedGradientFvPatchVectorField(tdpvf, iF),
107  traction_(tdpvf.traction_),
108  pressure_(tdpvf.pressure_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 void tractionDisplacementFvPatchVectorField::autoMap
115 (
116  const fvPatchFieldMapper& m
117 )
118 {
119  fixedGradientFvPatchVectorField::autoMap(m);
120  traction_.autoMap(m);
121  pressure_.autoMap(m);
122 }
123 
124 
125 void tractionDisplacementFvPatchVectorField::rmap
126 (
127  const fvPatchVectorField& ptf,
128  const labelList& addr
129 )
130 {
131  fixedGradientFvPatchVectorField::rmap(ptf, addr);
132 
133  const tractionDisplacementFvPatchVectorField& dmptf =
134  refCast<const tractionDisplacementFvPatchVectorField>(ptf);
135 
136  traction_.rmap(dmptf.traction_, addr);
137  pressure_.rmap(dmptf.pressure_, addr);
138 }
139 
140 
141 void tractionDisplacementFvPatchVectorField::updateCoeffs()
142 {
143  if (updated())
144  {
145  return;
146  }
147 
148  const dictionary& mechanicalProperties =
149  db().lookupObject<IOdictionary>("mechanicalProperties");
150 
151  const dictionary& thermalProperties =
152  db().lookupObject<IOdictionary>("thermalProperties");
153 
154 
155  const fvPatchField<scalar>& rho =
156  patch().lookupPatchField<volScalarField, scalar>("rho");
157 
158  const fvPatchField<scalar>& rhoE =
159  patch().lookupPatchField<volScalarField, scalar>("E");
160 
161  const fvPatchField<scalar>& nu =
162  patch().lookupPatchField<volScalarField, scalar>("nu");
163 
164  scalarField E(rhoE/rho);
165  scalarField mu(E/(2.0*(1.0 + nu)));
166  scalarField lambda(nu*E/((1.0 + nu)*(1.0 - 2.0*nu)));
167  scalarField threeK(E/(1.0 - 2.0*nu));
168 
169  if (mechanicalProperties.get<bool>("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  if (thermalProperties.get<bool>("thermalStress"))
189  {
190  const fvPatchField<scalar>& threeKalpha=
191  patch().lookupPatchField<volScalarField, scalar>("threeKalpha");
192 
193  const fvPatchField<scalar>& T =
194  patch().lookupPatchField<volScalarField, scalar>("T");
195 
196  gradient() += n*threeKalpha*T/twoMuLambda;
197  }
198 
199  fixedGradientFvPatchVectorField::updateCoeffs();
200 }
201 
202 
204 {
206  traction_.writeEntry("traction", os);
207  pressure_.writeEntry("pressure", os);
208  writeEntry("value", os);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
215 (
217  tractionDisplacementFvPatchVectorField
218 );
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 } // End namespace Foam
223 
224 // ************************************************************************* //
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:916
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:88
dictionary dict
fvPatchField< vector > fvPatchVectorField
volScalarField & rhoE
tractionDisplacementFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:55
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
OBJstream os(runTime.globalPath()/outputName)
const volScalarField & T
const dimensionedScalar mu
Atomic mass unit.
const std::string patch
OpenFOAM patch number as a std::string.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
volScalarField & nu
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157