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-2024 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class SolidThermophysicalTransportModel>
37 setZonesPatchFaces() const
38 {
39  if (!zoneCoordinateSystems_.size()) return;
40 
41  // Find all the patch faces adjacent to zones
42 
43  const fvMesh& mesh = this->thermo().mesh();
44  const fvBoundaryMesh& patches = mesh.boundary();
45  const labelList& own = mesh.faceOwner();
46 
47  zonesPatchFaces_.setSize(zoneCoordinateSystems_.size());
48  label zonei = 0;
49 
51  (
52  PtrDictionary<coordinateSystem>,
53  zoneCoordinateSystems_,
54  iter
55  )
56  {
57  zonesPatchFaces_[zonei].setSize(patches.size());
58 
59  const labelList& zoneCells = mesh.cellZones()[iter().name()];
60 
61  // Cell in zone
62  boolList cellInZone(mesh.nCells(), false);
63 
64  forAll(zoneCells, i)
65  {
66  cellInZone[zoneCells[i]] = true;
67  }
68 
70  {
71  const fvPatch& pp = patches[patchi];
72 
73  zonesPatchFaces_[zonei][patchi].setSize(pp.size());
74 
75  label nZonePatchFaces = 0;
76 
77  forAll(pp, patchFacei)
78  {
79  const label facei = pp.start() + patchFacei;
80 
81  if (cellInZone[own[facei]])
82  {
83  zonesPatchFaces_[zonei][patchi][nZonePatchFaces++] =
84  patchFacei;
85  }
86  }
87 
88  zonesPatchFaces_[zonei][patchi].setSize(nZonePatchFaces);
89  }
90 
91  zonei++;
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
98 template<class SolidThermophysicalTransportModel>
101 (
102  const alphaField& alpha,
103  const solidThermo& thermo
104 )
105 :
106  SolidThermophysicalTransportModel(typeName, alpha, thermo),
108  coordinateSystem_(coordinateSystem::New(thermo.mesh(), this->coeffDict())),
109  boundaryAligned_
110  (
111  this->coeffDict().template lookupOrDefault<Switch>
112  (
113  "boundaryAligned",
114  false
115  )
116  ),
117  aligned_(thermo.mesh().boundary().size(), true)
118 {
119  if (this->coeffDict().found("zones"))
120  {
121  const dictionary& zonesDict(this->coeffDict().subDict("zones"));
122 
123  Info<< " Reading coordinate system for zones:" << endl;
124 
125  forAllConstIter(dictionary, zonesDict, iter)
126  {
127  if (iter().isDict())
128  {
129  const word& zoneName = iter().keyword();
130  const dictionary& dict = iter().dict();
131 
132  Info<< " " << zoneName << endl;
133 
134  zoneCoordinateSystems_.insert
135  (
136  zoneName,
137  coordinateSystem::New(zoneName, dict).ptr()
138  );
139  }
140  }
141 
142  // Find all the patch faces adjacent to zones
143  setZonesPatchFaces();
144  }
145 
146  const fvMesh& mesh = thermo.mesh();
147  const fvBoundaryMesh& bMesh = mesh.boundary();
148 
149  bool aligned = true;
150 
152  {
153  if
154  (
155  !bMesh[patchi].coupled()
156  && returnReduce(bMesh[patchi].size(), sumOp<label>())
157  )
158  {
159  const vectorField n(bMesh[patchi].nf());
160  const vectorField nKappa(n & Kappa(patchi));
161 
162  // Calculate a patch alignment factor for Kappa
163  const scalar alignmentFactor =
164  gSum
165  (
166  mesh.magSf().boundaryField()[patchi]
167  *mag(nKappa - n*(nKappa & n))
168  )/gSum(mesh.magSf().boundaryField()[patchi]*mag(nKappa));
169 
170  if (alignmentFactor > 1e-3)
171  {
172  aligned_[patchi] = false;
173  aligned = false;
174 
175  Info<< " Kappa is not aligned with patch "
176  << bMesh[patchi].name()
177  << " by an alignment factor of " << alignmentFactor
178  << " (0=aligned, 1=unaligned)" << nl
179  << " heat-flux correction will be applied." << endl;
180  }
181  }
182  }
183 
184  if (!aligned && boundaryAligned_)
185  {
186  aligned_ = true;
187  aligned = true;
188 
189  Info<<
190  " boundaryAligned is set true, "
191  "boundary alignment of kappa will be enforced." << endl;
192  }
193 
194  // If Kappa is not aligned with any patch enable grad(T) caching
195  // because the patch grad(T) will be required for the heat-flux correction
196  if (!aligned)
197  {
198  mesh.solution().enableCache("grad(T)");
199  }
200 
201  Info << endl;
202 
203  mesh.schemes().setFluxRequired(thermo.T().name());
204 
205  this->printCoeffs(typeName);
206 }
207 
208 
209 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210 
211 template<class SolidThermophysicalTransportModel>
214 {
215 
216  return true;
217 }
218 
219 
220 template<class SolidThermophysicalTransportModel>
224 {
225  const solidThermo& thermo = this->thermo();
226  const fvMesh& mesh = thermo.mesh();
227 
228  const volVectorField& materialKappa = thermo.Kappa();
229 
231  (
233  (
234  "Kappa",
235  mesh,
236  dimensionedSymmTensor(materialKappa.dimensions(), Zero)
237  )
238  );
239  volSymmTensorField& Kappa = tKappa.ref();
240 
241  Kappa.primitiveFieldRef() =
242  coordinateSystem_.R(mesh.C()).transformDiagTensor(materialKappa);
243 
244  forAll(Kappa.boundaryField(), patchi)
245  {
246  Kappa.boundaryFieldRef()[patchi] =
247  coordinateSystem_.R(mesh.boundary()[patchi].Cf())
248  .transformDiagTensor(materialKappa.boundaryField()[patchi]);
249  }
250 
251  label zonei = 0;
252 
254  (
256  zoneCoordinateSystems_,
257  iter
258  )
259  {
260  const labelList& zoneCells = mesh.cellZones()[iter().name()];
261  const coordinateSystem& cs = iter();
262 
263  symmTensorField& KappaIf = Kappa;
264 
265  forAll(zoneCells, i)
266  {
267  const label celli = zoneCells[i];
268 
269  KappaIf[celli] =
271  (
272  mesh.C()[celli],
273  materialKappa[celli]
274  );
275  }
276 
277  forAll(zonesPatchFaces_[zonei], patchi)
278  {
279  symmTensorField& KappaPf = Kappa.boundaryFieldRef()[patchi];
280 
281  const vectorField& materialKappaPf =
282  materialKappa.boundaryField()[patchi];
283 
284  const vectorField& CPf = mesh.boundary()[patchi].Cf();
285 
286  forAll(zonesPatchFaces_[zonei][patchi], zonePatchFacei)
287  {
288  const label patchFacei =
289  zonesPatchFaces_[zonei][patchi][zonePatchFacei];
290 
291  KappaPf[patchFacei] =
293  (
294  CPf[patchFacei],
295  materialKappaPf[patchFacei]
296  );
297  }
298  }
299 
300  zonei++;
301  }
302 
303  return tKappa;
304 }
305 
306 
307 template<class SolidThermophysicalTransportModel>
311 (
312  const label patchi
313 ) const
314 {
315  const solidThermo& thermo = this->thermo();
316  const vectorField& CPf = thermo.mesh().boundary()[patchi].Cf();
317 
318  const vectorField& materialKappaPf(thermo.Kappa().boundaryField()[patchi]);
319 
320  tmp<symmTensorField> tKappa
321  (
322  coordinateSystem_.R(CPf).transformDiagTensor(materialKappaPf)
323  );
324  symmTensorField& KappaPf = tKappa.ref();
325 
326  label zonei = 0;
327 
329  (
330  PtrDictionary<coordinateSystem>,
331  zoneCoordinateSystems_,
332  iter
333  )
334  {
335  const coordinateSystem& cs = iter();
336 
337  forAll(zonesPatchFaces_[zonei][patchi], zonePatchFacei)
338  {
339  const label patchFacei =
340  zonesPatchFaces_[zonei][patchi][zonePatchFacei];
341 
342  KappaPf[patchFacei] =
343  cs.R().transformDiagTensor
344  (
345  CPf[patchFacei],
346  materialKappaPf[patchFacei]
347  );
348  }
349 
350  zonei++;
351  }
352 
353  return tKappa;
354 }
355 
356 
357 template<class SolidThermophysicalTransportModel>
361 {
363  return tmp<volScalarField>(nullptr);
364 }
365 
366 
367 template<class SolidThermophysicalTransportModel>
371 (
372  const label patchi
373 ) const
374 {
375  const vectorField n(this->thermo().mesh().boundary()[patchi].nf());
376  return n & Kappa(patchi) & n;
377 }
378 
379 
380 template<class SolidThermophysicalTransportModel>
384 {
385  const solidThermo& thermo = this->thermo();
386  const fvMesh& mesh = thermo.mesh();
387 
389  (
390  "q",
391  -fvm::laplacian(this->alpha()*Kappa(), thermo.T())().flux()/mesh.magSf()
392  );
393 }
394 
395 
396 template<class SolidThermophysicalTransportModel>
400 (
401  const label patchi
402 ) const
403 {
404  if (!aligned_[patchi])
405  {
406  tmp<volVectorField> gradT(fvc::grad(this->thermo().T()));
407 
408  const vectorField n(this->thermo().mesh().boundary()[patchi].nf());
409  const vectorField nKappa(n & Kappa(patchi));
410 
411  return
412  -(
413  this->alpha().boundaryField()[patchi]
414  *((nKappa - n*(nKappa & n)) & gradT().boundaryField()[patchi])
415  );
416  }
417  else
418  {
419  return tmp<scalarField>(nullptr);
420  }
421 }
422 
423 
424 template<class SolidThermophysicalTransportModel>
428 (
430 ) const
431 {
432  const solidThermo& thermo = this->thermo();
433  const volSymmTensorField Kappa(this->Kappa());
434  const surfaceVectorField& Sf = thermo.mesh().Sf();
435  const surfaceScalarField& magSf = thermo.mesh().magSf();
436 
437  // Return heat flux source as an implicit energy correction
438  // to the temperature gradient flux
439  return
440  -fvc::laplacian(this->alpha()*Kappa, thermo.T())
442  (
443  (Sf & fvc::interpolate(this->alpha()*Kappa/thermo.Cv()) & Sf)
444  /sqr(magSf),
445  e
446  );
447 }
448 
449 
450 template<class SolidThermophysicalTransportModel>
453 {
454  SolidThermophysicalTransportModel::predict();
455 
456  // Recalculate zonesPatchFaces if they have been deleted
457  // following mesh changes
458  if (!zonesPatchFaces_.size())
459  {
460  setZonesPatchFaces();
461  }
462 }
463 
464 
465 template<class SolidThermophysicalTransportModel>
468 {
469  return true;
470 }
471 
472 
473 template<class SolidThermophysicalTransportModel>
476 (
477  const polyTopoChangeMap& map
478 )
479 {
480  // Delete the cached zonesPatchFaces, will be re-created in predict
481  zonesPatchFaces_.clear();
482 }
483 
484 
485 template<class SolidThermophysicalTransportModel>
488 (
489  const polyMeshMap& map
490 )
491 {
492  // Delete the cached zonesPatchFaces, will be re-created in predict
493  zonesPatchFaces_.clear();
494 }
495 
496 
497 template<class SolidThermophysicalTransportModel>
500 (
501  const polyDistributionMap& map
502 )
503 {
504  // Delete the cached zonesPatchFaces, will be re-created in predict
505  zonesPatchFaces_.clear();
506 }
507 
508 
509 // ************************************************************************* //
bool found
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::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
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 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:162
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const volVectorField & C() const
Return cell centres.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1762
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 cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
Solid thermophysical transport model for anisotropic thermal conductivity.
Definition: anisotropic.H:121
virtual bool movePoints()
Update for mesh motion.
Definition: anisotropic.C:467
SolidThermophysicalTransportModel::alphaField alphaField
Definition: anisotropic.H:154
virtual tmp< scalarField > qCorr(const label patchi) const
Return the patch heat flux correction [W/m^2].
Definition: anisotropic.C:400
virtual void predict()
Correct the anisotropic viscosity.
Definition: anisotropic.C:452
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
Definition: anisotropic.C:488
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: anisotropic.C:500
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:383
anisotropic(const alphaField &alpha, const solidThermo &thermo)
Construct from solid thermophysical properties.
Definition: anisotropic.C:101
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
Definition: anisotropic.C:476
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: anisotropic.C:213
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: anisotropic.C:428
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:381
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
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
const doubleScalar e
Definition: doubleScalar.H:106
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:257
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 > &)
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:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
faceListList boundary(nPatches)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31