laserDTRM.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) 2017-2022 OpenCFD Ltd.
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 
28 #include "laserDTRM.H"
29 #include "fvmLaplacian.H"
30 #include "fvmSup.H"
32 #include "scatterModel.H"
33 #include "constants.H"
34 #include "unitConversion.H"
35 #include "interpolationCell.H"
36 #include "interpolationCellPoint.H"
37 #include "Random.H"
38 #include "OBJstream.H"
40 
41 using namespace Foam::constant;
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  namespace radiation
48  {
49  defineTypeNameAndDebug(laserDTRM, 0);
51  }
52 
54  (
55  Cloud<DTRMParticle>,
56  "DTRMCloud",
57  0
58  );
59 
60 } // End namespace Foam
61 
62 
64 Foam::radiation::laserDTRM::powerDistNames_
65 {
66  { powerDistributionMode::pdGaussian, "Gaussian" },
67  { powerDistributionMode::pdManual, "manual" },
68  { powerDistributionMode::pdUniform, "uniform" },
69  { powerDistributionMode::pdGaussianPeak, "GaussianPeak" },
70 };
71 
72 
73 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
74 
75 Foam::scalar Foam::radiation::laserDTRM::calculateIp(scalar r, scalar theta)
76 {
77  const scalar t = mesh_.time().value();
78  const scalar power = laserPower_->value(t);
79  switch (mode_)
80  {
81  case pdGaussianPeak:
82  {
83  return I0_*exp(-2.0*sqr(r)/sqr(sigma_));
84  break;
85  }
86  case pdGaussian:
87  {
88  scalar I0 = power/(mathematical::twoPi*sqr(sigma_));
89 
90  return I0*exp(-sqr(r)/2.0/sqr(sigma_));
91  break;
92  }
93  case pdManual:
94  {
95  return power*powerDistribution_()(theta, r);
96  break;
97  }
98  case pdUniform:
99  {
100  return power/(mathematical::pi*sqr(focalLaserRadius_));
101  break;
102  }
103  default:
104  {
106  << "Unhandled type " << powerDistNames_[mode_]
107  << abort(FatalError);
108  break;
109  }
110  }
111 
112  return 0;
113 }
114 
115 
116 Foam::tmp<Foam::volVectorField> Foam::radiation::laserDTRM::nHatfv
117 (
118  const volScalarField& alpha1,
119  const volScalarField& alpha2
120 ) const
121 {
122  const dimensionedScalar deltaN
123  (
124  "deltaN",
125  1e-7/cbrt(average(mesh_.V()))
126  );
127 
128  const volVectorField gradAlphaf
129  (
132  );
133 
134  // Face unit interface normal
135  return gradAlphaf/(mag(gradAlphaf)+ deltaN);
136 }
137 
138 
139 void Foam::radiation::laserDTRM::initialiseReflection()
140 {
141  if (found("reflectionModel"))
142  {
143  dictTable modelDicts(lookup("reflectionModel"));
144 
145  forAllConstIters(modelDicts, iter)
146  {
147  const phasePairKey& key = iter.key();
148 
149  reflections_.insert
150  (
151  key,
153  (
154  iter.val(),
155  mesh_
156  )
157  );
158  }
159 
160  reflectionSwitch_ = returnReduceOr(reflections_.size());
161  }
162 }
163 
164 
165 void Foam::radiation::laserDTRM::initialise()
166 {
167  // Initialise the DTRM particles
168  DTRMCloud_.clear();
169 
170  const scalar t = mesh_.time().value();
171  const vector lPosition = focalLaserPosition_->value(t);
172  const vector lDir = normalised(laserDirection_->value(t));
173 
174  DebugInfo
175  << "Laser position : " << lPosition << nl
176  << "Laser direction : " << lDir << endl;
177 
178  // Find a vector on the area plane. Normal to laser direction
179  vector rArea = Zero;
180  scalar magr = 0.0;
181 
182  {
183  Random rnd(1234);
184 
185  while (magr < VSMALL)
186  {
187  vector v = rnd.sample01<vector>();
188  rArea = v - (v & lDir)*lDir;
189  magr = mag(rArea);
190  }
191  }
192  rArea.normalise();
193 
194  scalar dr = focalLaserRadius_/ndr_;
195  scalar dTheta = mathematical::twoPi/ndTheta_;
196 
197  nParticles_ = ndr_*ndTheta_;
198 
199  switch (mode_)
200  {
201  case pdGaussianPeak:
202  {
203  I0_ = get<scalar>("I0");
204  sigma_ = get<scalar>("sigma");
205  break;
206  }
207  case pdGaussian:
208  {
209  sigma_ = get<scalar>("sigma");
210  break;
211  }
212  case pdManual:
213  {
214  powerDistribution_.reset
215  (
216  new interpolation2DTable<scalar>(*this)
217  );
218 
219  break;
220  }
221  case pdUniform:
222  {
223  break;
224  }
225  }
226 
227  // Count the number of missed positions
228  label nMissed = 0;
229 
230  // Target position
231  point p1 = vector::zero;
232 
233  // Seed DTRM particles
234  // TODO: currently only applicable to 3-D cases
235  point p0 = lPosition;
236  scalar power(0.0);
237  scalar area(0.0);
238  p1 = p0;
239  if (mesh_.nGeometricD() == 3)
240  {
241  for (label ri = 0; ri < ndr_; ri++)
242  {
243  scalar r1 = SMALL + dr*ri;
244 
245  scalar r2 = r1 + dr;
246 
247  scalar rP = ((r1 + r2)/2);
248 
249  // local radius on disk
250  vector localR = ((r1 + r2)/2)*rArea;
251 
252  // local final radius on disk
253  vector finalR = rP*rArea;
254 
255  scalar theta0 = 0.0;//dTheta/2.0;
256  for (label thetai = 0; thetai < ndTheta_; thetai++)
257  {
258  scalar theta1 = theta0 + SMALL + dTheta*thetai;
259 
260  scalar theta2 = theta1 + dTheta;
261 
262  scalar thetaP = (theta1 + theta2)/2.0;
263 
264  quaternion Q(lDir, thetaP);
265 
266  // Initial position on disk
267  vector initialPos = (Q.R() & localR);
268 
269  // Final position on disk
270  vector finalPos = (Q.R() & finalR);
271 
272  // Initial position
273  p0 = lPosition + initialPos;
274 
275  // calculate target point using new deviation rl
276  p1 = lPosition + finalPos + (0.5*maxTrackLength_*lDir);
277 
278  scalar Ip = calculateIp(rP, thetaP);
279 
280  scalar dAi = (sqr(r2) - sqr(r1))*(theta2 - theta1)/2.0;
281 
282  power += Ip*dAi;
283  area += dAi;
284 
285  label cellI = mesh_.findCell(p0);
286 
287  if (cellI != -1)
288  {
289  // Create a new particle
290  DTRMParticle* pPtr =
291  new DTRMParticle(mesh_, p0, p1, Ip, cellI, dAi, -1);
292 
293  // Add to cloud
294  DTRMCloud_.addParticle(pPtr);
295  }
296 
297  if (nMissed < 10 && returnReduceAnd(cellI < 0))
298  {
299  ++nMissed;
301  << "Cannot find owner cell for focalPoint at "
302  << p0 << endl;
303  }
304  }
305  }
306  }
307  else
308  {
310  << "Current functionality limited to 3-D cases"
311  << exit(FatalError);
312  }
313 
314  if (nMissed)
315  {
316  Info<< "Seeding missed " << nMissed << " locations" << endl;
317  }
318 
319  DebugInfo
320  << "Total Power in the laser : " << power << nl
321  << "Total Area in the laser : " << area << nl
322  << endl;
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
327 
328 Foam::radiation::laserDTRM::laserDTRM(const volScalarField& T)
329 :
330  radiationModel(typeName, T),
331  mode_(powerDistNames_.get("mode", *this)),
332  DTRMCloud_(mesh_, Foam::zero{}, "DTRMCloud"), // Empty cloud
333  nParticles_(0),
334  ndTheta_(get<label>("nTheta")),
335  ndr_(get<label>("nr")),
336  maxTrackLength_(mesh_.bounds().mag()),
337 
338  focalLaserPosition_
339  (
340  Function1<point>::New("focalLaserPosition", *this, &mesh_)
341  ),
342 
343  laserDirection_
344  (
345  Function1<vector>::New("laserDirection", *this, &mesh_)
346  ),
347 
348  focalLaserRadius_(get<scalar>("focalLaserRadius")),
349  qualityBeamLaser_
350  (
351  getOrDefault<scalar>("qualityBeamLaser", 0)
352  ),
353 
354  sigma_(0),
355  I0_(0),
356  laserPower_(Function1<scalar>::New("laserPower", *this, &mesh_)),
357  powerDistribution_(),
358 
359  reflectionSwitch_(false),
360 
361  alphaCut_(getOrDefault<scalar>("alphaCut", 0.5)),
362 
363  a_
364  (
365  IOobject
366  (
367  "a",
368  mesh_.time().timeName(),
369  mesh_,
372  ),
373  mesh_,
375  ),
376  e_
377  (
378  IOobject
379  (
380  "e",
381  mesh_.time().timeName(),
382  mesh_,
385  ),
386  mesh_,
388  ),
389  E_
390  (
391  IOobject
392  (
393  "E",
394  mesh_.time().timeName(),
395  mesh_,
398  ),
399  mesh_,
401  ),
402  Q_
403  (
404  IOobject
405  (
406  "Q",
407  mesh_.time().timeName(),
408  mesh_,
411  ),
412  mesh_,
414  )
415 {
416  initialiseReflection();
417 
418  initialise();
419 }
420 
421 
422 Foam::radiation::laserDTRM::laserDTRM
423 (
424  const dictionary& dict,
425  const volScalarField& T
426 )
427 :
428  radiationModel(typeName, dict, T),
429  mode_(powerDistNames_.get("mode", *this)),
430  DTRMCloud_(mesh_, Foam::zero{}, "DTRMCloud"), // Empty cloud
431  nParticles_(0),
432  ndTheta_(get<label>("nTheta")),
433  ndr_(get<label>("nr")),
434  maxTrackLength_(mesh_.bounds().mag()),
435 
436  focalLaserPosition_
437  (
438  Function1<point>::New("focalLaserPosition", *this, &mesh_)
439  ),
440  laserDirection_
441  (
442  Function1<vector>::New("laserDirection", *this, &mesh_)
443  ),
444 
445  focalLaserRadius_(get<scalar>("focalLaserRadius")),
446  qualityBeamLaser_
447  (
448  getOrDefault<scalar>("qualityBeamLaser", 0)
449  ),
450 
451  sigma_(0),
452  I0_(0),
453  laserPower_(Function1<scalar>::New("laserPower", *this, &mesh_)),
454  powerDistribution_(),
455 
456  reflectionSwitch_(false),
457 
458  alphaCut_(getOrDefault<scalar>("alphaCut", 0.5)),
459 
460  a_
461  (
462  IOobject
463  (
464  "a",
465  mesh_.time().timeName(),
466  mesh_,
469  ),
470  mesh_,
472  ),
473  e_
474  (
475  IOobject
476  (
477  "e",
478  mesh_.time().timeName(),
479  mesh_,
482  ),
483  mesh_,
485  ),
486  E_
487  (
488  IOobject
489  (
490  "E",
491  mesh_.time().timeName(),
492  mesh_,
495  ),
496  mesh_,
498  ),
499  Q_
500  (
501  IOobject
502  (
503  "Q",
504  mesh_.time().timeName(),
505  mesh_,
508  ),
509  mesh_,
511  )
512 {
513  initialiseReflection();
514  initialise();
515 }
516 
517 
518 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
519 
521 {
522  if (radiationModel::read())
523  {
524  return true;
525  }
526 
527  return false;
528 }
529 
530 Foam::label Foam::radiation::laserDTRM::nBands() const
531 {
532  return 1;
533 }
534 
535 
537 {
538  tmp<volScalarField> treflectingCells
539  (
540  new volScalarField
541  (
542  IOobject
543  (
544  "reflectingCellsVol",
545  mesh_.time().timeName(),
546  mesh_,
549  ),
550  mesh_,
551  dimensionedScalar("zero", dimless, -1)
552  )
553  );
554  volScalarField& reflectingCellsVol = treflectingCells.ref();
555 
556 
557  tmp<volVectorField> tnHat
558  (
559  new volVectorField
560  (
561  IOobject
562  (
563  "nHat",
564  mesh_.time().timeName(),
565  mesh_,
568  ),
569  mesh_,
571  )
572  );
573  volVectorField& nHat = tnHat.ref();
574 
575 
576  // Reset the field
577  Q_ == dimensionedScalar(Q_.dimensions(), Zero);
578 
579  a_ = absorptionEmission_->a();
580  e_ = absorptionEmission_->e();
581  E_ = absorptionEmission_->E();
582 
583  const interpolationCell<scalar> aInterp(a_);
584  const interpolationCell<scalar> eInterp(e_);
585  const interpolationCell<scalar> EInterp(E_);
586  const interpolationCell<scalar> TInterp(T_);
587 
588  labelField reflectingCells(mesh_.nCells(), -1);
589 
590  UPtrList<reflectionModel> reflectionUPtr;
591 
592  if (reflectionSwitch_)
593  {
594  reflectionUPtr.resize(reflections_.size());
595 
596  label reflectionModelId(0);
597  forAllIters(reflections_, iter1)
598  {
599  reflectionModel& model = iter1()();
600 
601  reflectionUPtr.set(reflectionModelId, &model);
602 
603  const volScalarField& alphaFrom =
604  mesh_.lookupObject<volScalarField>
605  (
606  IOobject::groupName("alpha", iter1.key().first())
607  );
608 
609  const volScalarField& alphaTo =
610  mesh_.lookupObject<volScalarField>
611  (
612  IOobject::groupName("alpha", iter1.key().second())
613  );
614 
615  const volVectorField nHatPhase(nHatfv(alphaFrom, alphaTo));
616 
617  const volScalarField gradAlphaf
618  (
619  fvc::grad(alphaFrom)
620  & fvc::grad(alphaTo)
621  );
622 
623  const volScalarField nearInterface(pos(alphaTo - alphaCut_));
624 
625  const volScalarField mask(nearInterface*gradAlphaf);
626 
627  forAll(alphaFrom, cellI)
628  {
629  if
630  (
631  nearInterface[cellI]
632  && mag(nHatPhase[cellI]) > 0.99
633  && mask[cellI] < 0
634  )
635  {
636  reflectingCells[cellI] = reflectionModelId;
637  reflectingCellsVol[cellI] = reflectionModelId;
638  if (mag(nHat[cellI]) == 0.0)
639  {
640  nHat[cellI] += nHatPhase[cellI];
641  }
642  }
643  }
644  reflectionModelId++;
645  }
646  }
647 
648  interpolationCellPoint<vector> nHatInterp(nHat);
649 
650  DTRMParticle::trackingData td
651  (
652  DTRMCloud_,
653  aInterp,
654  eInterp,
655  EInterp,
656  TInterp,
657  nHatInterp,
658  reflectingCells,
659  reflectionUPtr,
660  Q_
661  );
662 
663  Info<< "Move particles..."
664  << returnReduce(DTRMCloud_.size(), sumOp<label>()) << endl;
665 
666  DTRMCloud_.move(DTRMCloud_, td, mesh_.time().deltaTValue());
667 
668  // Normalize by cell volume
669  Q_.primitiveFieldRef() /= mesh_.V();
670 
671  if (debug)
672  {
673  Info<< "Final number of particles..."
674  << returnReduce(DTRMCloud_.size(), sumOp<label>()) << endl;
675 
676  pointField lines(2*DTRMCloud_.size());
677  {
678  label i = 0;
679  for (const DTRMParticle& p : DTRMCloud_)
680  {
681  lines[i] = p.position();
682  lines[i+1] = p.p0();
683  i += 2;
684  }
685  }
686 
688 
689  if (Pstream::master())
690  {
691  OBJstream os(type() + ":particlePath.obj");
692 
693  for (label pointi = 0; pointi < lines.size(); pointi += 2)
694  {
695  os.writeLine(lines[pointi], lines[pointi+1]);
696  }
697  }
698 
699  scalar totalQ = gSum(Q_.primitiveFieldRef()*mesh_.V());
700  Info << "Total energy absorbed [W]: " << totalQ << endl;
701 
702  if (mesh_.time().writeTime())
703  {
704  reflectingCellsVol.write();
705  nHat.write();
706  }
707  }
708 
709  // Clear and initialise the cloud
710  // NOTE: Possible to reset original particles, but this requires
711  // data transfer for the cloud in differet processors.
712  initialise();
713 }
714 
715 
717 {
719  (
720  IOobject
721  (
722  "zero",
723  mesh_.time().timeName(),
724  mesh_,
728  ),
729  mesh_,
731  );
732 }
733 
734 
737 {
738  return Q_.internalField();
739 }
740 
741 
742 // ************************************************************************* //
Different types of constants.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
bool read()
Read radiation properties dictionary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
static autoPtr< reflectionModel > New(const dictionary &dict, const fvMesh &mesh)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
DTRMParticle(const polyMesh &mesh, const vector &position, const vector &targetPosition, const scalar I, const label cellI, const scalar dA, const label transmissiveId)
Construct from components, with searching for tetFace and.
Calculate the matrix for the laplacian of the field.
const volScalarField & alpha2
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const dimensionSet dimless
Dimensionless.
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual label nBands() const
Number of bands for this radiation model.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar exp(const dimensionedScalar &ds)
virtual bool read()=0
Read radiationProperties dictionary.
constexpr scalar twoPi(2 *M_PI)
const wordList area
Standard area field types (scalar, vector, tensor, etc)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar cbrt(const dimensionedScalar &ds)
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
constexpr scalar pi(M_PI)
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
const dimensionSet dimPower
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
void calculate()
Solve radiation equation(s)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar pow3(const dimensionedScalar &ds)
vector point
Point is a vector.
Definition: point.H:37
static void gatherInplaceOp(List< Type > &fld, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking, const label comm=UPstream::worldComm)
Inplace collect data in processor order on master (in serial: a no-op).
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
dimensionedScalar pow4(const dimensionedScalar &ds)
Nothing to be read.
Automatically write from objectRegistry::writeObject()
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
#define addToRadiationRunTimeSelectionTables(model)
bool found
defineTemplateTypeNameAndDebugWithName(psiReactionsSensitivityAnalysisFunctionObject, "psiReactionsSensitivityAnalysis", 0)
Do not request registration (bool: false)
const volScalarField & p0
Definition: EEqn.H:36
Calculate the finiteVolume matrix for implicit and explicit sources.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1