Go to the documentation of this file.
63 Foam::radiation::laserDTRM::powerDistNames_
65 { powerDistributionMode::pdGaussian,
"Gaussian" },
66 { powerDistributionMode::pdManual,
"manual" },
67 { powerDistributionMode::pdUniform,
"uniform" },
68 { powerDistributionMode::pdGaussianPeak,
"GaussianPeak" },
74 Foam::scalar Foam::radiation::laserDTRM::calculateIp(scalar r, scalar theta)
76 const scalar t = mesh_.time().value();
77 const scalar power = laserPower_->value(t);
94 return power*powerDistribution_()(theta, r);
105 <<
"Unhandled type " << powerDistNames_[mode_]
134 return gradAlphaf/(
mag(gradAlphaf)+ deltaN);
138 void Foam::radiation::laserDTRM::initialiseReflection()
140 if (
found(
"reflectionModel"))
142 dictTable modelDicts(lookup(
"reflectionModel"));
146 const phasePairKey&
key = iter.key();
159 if (reflections_.size())
161 reflectionSwitch_ =
true;
164 reflectionSwitch_ =
returnReduce(reflectionSwitch_, orOp<bool>());
169 void Foam::radiation::laserDTRM::initialise()
174 const scalar t = mesh_.time().value();
175 const vector lPosition = focalLaserPosition_->value(t);
179 <<
"Laser position : " << lPosition <<
nl
180 <<
"Laser direction : " << lDir <<
endl;
189 while (magr < VSMALL)
192 rArea = v - (v & lDir)*lDir;
198 scalar dr = focalLaserRadius_/ndr_;
201 nParticles_ = ndr_*ndTheta_;
207 I0_ = get<scalar>(
"I0");
208 sigma_ = get<scalar>(
"sigma");
213 sigma_ = get<scalar>(
"sigma");
218 powerDistribution_.reset
220 new interpolation2DTable<scalar>(*
this)
243 if (mesh_.nGeometricD() == 3)
245 for (label ri = 0; ri < ndr_; ri++)
247 scalar r1 = SMALL + dr*ri;
251 scalar rP = ((r1 + r2)/2);
254 vector localR = ((r1 + r2)/2)*rArea;
260 for (label thetai = 0; thetai < ndTheta_; thetai++)
262 scalar theta1 = theta0 + SMALL + dTheta*thetai;
264 scalar theta2 = theta1 + dTheta;
266 scalar thetaP = (theta1 + theta2)/2.0;
268 quaternion Q(lDir, thetaP);
271 vector initialPos = (Q.R() & localR);
274 vector finalPos = (Q.R() & finalR);
277 p0 = lPosition + initialPos;
280 p1 = lPosition + finalPos + (0.5*maxTrackLength_*lDir);
282 scalar Ip = calculateIp(rP, thetaP);
284 scalar dAi = (
sqr(r2) -
sqr(r1))*(theta2 - theta1)/2.0;
289 label cellI = mesh_.findCell(
p0);
295 new DTRMParticle(mesh_,
p0, p1, Ip, cellI, dAi, -1);
298 DTRMCloud_.addParticle(pPtr);
306 <<
"Cannot find owner cell for focalPoint at "
316 <<
"Current functionality limited to 3-D cases"
322 Info<<
"Seeding missed " << nMissed <<
" locations" <<
endl;
326 <<
"Total Power in the laser : " << power <<
nl
327 <<
"Total Area in the laser : " <<
area <<
nl
336 radiationModel(typeName,
T),
337 mode_(powerDistNames_.
get(
"mode", *this)),
338 DTRMCloud_(mesh_,
"DTRMCloud",
IDLList<DTRMParticle>()),
340 ndTheta_(
get<label>(
"nTheta")),
341 ndr_(
get<label>(
"nr")),
342 maxTrackLength_(mesh_.bounds().
mag()),
346 Function1<
point>::
New(
"focalLaserPosition", *this, &mesh_)
351 Function1<
vector>::
New(
"laserDirection", *this, &mesh_)
354 focalLaserRadius_(
get<scalar>(
"focalLaserRadius")),
357 getOrDefault<scalar>(
"qualityBeamLaser", 0)
362 laserPower_(Function1<scalar>::
New(
"laserPower", *this, &mesh_)),
363 powerDistribution_(),
365 reflectionSwitch_(false),
367 alphaCut_(getOrDefault<scalar>(
"alphaCut", 0.5)),
422 initialiseReflection();
428 Foam::radiation::laserDTRM::laserDTRM
430 const dictionary&
dict,
434 radiationModel(typeName,
dict,
T),
435 mode_(powerDistNames_.
get(
"mode", *this)),
436 DTRMCloud_(mesh_,
"DTRMCloud",
IDLList<DTRMParticle>()),
438 ndTheta_(
get<label>(
"nTheta")),
439 ndr_(
get<label>(
"nr")),
440 maxTrackLength_(mesh_.bounds().
mag()),
444 Function1<
point>::
New(
"focalLaserPosition", *this, &mesh_)
448 Function1<
vector>::
New(
"laserDirection", *this, &mesh_)
451 focalLaserRadius_(
get<scalar>(
"focalLaserRadius")),
454 getOrDefault<scalar>(
"qualityBeamLaser", 0)
459 laserPower_(Function1<scalar>::
New(
"laserPower", *this, &mesh_)),
460 powerDistribution_(),
462 reflectionSwitch_(false),
464 alphaCut_(getOrDefault<scalar>(
"alphaCut", 0.5)),
519 initialiseReflection();
544 tmp<volScalarField> treflectingCells
550 "reflectingCellsVol",
551 mesh_.time().timeName(),
563 tmp<volVectorField> tnHat
570 mesh_.time().timeName(),
585 a_ = absorptionEmission_->a();
586 e_ = absorptionEmission_->e();
587 E_ = absorptionEmission_->E();
589 const interpolationCell<scalar> aInterp(a_);
590 const interpolationCell<scalar> eInterp(e_);
591 const interpolationCell<scalar> EInterp(E_);
592 const interpolationCell<scalar> TInterp(T_);
594 labelField reflectingCells(mesh_.nCells(), -1);
596 UPtrList<reflectionModel> reflectionUPtr;
598 if (reflectionSwitch_)
600 reflectionUPtr.resize(reflections_.size());
602 label reflectionModelId(0);
605 reflectionModel& model = iter1()();
607 reflectionUPtr.set(reflectionModelId, &model);
638 &&
mag(nHatPhase[cellI]) > 0.99
642 reflectingCells[cellI] = reflectionModelId;
643 reflectingCellsVol[cellI] = reflectionModelId;
644 if (
mag(nHat[cellI]) == 0.0)
646 nHat[cellI] += nHatPhase[cellI];
654 interpolationCellPoint<vector> nHatInterp(nHat);
656 DTRMParticle::trackingData td
669 Info<<
"Move particles..."
672 DTRMCloud_.move(DTRMCloud_, td, mesh_.time().deltaTValue());
675 Q_.primitiveFieldRef() /= mesh_.V();
679 Info<<
"Final number of particles..."
682 OFstream osRef(
type() +
":particlePath.obj");
688 DynamicList<point> positionsMyProc;
689 DynamicList<point> p0MyProc;
691 for (
const DTRMParticle&
p : DTRMCloud_)
693 positionsMyProc.append(
p.position());
694 p0MyProc.append(
p.p0());
715 osRef <<
"l " << vertI-1 <<
' ' << vertI <<
nl;
721 scalar totalQ =
gSum(Q_.primitiveFieldRef()*mesh_.V());
722 Info <<
"Total energy absorbed [W]: " << totalQ <<
endl;
724 if (mesh_.time().outputTime())
726 reflectingCellsVol.
write();
745 mesh_.time().timeName(),
760 return Q_.internalField();
vectorField pointField
pointField is a vectorField.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
const volScalarField & alpha2
Different types of constants.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
ILList< DLListBase, T > IDLList
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Unit conversion functions.
Ostream & endl(Ostream &os)
Type gSum(const FieldField< Field, Type > &f)
const volScalarField & alpha1
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
Field< label > labelField
Specialisation of Field<T> for label.
dimensionedScalar pow3(const dimensionedScalar &ds)
Calculate the matrix for the laplacian of the field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define forAllIters(container, iter)
constexpr scalar twoPi(2 *M_PI)
const dimensionSet dimPower
virtual tmp< volScalarField > Rp() const
Macros for easy insertion into run-time selection tables.
virtual bool write(const token &tok)=0
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
errorManip< error > abort(error &err)
Vector< scalar > vector
A scalar version of the templated Vector.
Calculate the matrix for implicit and explicit sources.
errorManipArg< error, int > exit(error &err, const int errNo=1)
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual label nBands() const
#define FatalErrorInFunction
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
static rangeType allProcs(const label communicator=worldComm)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static int myProcNo(const label communicator=worldComm)
constexpr scalar pi(M_PI)
forAllConstIters(mixture.phases(), phase)
Internal & ref(const bool updateAccessTime=true)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
const dimensionedScalar e
static tmp< T > New(Args &&... args)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
static autoPtr< reflectionModel > New(const dictionary &dict, const fvMesh &mesh)
#define addToRadiationRunTimeSelectionTables(model)
static word groupName(StringType base, const word &group)
const volScalarField & p0
vector point
Point is a vector.
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
static label nProcs(const label communicator=worldComm)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimless
defineTemplateTypeNameAndDebugWithName(psiReactionsSensitivityAnalysisFunctionObject, "psiReactionsSensitivityAnalysis", 0)
dimensionedScalar pos(const dimensionedScalar &ds)