solidDisplacementThermo.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) 2019-2020 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 "fvMesh.H"
28 #include "fvmLaplacian.H"
29 #include "fvcSnGrad.H"
30 #include "surfaceInterpolate.H"
31 
32 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(solidDisplacementThermo, 0);
37 }
38 
39 
40 void Foam::solidDisplacementThermo::readProperty(volScalarField& prop) const
41 {
42  const dictionary& propDict(subDict(prop.name()));
43  const word propType(propDict.lookup("type"));
44 
45  if (propType == "uniform")
46  {
47  prop == dimensionedScalar
48  (
49  prop.name(),
50  prop.dimensions(),
51  propDict.lookup<scalar>("value")
52  );
53  }
54  else if (propType == "field")
55  {
56  const volScalarField propField
57  (
58  IOobject
59  (
60  prop.name(),
61  prop.mesh().time().timeName(0),
62  prop.mesh(),
65  ),
66  prop.mesh()
67  );
68  prop == propField;
69  }
70  else
71  {
73  << "Valid type entries are uniform or field for " << prop.name()
74  << abort(FatalError);
75  }
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 (
83  const fvMesh& mesh,
84  const word& phaseName
85 )
86 :
87  solidThermo::composite(mesh, phaseName),
88  planeStress_(lookup("planeStress")),
89  thermalStress_(lookup("thermalStress")),
90  Cp_
91  (
92  IOobject
93  (
94  phasePropertyName("Cp"),
95  mesh.time().timeName(),
96  mesh,
97  IOobject::NO_READ,
98  IOobject::NO_WRITE
99  ),
100  mesh,
102  ),
103  kappa_
104  (
105  IOobject
106  (
107  phasePropertyName("kappa"),
108  mesh.time().timeName(),
109  mesh,
110  IOobject::NO_READ,
111  IOobject::NO_WRITE
112  ),
113  mesh,
114  Cp_.dimensions()*dimensionSet(1, -1, -1, 0, 0)
115  ),
116  E_
117  (
118  IOobject
119  (
120  phasePropertyName("E"),
121  mesh.time().timeName(),
122  mesh,
123  IOobject::NO_READ,
124  IOobject::NO_WRITE
125  ),
126  mesh,
128  ),
129  nu_
130  (
131  IOobject
132  (
133  phasePropertyName("nu"),
134  mesh.time().timeName(),
135  mesh,
136  IOobject::NO_READ,
137  IOobject::NO_WRITE
138  ),
139  mesh,
140  dimless
141  ),
142  alphav_
143  (
144  IOobject
145  (
146  phasePropertyName("alphav"),
147  mesh.time().timeName(),
148  mesh,
149  IOobject::NO_READ,
150  IOobject::NO_WRITE
151  ),
152  mesh,
154  )
155 {
156  readProperty(rho_);
157  readProperty(Cp_);
158  readProperty(kappa_);
159  readProperty(E_);
160  readProperty(nu_);
161  readProperty(alphav_);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
166 
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  return rho_;
176 }
177 
178 
180 (
181  const label patchi
182 ) const
183 {
184  return rho_.boundaryField()[patchi];
185 }
186 
187 
189 {
190  return E_;
191 }
192 
193 
195 (
196  const label patchi
197 ) const
198 {
199  return E_.boundaryField()[patchi];
200 }
201 
202 
204 {
205  return nu_;
206 }
207 
208 
210 (
211  const label patchi
212 ) const
213 {
214  return nu_.boundaryField()[patchi];
215 }
216 
217 
219 {
220  return alphav_;
221 }
222 
223 
225 (
226  const label patchi
227 ) const
228 {
229  return alphav_.boundaryField()[patchi];
230 }
231 
232 
234 {
236  return rho_;
237 }
238 
239 
241 {
243  return rho_;
244 }
245 
246 
248 (
249  const scalarField& T,
250  const labelList& cells
251 ) const
252 {
254  return tmp<scalarField>(nullptr);
255 }
256 
258 (
259  const volScalarField& p,
260  const volScalarField& T
261 ) const
262 {
264  return tmp<volScalarField>(nullptr);
265 }
266 
267 
269 (
270  const scalarField& T,
271  const label patchi
272 ) const
273 {
275  return tmp<scalarField>(nullptr);
276 }
277 
278 
280 {
282  return tmp<volScalarField>(nullptr);
283 }
284 
285 
287 (
288  const volScalarField& p,
289  const volScalarField& T
290 ) const
291 {
293  return tmp<volScalarField>(nullptr);
294 }
295 
296 
298 (
299  const scalarField& T,
300  const label patchi
301 ) const
302 {
304  return tmp<scalarField>(nullptr);
305 }
306 
307 
309 (
310  const scalarField& T,
311  const labelList& cells
312 ) const
313 {
315  return tmp<scalarField>(nullptr);
316 }
317 
318 
320 {
322  return tmp<volScalarField>(nullptr);
323 }
324 
325 
327 (
328  const volScalarField& p,
329  const volScalarField& T
330 ) const
331 {
333  return tmp<volScalarField>(nullptr);
334 }
335 
336 
338 (
339  const scalarField& T,
340  const label patchi
341 ) const
342 {
344  return tmp<scalarField>(nullptr);
345 }
346 
347 
349 (
350  const scalarField& T,
351  const labelList& cells
352 ) const
353 {
355  return tmp<scalarField>(nullptr);
356 }
357 
358 
360 {
362  return tmp<volScalarField>(nullptr);
363 }
364 
365 
367 (
368  const volScalarField& h,
369  const volScalarField& p,
370  const volScalarField& T0
371 ) const
372 {
374  return tmp<volScalarField>(nullptr);
375 }
376 
377 
379 (
380  const scalarField& he,
381  const scalarField& T0,
382  const labelList& cells
383 ) const
384 {
386  return tmp<scalarField>(nullptr);
387 }
388 
389 
391 (
392  const scalarField& he,
393  const scalarField& T0,
394  const label patchi
395 ) const
396 {
398  return tmp<scalarField>(nullptr);
399 }
400 
401 
403 {
404  return Cp_;
405 }
406 
407 
409 (
410  const scalarField& T,
411  const label patchi
412 ) const
413 {
414  return Cp_.boundaryField()[patchi];
415 }
416 
417 
419 {
420  return Cp_;
421 }
422 
423 
425 (
426  const scalarField& T,
427  const label patchi
428 ) const
429 {
430  return Cp_.boundaryField()[patchi];
431 }
432 
433 
435 {
436  return Cp_;
437 }
438 
439 
441 (
442  const scalarField& T,
443  const label patchi
444 ) const
445 {
446  return Cp_.boundaryField()[patchi];
447 }
448 
449 
451 {
452  return kappa_;
453 }
454 
455 
457 (
458  const label patchi
459 ) const
460 {
461  return kappa_.boundaryField()[patchi];
462 }
463 
464 
466 {
468  return tmp<volScalarField>(nullptr);
469 }
470 
471 
473 (
474  const label patchi
475 ) const
476 {
478  return tmp<scalarField>(nullptr);
479 }
480 
481 
483 (
484  const volScalarField& alphat
485 ) const
486 {
488  return tmp<volScalarField>(nullptr);
489 }
490 
491 
493 (
494  const scalarField& alphat,
495  const label patchi
496 ) const
497 {
499  return tmp<scalarField>(nullptr);
500 }
501 
502 
504 (
505  const volScalarField& alphat
506 ) const
507 {
509  return tmp<volScalarField>(nullptr);
510 }
511 
512 
514 (
515  const scalarField& alphat,
516  const label patchi
517 ) const
518 {
520  return tmp<scalarField>(nullptr);
521 }
522 
523 
525 {
527  return tmp<volVectorField>(nullptr);
528 }
529 
530 
532 (
533  const label patchi
534 ) const
535 {
537  return tmp<vectorField>(nullptr);
538 }
539 
540 
542 (
543  const label patchi
544 ) const
545 {
547  return tmp<symmTensorField>(nullptr);
548 }
549 
550 
552 {}
553 
554 
556 {
557  return regIOobject::read();
558 }
559 
560 
563 {
564  return -fvc::interpolate(kappa_)*fvc::snGrad(T());
565 }
566 
567 
570 {
571  return -fvm::laplacian(kappa_, T);
572 }
573 
574 
575 // ************************************************************************* //
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
virtual bool read()
Read object.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual ~solidDisplacementThermo()
Destructor.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionSet dimPressure
Calculate the matrix for the laplacian of the field.
Calculate the snGrad of the given volField.
const dimensionSet dimless
virtual tmp< symmTensorField > KappaLocal(const label patchi) const
Anisotropic thermal conductivity for patch.
virtual tmp< surfaceScalarField > q() const
Return the heat flux.
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [W/m/K].
const dimensionedScalar h
Planck constant.
volScalarField rho_
Density field [kg/m^3].
Definition: solidThermo.H:173
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
solidDisplacementThermo(const fvMesh &, const word &phaseName=word::null)
Construct from mesh and phase name.
const cellShapeList & cells
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:449
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
virtual const volScalarField & alphav() const
Volumetric thermal expansion coefficient [1/T].
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
const dimensionSet dimEnergy
virtual const volScalarField & E() const
Youngs modulus [Pa].
const dimensionSet dimMass
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual const volScalarField & nu() const
Poisson&#39;s ratio [].
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
virtual void correct()
Update properties.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool read()
Read thermophysicalProperties dictionary.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal turbulent diffusivity for temperature.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
const dimensionSet dimTemperature
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22