solidDisplacement.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) 2023 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 
26 #include "solidDisplacement.H"
27 #include "fvcGrad.H"
28 #include "fvcDiv.H"
29 #include "fvcLaplacian.H"
30 #include "fvmD2dt2.H"
31 #include "fvmLaplacian.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace solvers
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
47 
49 {
51 }
52 
53 
55 {
56  solid::read();
57 
58  nCorr = pimple.dict().lookupOrDefault<int>("nCorrectors", 1);
59  convergenceTolerance = pimple.dict().lookupOrDefault<scalar>("D", 0);
60  pimple.dict().lookup("compactNormalStress") >> compactNormalStress;
61  accFac = pimple.dict().lookupOrDefault<scalar>("accelerationFactor", 1);
62 
63  return true;
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 :
71  solid
72  (
73  mesh,
75  ),
76 
77  thermo_(refCast<solidDisplacementThermo>(solid::thermo_)),
78 
79  compactNormalStress(pimple.dict().lookup("compactNormalStress")),
80 
81  D_
82  (
83  IOobject
84  (
85  "D",
86  runTime.name(),
87  mesh,
88  IOobject::MUST_READ,
89  IOobject::AUTO_WRITE
90  ),
91  mesh
92  ),
93 
94  E(thermo_.E()),
95  nu(thermo_.nu()),
96 
97  mu(E/(2*(1 + nu))),
98 
99  lambda
100  (
101  thermo_.planeStress()
102  ? nu*E/((1 + nu)*(1 - nu))
103  : nu*E/((1 + nu)*(1 - 2*nu))
104  ),
105 
106  threeK
107  (
108  thermo_.planeStress()
109  ? E/(1 - nu)
110  : E/(1 - 2*nu)
111  ),
112 
113  threeKalpha("threeKalpha", threeK*thermo_.alphav()),
114 
115  sigmaD
116  (
117  IOobject
118  (
119  "sigmaD",
120  runTime.name(),
121  mesh
122  ),
123  mu*twoSymm(fvc::grad(D_)) + lambda*(I*tr(fvc::grad(D_)))
124  ),
125 
126  divSigmaExp
127  (
128  IOobject
129  (
130  "divSigmaExp",
131  runTime.name(),
132  mesh
133  ),
134  fvc::div(sigmaD)
135  - (
136  compactNormalStress
137  ? fvc::laplacian(2*mu + lambda, D_, "laplacian(DD,D)")
138  : fvc::div((2*mu + lambda)*fvc::grad(D_), "div(sigmaD)")
139  )
140  ),
141 
142  thermo(thermo_),
143  D(D_)
144 {
146 
147  // Read the controls
148  read();
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153 
155 {}
156 
157 
158 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
159 
161 {
162  if (thermo.thermalStress())
163  {
165  }
166 }
167 
168 
170 {
171  if (thermo.thermalStress())
172  {
174  }
175 }
176 
177 
179 {
180  volVectorField& D(D_);
181 
182  const volScalarField& rho = thermo_.rho();
183 
184  int iCorr = 0;
185  scalar initialResidual = 0;
186 
187  {
188  {
189  fvVectorMatrix DEqn
190  (
191  fvm::d2dt2(rho, D)
192  ==
193  fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
194  + divSigmaExp
195  + rho*fvModels().d2dt2(D)
196  );
197 
198  if (thermo.thermalStress())
199  {
200  DEqn += fvc::grad(threeKalpha*T);
201  }
202 
203  fvConstraints().constrain(DEqn);
204 
205  initialResidual = DEqn.solve().max().initialResidual();
206 
207  // For steady-state optionally accelerate the solution
208  // by over-relaxing the displacement
209  if (mesh.schemes().steady() && accFac > 1)
210  {
211  D += (accFac - 1)*(D - D.oldTime());
212  }
213 
214  if (!compactNormalStress)
215  {
216  divSigmaExp = fvc::div(DEqn.flux());
217  }
218  }
219 
220  const volTensorField gradD(fvc::grad(D));
221  sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD);
222 
223  if (compactNormalStress)
224  {
225  divSigmaExp = fvc::div
226  (
227  sigmaD - (2*mu + lambda)*gradD,
228  "div(sigmaD)"
229  );
230  }
231  else
232  {
233  divSigmaExp += fvc::div(sigmaD);
234  }
235 
236  } while (initialResidual > convergenceTolerance && ++iCorr < nCorr);
237 }
238 
239 
241 {
242  if (thermo.thermalStress())
243  {
245  }
246 }
247 
248 
250 {
251  if (runTime.writeTime())
252  {
254  (
255  IOobject
256  (
257  "sigma",
258  runTime.name(),
259  mesh
260  ),
261  sigmaD
262  );
263 
264  if (thermo.thermalStress())
265  {
266  sigma = sigma - I*(threeKalpha*thermo.T());
267  }
268 
269  volScalarField sigmaEq
270  (
271  IOobject
272  (
273  "sigmaEq",
274  runTime.name(),
275  mesh
276  ),
277  sqrt((3.0/2.0)*magSqr(dev(sigma)))
278  );
279 
280  Info<< "Max sigmaEq = " << max(sigmaEq).value()
281  << endl;
282 
283  sigma.write();
284  sigmaEq.write();
285  }
286 }
287 
288 
289 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:964
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1762
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
virtual bool modified() const
Return true if the object's file (or files for objectRegistry)
virtual bool write(const bool write=true) const
Write using setting from DB.
Fundamental solid thermodynamic properties.
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
Solver module for steady or transient segregated finite-volume solution of linear-elastic,...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
virtual ~solidDisplacement()
Destructor.
virtual void prePredictor()
Called at the beginning of the PIMPLE loop.
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
virtual bool dependenciesModified() const
Return true if the solver's dependencies have been modified.
solidDisplacement(fvMesh &mesh)
Construct from region mesh.
const volVectorField & D
Reference to the Displacement field.
virtual void pressureCorrector()
Construct and solve the displacement equation to obtain the stress.
virtual void postCorrector()
Correct the thermophysical transport modelling.
virtual bool read()
Read controls.
Solver module for thermal transport in solid domains and regions for conjugate heat transfer,...
Definition: solid.H:56
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
Definition: solid.C:218
virtual void prePredictor()
Called at the beginning of the PIMPLE loop.
Definition: solid.C:205
virtual bool dependenciesModified() const
Return true if the solver's dependencies have been modified.
Definition: solid.C:79
virtual void postCorrector()
Correct the thermophysical transport modelling.
Definition: solid.C:250
virtual bool read()
Read controls.
Definition: solid.C:85
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the laplacian of the given field.
Calculate the matrix for the second-order temporal derivative.
Calculate the matrix for the laplacian of the field.
dimensionedScalar lambda(viscosity->lookup("lambda"))
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< VolField< Type > > d2dt2(const VolField< Type > &vf)
Definition: fvcD2dt2.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &vf)
Definition: fvmD2dt2.C:46
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
messageStream Info
static const Identity< scalar > I
Definition: Identity.H:93
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31