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-2020 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"
35 #include "unitConversion.H"
36 #include "interpolationCell.H"
37 #include "interpolationCellPoint.H"
38 #include "Random.H"
39 
40 using namespace Foam::constant;
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  namespace radiation
47  {
48  defineTypeNameAndDebug(laserDTRM, 0);
50  }
51 
53  (
54  Cloud<DTRMParticle>,
55  "DTRMCloud",
56  0
57  );
58 
59 } // End namespace Foam
60 
61 
63 Foam::radiation::laserDTRM::powerDistNames_
64 {
65  { powerDistributionMode::pdGaussian, "Gaussian" },
66  { powerDistributionMode::pdManual, "manual" },
67  { powerDistributionMode::pdUniform, "uniform" },
68  { powerDistributionMode::pdGaussianPeak, "GaussianPeak" },
69 };
70 
71 
72 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
73 
74 Foam::scalar Foam::radiation::laserDTRM::calculateIp(scalar r, scalar theta)
75 {
76  const scalar t = mesh_.time().value();
77  const scalar power = laserPower_->value(t);
78  switch (mode_)
79  {
80  case pdGaussianPeak:
81  {
82  return I0_*exp(-2.0*sqr(r)/sqr(sigma_));
83  break;
84  }
85  case pdGaussian:
86  {
87  scalar I0 = power/(mathematical::twoPi*sqr(sigma_));
88 
89  return I0*exp(-sqr(r)/2.0/sqr(sigma_));
90  break;
91  }
92  case pdManual:
93  {
94  return power*powerDistribution_()(theta, r);
95  break;
96  }
97  case pdUniform:
98  {
99  return power/(mathematical::pi*sqr(focalLaserRadius_));
100  break;
101  }
102  default:
103  {
105  << "Unhandled type " << powerDistNames_[mode_]
106  << abort(FatalError);
107  break;
108  }
109  }
110 
111  return 0;
112 }
113 
114 
115 Foam::tmp<Foam::volVectorField> Foam::radiation::laserDTRM::nHatfv
116 (
117  const volScalarField& alpha1,
118  const volScalarField& alpha2
119 ) const
120 {
121  const dimensionedScalar deltaN
122  (
123  "deltaN",
124  1e-7/cbrt(average(mesh_.V()))
125  );
126 
127  const volVectorField gradAlphaf
128  (
131  );
132 
133  // Face unit interface normal
134  return gradAlphaf/(mag(gradAlphaf)+ deltaN);
135 }
136 
137 
138 void Foam::radiation::laserDTRM::initialiseReflection()
139 {
140  if (found("reflectionModel"))
141  {
142  dictTable modelDicts(lookup("reflectionModel"));
143 
144  forAllConstIters(modelDicts, iter)
145  {
146  const phasePairKey& key = iter.key();
147 
148  reflections_.insert
149  (
150  key,
152  (
153  iter.val(),
154  mesh_
155  )
156  );
157  }
158 
159  if (reflections_.size())
160  {
161  reflectionSwitch_ = true;
162  }
163 
164  reflectionSwitch_ = returnReduce(reflectionSwitch_, orOp<bool>());
165  }
166 }
167 
168 
169 void Foam::radiation::laserDTRM::initialise()
170 {
171  // Initialise the DTRM particles
172  DTRMCloud_.clear();
173 
174  const scalar t = mesh_.time().value();
175  const vector lPosition = focalLaserPosition_->value(t);
176  const vector lDir = normalised(laserDirection_->value(t));
177 
178  DebugInfo
179  << "Laser position : " << lPosition << nl
180  << "Laser direction : " << lDir << endl;
181 
182  // Find a vector on the area plane. Normal to laser direction
183  vector rArea = Zero;
184  scalar magr = 0.0;
185 
186  {
187  Random rnd(1234);
188 
189  while (magr < VSMALL)
190  {
191  vector v = rnd.sample01<vector>();
192  rArea = v - (v & lDir)*lDir;
193  magr = mag(rArea);
194  }
195  }
196  rArea.normalise();
197 
198  scalar dr = focalLaserRadius_/ndr_;
199  scalar dTheta = mathematical::twoPi/ndTheta_;
200 
201  nParticles_ = ndr_*ndTheta_;
202 
203  switch (mode_)
204  {
205  case pdGaussianPeak:
206  {
207  I0_ = get<scalar>("I0");
208  sigma_ = get<scalar>("sigma");
209  break;
210  }
211  case pdGaussian:
212  {
213  sigma_ = get<scalar>("sigma");
214  break;
215  }
216  case pdManual:
217  {
218  powerDistribution_.reset
219  (
220  new interpolation2DTable<scalar>(*this)
221  );
222 
223  break;
224  }
225  case pdUniform:
226  {
227  break;
228  }
229  }
230 
231  // Count the number of missed positions
232  label nMissed = 0;
233 
234  // Target position
235  point p1 = vector::zero;
236 
237  // Seed DTRM particles
238  // TODO: currently only applicable to 3-D cases
239  point p0 = lPosition;
240  scalar power(0.0);
241  scalar area(0.0);
242  p1 = p0;
243  if (mesh_.nGeometricD() == 3)
244  {
245  for (label ri = 0; ri < ndr_; ri++)
246  {
247  scalar r1 = SMALL + dr*ri;
248 
249  scalar r2 = r1 + dr;
250 
251  scalar rP = ((r1 + r2)/2);
252 
253  // local radius on disk
254  vector localR = ((r1 + r2)/2)*rArea;
255 
256  // local final radius on disk
257  vector finalR = rP*rArea;
258 
259  scalar theta0 = 0.0;//dTheta/2.0;
260  for (label thetai = 0; thetai < ndTheta_; thetai++)
261  {
262  scalar theta1 = theta0 + SMALL + dTheta*thetai;
263 
264  scalar theta2 = theta1 + dTheta;
265 
266  scalar thetaP = (theta1 + theta2)/2.0;
267 
268  quaternion Q(lDir, thetaP);
269 
270  // Initial position on disk
271  vector initialPos = (Q.R() & localR);
272 
273  // Final position on disk
274  vector finalPos = (Q.R() & finalR);
275 
276  // Initial position
277  p0 = lPosition + initialPos;
278 
279  // calculate target point using new deviation rl
280  p1 = lPosition + finalPos + (0.5*maxTrackLength_*lDir);
281 
282  scalar Ip = calculateIp(rP, thetaP);
283 
284  scalar dAi = (sqr(r2) - sqr(r1))*(theta2 - theta1)/2.0;
285 
286  power += Ip*dAi;
287  area += dAi;
288 
289  label cellI = mesh_.findCell(p0);
290 
291  if (cellI != -1)
292  {
293  // Create a new particle
294  DTRMParticle* pPtr =
295  new DTRMParticle(mesh_, p0, p1, Ip, cellI, dAi, -1);
296 
297  // Add to cloud
298  DTRMCloud_.addParticle(pPtr);
299  }
300 
301  if (returnReduce(cellI, maxOp<label>()) == -1)
302  {
303  if (++nMissed <= 10)
304  {
306  << "Cannot find owner cell for focalPoint at "
307  << p0 << endl;
308  }
309  }
310  }
311  }
312  }
313  else
314  {
316  << "Current functionality limited to 3-D cases"
317  << exit(FatalError);
318  }
319 
320  if (nMissed)
321  {
322  Info<< "Seeding missed " << nMissed << " locations" << endl;
323  }
324 
325  DebugInfo
326  << "Total Power in the laser : " << power << nl
327  << "Total Area in the laser : " << area << nl
328  << endl;
329 }
330 
331 
332 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
333 
334 Foam::radiation::laserDTRM::laserDTRM(const volScalarField& T)
335 :
336  radiationModel(typeName, T),
337  mode_(powerDistNames_.get("mode", *this)),
338  DTRMCloud_(mesh_, "DTRMCloud", IDLList<DTRMParticle>()),
339  nParticles_(0),
340  ndTheta_(get<label>("nTheta")),
341  ndr_(get<label>("nr")),
342  maxTrackLength_(mesh_.bounds().mag()),
343 
344  focalLaserPosition_
345  (
346  Function1<point>::New("focalLaserPosition", *this, &mesh_)
347  ),
348 
349  laserDirection_
350  (
351  Function1<vector>::New("laserDirection", *this, &mesh_)
352  ),
353 
354  focalLaserRadius_(get<scalar>("focalLaserRadius")),
355  qualityBeamLaser_
356  (
357  getOrDefault<scalar>("qualityBeamLaser", 0)
358  ),
359 
360  sigma_(0),
361  I0_(0),
362  laserPower_(Function1<scalar>::New("laserPower", *this, &mesh_)),
363  powerDistribution_(),
364 
365  reflectionSwitch_(false),
366 
367  alphaCut_(getOrDefault<scalar>("alphaCut", 0.5)),
368 
369  a_
370  (
371  IOobject
372  (
373  "a",
374  mesh_.time().timeName(),
375  mesh_,
376  IOobject::NO_READ,
377  IOobject::NO_WRITE
378  ),
379  mesh_,
381  ),
382  e_
383  (
384  IOobject
385  (
386  "e",
387  mesh_.time().timeName(),
388  mesh_,
389  IOobject::NO_READ,
390  IOobject::NO_WRITE
391  ),
392  mesh_,
394  ),
395  E_
396  (
397  IOobject
398  (
399  "E",
400  mesh_.time().timeName(),
401  mesh_,
402  IOobject::NO_READ,
403  IOobject::NO_WRITE
404  ),
405  mesh_,
407  ),
408  Q_
409  (
410  IOobject
411  (
412  "Q",
413  mesh_.time().timeName(),
414  mesh_,
415  IOobject::NO_READ,
416  IOobject::AUTO_WRITE
417  ),
418  mesh_,
420  )
421 {
422  initialiseReflection();
423 
424  initialise();
425 }
426 
427 
428 Foam::radiation::laserDTRM::laserDTRM
429 (
430  const dictionary& dict,
431  const volScalarField& T
432 )
433 :
434  radiationModel(typeName, dict, T),
435  mode_(powerDistNames_.get("mode", *this)),
436  DTRMCloud_(mesh_, "DTRMCloud", IDLList<DTRMParticle>()),
437  nParticles_(0),
438  ndTheta_(get<label>("nTheta")),
439  ndr_(get<label>("nr")),
440  maxTrackLength_(mesh_.bounds().mag()),
441 
442  focalLaserPosition_
443  (
444  Function1<point>::New("focalLaserPosition", *this, &mesh_)
445  ),
446  laserDirection_
447  (
448  Function1<vector>::New("laserDirection", *this, &mesh_)
449  ),
450 
451  focalLaserRadius_(get<scalar>("focalLaserRadius")),
452  qualityBeamLaser_
453  (
454  getOrDefault<scalar>("qualityBeamLaser", 0)
455  ),
456 
457  sigma_(0),
458  I0_(0),
459  laserPower_(Function1<scalar>::New("laserPower", *this, &mesh_)),
460  powerDistribution_(),
461 
462  reflectionSwitch_(false),
463 
464  alphaCut_(getOrDefault<scalar>("alphaCut", 0.5)),
465 
466  a_
467  (
468  IOobject
469  (
470  "a",
471  mesh_.time().timeName(),
472  mesh_,
473  IOobject::NO_READ,
474  IOobject::NO_WRITE
475  ),
476  mesh_,
478  ),
479  e_
480  (
481  IOobject
482  (
483  "e",
484  mesh_.time().timeName(),
485  mesh_,
486  IOobject::NO_READ,
487  IOobject::NO_WRITE
488  ),
489  mesh_,
491  ),
492  E_
493  (
494  IOobject
495  (
496  "E",
497  mesh_.time().timeName(),
498  mesh_,
499  IOobject::NO_READ,
500  IOobject::NO_WRITE
501  ),
502  mesh_,
504  ),
505  Q_
506  (
507  IOobject
508  (
509  "Q",
510  mesh_.time().timeName(),
511  mesh_,
512  IOobject::NO_READ,
513  IOobject::AUTO_WRITE
514  ),
515  mesh_,
517  )
518 {
519  initialiseReflection();
520  initialise();
521 }
522 
523 
524 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
525 
527 {
528  if (radiationModel::read())
529  {
530  return true;
531  }
532 
533  return false;
534 }
535 
536 Foam::label Foam::radiation::laserDTRM::nBands() const
537 {
538  return 1;
539 }
540 
541 
543 {
544  tmp<volScalarField> treflectingCells
545  (
546  new volScalarField
547  (
548  IOobject
549  (
550  "reflectingCellsVol",
551  mesh_.time().timeName(),
552  mesh_,
555  ),
556  mesh_,
557  dimensionedScalar("zero", dimless, -1)
558  )
559  );
560  volScalarField& reflectingCellsVol = treflectingCells.ref();
561 
562 
563  tmp<volVectorField> tnHat
564  (
565  new volVectorField
566  (
567  IOobject
568  (
569  "nHat",
570  mesh_.time().timeName(),
571  mesh_,
574  ),
575  mesh_,
577  )
578  );
579  volVectorField& nHat = tnHat.ref();
580 
581 
582  // Reset the field
583  Q_ == dimensionedScalar(Q_.dimensions(), Zero);
584 
585  a_ = absorptionEmission_->a();
586  e_ = absorptionEmission_->e();
587  E_ = absorptionEmission_->E();
588 
589  const interpolationCell<scalar> aInterp(a_);
590  const interpolationCell<scalar> eInterp(e_);
591  const interpolationCell<scalar> EInterp(E_);
592  const interpolationCell<scalar> TInterp(T_);
593 
594  labelField reflectingCells(mesh_.nCells(), -1);
595 
596  UPtrList<reflectionModel> reflectionUPtr;
597 
598  if (reflectionSwitch_)
599  {
600  reflectionUPtr.resize(reflections_.size());
601 
602  label reflectionModelId(0);
603  forAllIters(reflections_, iter1)
604  {
605  reflectionModel& model = iter1()();
606 
607  reflectionUPtr.set(reflectionModelId, &model);
608 
609  const volScalarField& alphaFrom =
610  mesh_.lookupObject<volScalarField>
611  (
612  IOobject::groupName("alpha", iter1.key().first())
613  );
614 
615  const volScalarField& alphaTo =
616  mesh_.lookupObject<volScalarField>
617  (
618  IOobject::groupName("alpha", iter1.key().second())
619  );
620 
621  const volVectorField nHatPhase(nHatfv(alphaFrom, alphaTo));
622 
623  const volScalarField gradAlphaf
624  (
625  fvc::grad(alphaFrom)
626  & fvc::grad(alphaTo)
627  );
628 
629  const volScalarField nearInterface(pos(alphaTo - alphaCut_));
630 
631  const volScalarField mask(nearInterface*gradAlphaf);
632 
633  forAll(alphaFrom, cellI)
634  {
635  if
636  (
637  nearInterface[cellI]
638  && mag(nHatPhase[cellI]) > 0.99
639  && mask[cellI] < 0
640  )
641  {
642  reflectingCells[cellI] = reflectionModelId;
643  reflectingCellsVol[cellI] = reflectionModelId;
644  if (mag(nHat[cellI]) == 0.0)
645  {
646  nHat[cellI] += nHatPhase[cellI];
647  }
648  }
649  }
650  reflectionModelId++;
651  }
652  }
653 
654  interpolationCellPoint<vector> nHatInterp(nHat);
655 
656  DTRMParticle::trackingData td
657  (
658  DTRMCloud_,
659  aInterp,
660  eInterp,
661  EInterp,
662  TInterp,
663  nHatInterp,
664  reflectingCells,
665  reflectionUPtr,
666  Q_
667  );
668 
669  Info<< "Move particles..."
670  << returnReduce(DTRMCloud_.size(), sumOp<label>()) << endl;
671 
672  DTRMCloud_.move(DTRMCloud_, td, mesh_.time().deltaTValue());
673 
674  // Normalize by cell volume
675  Q_.primitiveFieldRef() /= mesh_.V();
676 
677  if (debug)
678  {
679  Info<< "Final number of particles..."
680  << returnReduce(DTRMCloud_.size(), sumOp<label>()) << endl;
681 
682  OFstream osRef(type() + ":particlePath.obj");
683  label vertI = 0;
684 
685  List<pointField> positions(Pstream::nProcs());
686  List<pointField> p0(Pstream::nProcs());
687 
688  DynamicList<point> positionsMyProc;
689  DynamicList<point> p0MyProc;
690 
691  for (const DTRMParticle& p : DTRMCloud_)
692  {
693  positionsMyProc.append(p.position());
694  p0MyProc.append(p.p0());
695  }
696 
697  positions[Pstream::myProcNo()].transfer(positionsMyProc);
698  p0[Pstream::myProcNo()].transfer(p0MyProc);
699 
700  Pstream::gatherList(positions);
701  Pstream::scatterList(positions);
704 
705  for (const int proci : Pstream::allProcs())
706  {
707  const pointField& pos = positions[proci];
708  const pointField& pfinal = p0[proci];
709  forAll(pos, i)
710  {
711  meshTools::writeOBJ(osRef, pos[i]);
712  vertI++;
713  meshTools::writeOBJ(osRef, pfinal[i]);
714  vertI++;
715  osRef << "l " << vertI-1 << ' ' << vertI << nl;
716  }
717  }
718 
719  osRef.flush();
720 
721  scalar totalQ = gSum(Q_.primitiveFieldRef()*mesh_.V());
722  Info << "Total energy absorbed [W]: " << totalQ << endl;
723 
724  if (mesh_.time().outputTime())
725  {
726  reflectingCellsVol.write();
727  nHat.write();
728  }
729  }
730 
731  // Clear and initialise the cloud
732  // NOTE: Possible to reset original particles, but this requires
733  // data transfer for the cloud in differet processors.
734  initialise();
735 }
736 
737 
739 {
741  (
742  IOobject
743  (
744  "zero",
745  mesh_.time().timeName(),
746  mesh_,
749  false
750  ),
751  mesh_,
753  );
754 }
755 
756 
759 {
760  return Q_.internalField();
761 }
762 
763 
764 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:53
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
scatterModel.H
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Definition: meshTools.C:196
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
Foam::constant
Different types of constants.
Definition: atomicConstants.C:31
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Definition: gatherScatterList.C:209
Foam::IDLList
ILList< DLListBase, T > IDLList
Definition: IDLList.H:39
interpolationCell.H
Foam::glTF::key
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
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:587
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:254
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
interpolationCellPoint.H
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:46
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:93
Foam::radiation::radiationModel::read
virtual bool read()=0
Definition: radiationModel.C:203
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:82
Foam::Info
messageStream Info
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::radiation::laserDTRM::calculate
void calculate()
Foam::radiation::laserDTRM::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
forAllIters
#define forAllIters(container, iter)
Definition: stdFoam.H:264
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
Foam::dimPower
const dimensionSet dimPower
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::radiation::laserDTRM::Rp
virtual tmp< volScalarField > Rp() const
Foam::FatalError
error FatalError
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::Ostream::write
virtual bool write(const token &tok)=0
Foam
Definition: atmBoundaryLayer.C:26
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:47
Random.H
constants.H
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
laserDTRM.H
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:480
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::radiation::laserDTRM::nBands
virtual label nBands() const
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Definition: gatherScatterList.C:46
Foam::UPstream::allProcs
static rangeType allProcs(const label communicator=worldComm)
Definition: UPstream.H:504
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
DebugInfo
#define DebugInfo
Definition: messageStream.H:436
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Definition: UPstream.H:459
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:424
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
absorptionEmissionModel.H
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Definition: GeometricField.C:742
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Foam::radiation::laserDTRM::read
bool read()
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::fieldTypes::area
const wordList area
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Form zero
Definition: VectorSpace.H:111
Foam::radiation::reflectionModel::New
static autoPtr< reflectionModel > New(const dictionary &dict, const fvMesh &mesh)
addToRadiationRunTimeSelectionTables
#define addToRadiationRunTimeSelectionTables(model)
Definition: radiationModel.H:285
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::point
vector point
Point is a vector.
Definition: point.H:37
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:148
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Definition: UPstream.H:441
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:321
Foam::dimless
const dimensionSet dimless
Definition: dimensionSets.C:182
Foam::defineTemplateTypeNameAndDebugWithName
defineTemplateTypeNameAndDebugWithName(psiReactionsSensitivityAnalysisFunctionObject, "psiReactionsSensitivityAnalysis", 0)
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:170