anisotropic.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) 2022 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 "anisotropic.H"
27 #include "fvmLaplacian.H"
28 #include "fvcLaplacian.H"
29 #include "fvcSnGrad.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace solidThermophysicalTransportModels
37 {
40  (
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::solidThermophysicalTransportModels::anisotropic::
52 setZonesPatchFaces() const
53 {
54  if (!zoneCoordinateSystems_.size()) return;
55 
56  // Find all the patch faces adjacent to zones
57 
58  const fvMesh& mesh = thermo().mesh();
59  const fvBoundaryMesh& patches = mesh.boundary();
60  const labelList& own = mesh.faceOwner();
61 
62  zonesPatchFaces_.setSize(zoneCoordinateSystems_.size());
63  label zonei = 0;
64 
66  (
67  PtrDictionary<coordinateSystem>,
68  zoneCoordinateSystems_,
69  iter
70  )
71  {
72  zonesPatchFaces_[zonei].setSize(patches.size());
73 
74  const labelList& zoneCells = mesh.cellZones()[iter().name()];
75 
76  // Cell in zone
77  boolList cellInZone(mesh.nCells(), false);
78 
79  forAll(zoneCells, i)
80  {
81  cellInZone[zoneCells[i]] = true;
82  }
83 
85  {
86  const fvPatch& pp = patches[patchi];
87 
88  zonesPatchFaces_[zonei][patchi].setSize(pp.size());
89 
90  label nZonePatchFaces = 0;
91 
92  forAll(pp, patchFacei)
93  {
94  const label facei = pp.start() + patchFacei;
95 
96  if (cellInZone[own[facei]])
97  {
98  zonesPatchFaces_[zonei][patchi][nZonePatchFaces++] =
99  patchFacei;
100  }
101  }
102 
103  zonesPatchFaces_[zonei][patchi].setSize(nZonePatchFaces);
104  }
105 
106  zonei++;
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
114 (
115  const solidThermo& thermo
116 )
117 :
119  UpdateableMeshObject(*this, thermo.mesh()),
120  coordinateSystem_(coordinateSystem::New(thermo.mesh(), coeffDict())),
121  boundaryAligned_
122  (
123  coeffDict().lookupOrDefault<Switch>("boundaryAligned", false)
124  ),
125  aligned_(thermo.mesh().boundary().size(), true)
126 {
127  if (coeffDict().found("zones"))
128  {
129  const dictionary& zonesDict(coeffDict().subDict("zones"));
130 
131  Info<< " Reading coordinate system for zones:" << endl;
132 
133  forAllConstIter(dictionary, zonesDict, iter)
134  {
135  if (iter().isDict())
136  {
137  const word& zoneName = iter().keyword();
138  const dictionary& dict = iter().dict();
139 
140  Info<< " " << zoneName << endl;
141 
142  zoneCoordinateSystems_.insert
143  (
144  zoneName,
145  coordinateSystem::New(zoneName, dict).ptr()
146  );
147  }
148  }
149 
150  // Find all the patch faces adjacent to zones
151  setZonesPatchFaces();
152  }
153 
154  const fvMesh& mesh = thermo.mesh();
155  const fvBoundaryMesh& bMesh = mesh.boundary();
156 
157  bool aligned = true;
158 
160  {
161  if
162  (
163  !bMesh[patchi].coupled()
164  && returnReduce(bMesh[patchi].size(), sumOp<label>())
165  )
166  {
167  const vectorField n(bMesh[patchi].nf());
168  const vectorField nKappa(n & Kappa(patchi));
169 
170  // Calculate a patch alignment factor for Kappa
171  const scalar alignmentFactor =
172  gSum
173  (
174  mesh.magSf().boundaryField()[patchi]
175  *mag(nKappa - n*(nKappa & n))
176  )/gSum(mesh.magSf().boundaryField()[patchi]*mag(nKappa));
177 
178  if (alignmentFactor > 1e-3)
179  {
180  aligned_[patchi] = false;
181  aligned = false;
182 
183  Info<< " Kappa is not aligned with patch "
184  << bMesh[patchi].name()
185  << " by an alignment factor of " << alignmentFactor
186  << " (0=aligned, 1=unaligned)" << nl
187  << " heat-flux correction will be applied." << endl;
188  }
189  }
190  }
191 
192  if (!aligned && boundaryAligned_)
193  {
194  aligned_ = true;
195  aligned = true;
196 
197  Info<<
198  " boundaryAligned is set true, "
199  "boundary alignment of kappa will be enforced." << endl;
200  }
201 
202  // If Kappa is not aligned with any patch enable grad(T) caching
203  // because the patch grad(T) will be required for the heat-flux correction
204  if (!aligned)
205  {
206  mesh.solution().enableCache("grad(T)");
207  }
208 
209  Info << endl;
210 
211  mesh.schemes().setFluxRequired(thermo.T().name());
212 
213  this->printCoeffs(typeName);
214 }
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 {
221 
222  return true;
223 }
224 
225 
227 Foam::solidThermophysicalTransportModels::anisotropic::Kappa() const
228 {
229  const solidThermo& thermo = this->thermo();
230  const fvMesh& mesh = thermo.mesh();
231 
232  const volVectorField& materialKappa = thermo.Kappa();
233 
235  (
237  (
238  "Kappa",
239  mesh,
240  dimensionedSymmTensor(materialKappa.dimensions(), Zero)
241  )
242  );
243  volSymmTensorField& Kappa = tKappa.ref();
244 
245  Kappa.primitiveFieldRef() =
246  coordinateSystem_.R(mesh.C()).transformDiagTensor(materialKappa);
247 
248  forAll(Kappa.boundaryField(), patchi)
249  {
250  Kappa.boundaryFieldRef()[patchi] =
251  coordinateSystem_.R(mesh.boundary()[patchi].Cf())
252  .transformDiagTensor(materialKappa.boundaryField()[patchi]);
253  }
254 
255  label zonei = 0;
256 
258  (
260  zoneCoordinateSystems_,
261  iter
262  )
263  {
264  const labelList& zoneCells = mesh.cellZones()[iter().name()];
265  const coordinateSystem& cs = iter();
266 
267  symmTensorField& KappaIf = Kappa;
268 
269  forAll(zoneCells, i)
270  {
271  const label celli = zoneCells[i];
272 
273  KappaIf[celli] =
275  (
276  mesh.C()[celli],
277  materialKappa[celli]
278  );
279  }
280 
281  forAll(zonesPatchFaces_[zonei], patchi)
282  {
283  symmTensorField& KappaPf = Kappa.boundaryFieldRef()[patchi];
284 
285  const vectorField& materialKappaPf =
286  materialKappa.boundaryField()[patchi];
287 
288  const vectorField& CPf = mesh.boundary()[patchi].Cf();
289 
290  forAll(zonesPatchFaces_[zonei][patchi], zonePatchFacei)
291  {
292  const label patchFacei =
293  zonesPatchFaces_[zonei][patchi][zonePatchFacei];
294 
295  KappaPf[patchFacei] =
297  (
298  CPf[patchFacei],
299  materialKappaPf[patchFacei]
300  );
301  }
302  }
303 
304  zonei++;
305  }
306 
307  return tKappa;
308 }
309 
310 
312 Foam::solidThermophysicalTransportModels::anisotropic::Kappa
313 (
314  const label patchi
315 ) const
316 {
317  const solidThermo& thermo = this->thermo();
318  const vectorField& CPf = thermo.mesh().boundary()[patchi].Cf();
319 
320  const vectorField& materialKappaPf(thermo.Kappa().boundaryField()[patchi]);
321 
322  tmp<symmTensorField> tKappa
323  (
324  coordinateSystem_.R(CPf).transformDiagTensor(materialKappaPf)
325  );
326  symmTensorField& KappaPf = tKappa.ref();
327 
328  label zonei = 0;
329 
331  (
332  PtrDictionary<coordinateSystem>,
333  zoneCoordinateSystems_,
334  iter
335  )
336  {
337  const coordinateSystem& cs = iter();
338 
339  forAll(zonesPatchFaces_[zonei][patchi], zonePatchFacei)
340  {
341  const label patchFacei =
342  zonesPatchFaces_[zonei][patchi][zonePatchFacei];
343 
344  KappaPf[patchFacei] =
345  cs.R().transformDiagTensor
346  (
347  CPf[patchFacei],
348  materialKappaPf[patchFacei]
349  );
350  }
351 
352  zonei++;
353  }
354 
355  return tKappa;
356 }
357 
358 
361 {
363  return tmp<volScalarField>(nullptr);
364 }
365 
366 
369 (
370  const label patchi
371 ) const
372 {
373  const vectorField n(thermo().mesh().boundary()[patchi].nf());
374  return n & Kappa(patchi) & n;
375 }
376 
377 
380 {
381  const solidThermo& thermo = this->thermo();
382  const fvMesh& mesh = thermo.mesh();
383 
385  (
386  "q",
387  -fvm::laplacian(Kappa(), thermo.T())().flux()/mesh.magSf()
388  );
389 }
390 
391 
394 (
395  const label patchi
396 ) const
397 {
398  if (!aligned_[patchi])
399  {
401 
402  const vectorField n(thermo().mesh().boundary()[patchi].nf());
403  const vectorField nKappa(n & Kappa(patchi));
404 
405  return -(nKappa - n*(nKappa & n)) & gradT().boundaryField()[patchi];
406  }
407  else
408  {
409  return tmp<scalarField>(nullptr);
410  }
411 }
412 
413 
416 (
418 ) const
419 {
420  const solidThermo& thermo = this->thermo();
421  const volSymmTensorField Kappa(this->Kappa());
422  const surfaceVectorField& Sf = thermo.mesh().Sf();
423  const surfaceScalarField& magSf = thermo.mesh().magSf();
424 
425  // Return heat flux source as an implicit energy correction
426  // to the temperature gradient flux
427  return
428  -fvc::laplacian(Kappa, thermo.T())
430  (
431  (Sf & fvc::interpolate(Kappa/thermo.Cv()) & Sf)/sqr(magSf),
432  e
433  );
434 }
435 
436 
438 {
440 
441  // Recalculate zonesPatchFaces if they have been deleted
442  // following mesh changes
443  if (!zonesPatchFaces_.size())
444  {
445  setZonesPatchFaces();
446  }
447 }
448 
449 
451 {
452  return true;
453 }
454 
455 
457 (
458  const polyTopoChangeMap& map
459 )
460 {
461  // Delete the cached zonesPatchFaces, will be re-created in predict
462  zonesPatchFaces_.clear();
463 }
464 
465 
467 (
468  const polyMeshMap& map
469 )
470 {
471  // Delete the cached zonesPatchFaces, will be re-created in predict
472  zonesPatchFaces_.clear();
473 }
474 
475 
477 (
478  const polyDistributionMap& map
479 )
480 {
481  // Delete the cached zonesPatchFaces, will be re-created in predict
482  zonesPatchFaces_.clear();
483 }
484 
485 
486 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
void insert(const word &, T *)
Add at head of dictionary.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Internal & ref()
Return a reference to the dimensioned internal field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const word & name() const
Return name.
Definition: IOobject.H:310
void setSize(const label)
Reset size of List.
Definition: List.C:281
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const volScalarField & T() const =0
Temperature [K].
virtual symmTensor transformDiagTensor(const vector &p, const vector &v) const =0
Transform diagTensor masquerading as a vector using transformation.
Base class for other coordinate system specifications.
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
const coordinateRotation & R() const
Return const reference to co-ordinate rotation.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:952
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const volVectorField & C() const
Return cell centres.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
const fvSolution & solution() const
Return the fvSchemes.
Definition: fvMesh.C:1759
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:451
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:56
Abstract base class for solid thermophysical transport models.
virtual void printCoeffs(const word &type)
Print model coefficients.
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
virtual void predict()=0
Predict the thermophysical transport coefficients if possible.
virtual const solidThermo & thermo() const
Access function to solid thermophysical properties.
Solid thermophysical transport model for anisotropic thermal conductivity.
Definition: anisotropic.H:120
virtual bool movePoints()
Update for mesh motion.
Definition: anisotropic.C:450
virtual tmp< scalarField > qCorr(const label patchi) const
Return the patch heat flux correction [W/m^2].
Definition: anisotropic.C:394
virtual void predict()
Correct the anisotropic viscosity.
Definition: anisotropic.C:437
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
Definition: anisotropic.C:467
anisotropic(const solidThermo &thermo)
Construct from solid thermophysical properties.
Definition: anisotropic.C:114
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: anisotropic.C:477
virtual tmp< volScalarField > kappa() const
Thermal conductivity [W/m/K].
Definition: anisotropic.C:360
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: anisotropic.C:379
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
Definition: anisotropic.C:457
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: anisotropic.C:219
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: anisotropic.C:416
void enableCache(const word &name) const
Enable caching of the given field.
Definition: solution.C:234
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
Calculate the laplacian of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
const polyBoundaryMesh & bMesh
label patchi
const fvPatchList & patches
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
Definition: fvmLaplacian.C:340
addToRunTimeSelectionTable(solidThermophysicalTransportModel, anisotropic, dictionary)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
const doubleScalar e
Definition: doubleScalar.H:105
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
faceListList boundary(nPatches)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31