PLaSK library
Loading...
Searching...
No Matches
plask::optical::modal::ExpansionBessel Struct Referenceabstract

#include <solvers/optical/modal/bessel/expansioncyl.hpp>

Inheritance diagram for plask::optical::modal::ExpansionBessel:
[legend]
Collaboration diagram for plask::optical::modal::ExpansionBessel:
[legend]

Classes

struct  Integrals
 Matrices with computed integrals necessary to construct RE and RH matrices. More...
 
struct  Segment
 Integration segment data. More...
 

Public Member Functions

 ExpansionBessel (BesselSolverCyl *solver)
 Create new expansion.
 
virtual ~ExpansionBessel ()
 
void init1 ()
 Init expansion.
 
virtual void init2 ()=0
 Perform m-specific initialization.
 
void init3 ()
 Estimate required integration order.
 
virtual void reset ()
 Free allocated memory.
 
bool diagonalQE (size_t l) const override
 
size_t matrixSize () const override
 Return size of the expansion matrix (equal to the number of expansion coefficients)
 
void prepareField () override
 Prepare for computatiations of the fields.
 
void cleanupField () override
 Cleanup after computatiations of the fields.
 
LazyData< Vec< 3, dcomplex > > getField (size_t layer, const shared_ptr< const typename LevelsAdapter::Level > &level, const cvector &E, const cvector &H) override
 Compute electric or magnetic field on dst_mesh at certain level.
 
LazyData< Tensor3< dcomplex > > getMaterialEps (size_t layer, const shared_ptr< const typename LevelsAdapter::Level > &level, InterpolationMethod interp=INTERPOLATION_DEFAULT) override
 Get epsilons index back from expansion.
 
double integratePoyntingVert (const cvector &E, const cvector &H) override
 Compute vertical component of the Poynting vector for specified fields.
 
double integrateField (WhichField field, size_t layer, const cmatrix &TE, const cmatrix &TH, const std::function< std::pair< dcomplex, dcomplex >(size_t, size_t)> &vertical) override
 Compute ½ En·conj(Em) or ½ Hn·conj(Hm)
 
void beforeGetEpsilon () override
 Prepare retrieval of refractive index.
 
void afterGetEpsilon () override
 Finish retrieval of refractive index.
 
unsigned getM () const
 
void setM (unsigned n)
 
size_t idxs (size_t i)
 Get $ X_s $ index.
 
size_t idxp (size_t i)
 Get $ X_p $ index.
 
cmatrix epsV_k (size_t layer)
 
cmatrix epsTss (size_t layer)
 
cmatrix epsTsp (size_t layer)
 
cmatrix epsTps (size_t layer)
 
cmatrix epsTpp (size_t layer)
 
virtual std::vector< doublegetKpts ()
 Return expansion wavevectors.
 
- Public Member Functions inherited from plask::optical::modal::Expansion
 Expansion (ModalBase *solver)
 
virtual ~Expansion ()
 
TempMatrix getTempMatrix ()
 
double getLam0 () const
 Get lam0.
 
void setLam0 (double lam)
 Set lam0.
 
void clearLam0 ()
 Clear lam0.
 
dcomplex getK0 () const
 Get current k0.
 
void setK0 (dcomplex k)
 Set current k0.
 
void computeIntegrals ()
 Compute all expansion coefficients.
 
virtual bool diagonalQE (size_t PLASK_UNUSED(l)) const
 Tell if matrix for i-th layer is diagonal.
 
virtual void getMatrices (size_t layer, cmatrix &RE, cmatrix &RH)=0
 Get RE anf RH matrices.
 
virtual void getDiagonalEigenvectors (cmatrix &Te, cmatrix &Te1, const cmatrix &RE, const cdiagonal &gamma)
 Get eigenvectors with some physical meaning when the layer is diagonal.
 
void initField (WhichField which, InterpolationMethod method)
 Prepare for computations of the fields.
 
double getModeFlux (size_t n, const cmatrix &TE, const cmatrix &TH)
 Compute vertical component of the Poynting vector for fields specified as matrix column.
 

Public Attributes

int m
 Angular dependency index.
 
bool initialized
 Expansion is initialized.
 
bool m_changed
 m has changed and init2 must be called
 
OrderedAxis rbounds
 Horizontal axis with separate integration intervals.
 
std::vector< doublekpts
 Argument coefficients for Bessel expansion base (zeros of Bessel function for finite domain)
 
shared_ptr< RectangularMesh< 2 > > mesh
 Mesh for getting material data.
 
- Public Attributes inherited from plask::optical::modal::Expansion
WhichField which_field
 
InterpolationMethod field_interpolation
 
ModalBasesolver
 Solver which performs calculations (and is the interface to the outside world)
 
dcomplex k0
 Frequency for which the actual computations are performed.
 
double lam0
 Material parameters wavelength.
 
LazyData< doubletemperature
 Obtained temperature.
 
bool gain_connected
 Flag indicating if the gain is connected.
 
bool epsilon_connected
 Flag indicating if inEpsilon is connected.
 
LazyData< Tensor2< double > > gain
 Obtained gain.
 
LazyData< Tensor3< dcomplex > > epsilons
 Obtained epsilons.
 
LazyData< doublecarriers
 Carriers concentration.
 

Protected Member Functions

void beforeLayersIntegrals (dcomplex lam, dcomplex glam) override
 
Tensor3< dcomplex > getEps (size_t layer, size_t ri, double r, double matz, double lam, double glam)
 
void layerIntegrals (size_t layer, double lam, double glam) override
 Compute integrals for RE and RH matrices.
 
virtual void integrateParams (Integrals &integrals, const dcomplex *datap, const dcomplex *datar, const dcomplex *dataz, dcomplex datap0, dcomplex datar0, dcomplex dataz0)=0
 
virtual double fieldFactor (size_t i)=0
 
virtual cmatrix getHzMatrix (const cmatrix &Bz, cmatrix &Hz)=0
 
- Protected Member Functions inherited from plask::optical::modal::Expansion
virtual void beforeLayersIntegrals (dcomplex PLASK_UNUSED(lam), dcomplex PLASK_UNUSED(glam))
 Method called before layer integrals are computed.
 
virtual void afterLayersIntegrals ()
 Method called after layer integrals are computed.
 

Protected Attributes

std::vector< Segmentsegments
 Integration segments.
 
std::vector< booldiagonals
 Information if the layer is diagonal.
 
std::vector< Integralslayers_integrals
 Computed integrals.
 
- Protected Attributes inherited from plask::optical::modal::Expansion
TempMatrixPool temporary
 

Additional Inherited Members

- Public Types inherited from plask::optical::modal::Expansion
enum  Component { E_UNSPECIFIED = 0 , E_TRAN = 1 , E_LONG = 2 }
 Specified component in polarization or symmetry. More...
 

Detailed Description

Definition at line 27 of file expansioncyl.hpp.

Constructor & Destructor Documentation

◆ ExpansionBessel()

plask::optical::modal::ExpansionBessel::ExpansionBessel ( BesselSolverCyl solver)

Create new expansion.

Parameters
solversolver which performs calculations

Definition at line 32 of file expansioncyl.cpp.

◆ ~ExpansionBessel()

virtual plask::optical::modal::ExpansionBessel::~ExpansionBessel ( )
inlinevirtual

Definition at line 50 of file expansioncyl.hpp.

Member Function Documentation

◆ afterGetEpsilon()

void plask::optical::modal::ExpansionBessel::afterGetEpsilon ( )
overridevirtual

Finish retrieval of refractive index.

Reimplemented from plask::optical::modal::Expansion.

Definition at line 162 of file expansioncyl.cpp.

◆ beforeGetEpsilon()

void plask::optical::modal::ExpansionBessel::beforeGetEpsilon ( )
overridevirtual

Prepare retrieval of refractive index.

Reimplemented from plask::optical::modal::Expansion.

Definition at line 150 of file expansioncyl.cpp.

◆ beforeLayersIntegrals()

void plask::optical::modal::ExpansionBessel::beforeLayersIntegrals ( dcomplex  lam,
dcomplex  glam 
)
overrideprotected

Definition at line 145 of file expansioncyl.cpp.

◆ cleanupField()

void plask::optical::modal::ExpansionBessel::cleanupField ( )
overridevirtual

Cleanup after computatiations of the fields.

Reimplemented from plask::optical::modal::Expansion.

Definition at line 400 of file expansioncyl.cpp.

◆ diagonalQE()

bool plask::optical::modal::ExpansionBessel::diagonalQE ( size_t  l) const
inlineoverride

Definition at line 64 of file expansioncyl.hpp.

◆ epsTpp()

cmatrix plask::optical::modal::ExpansionBessel::epsTpp ( size_t  layer)

Definition at line 387 of file expansioncyl.cpp.

◆ epsTps()

cmatrix plask::optical::modal::ExpansionBessel::epsTps ( size_t  layer)

Definition at line 380 of file expansioncyl.cpp.

◆ epsTsp()

cmatrix plask::optical::modal::ExpansionBessel::epsTsp ( size_t  layer)

Definition at line 373 of file expansioncyl.cpp.

◆ epsTss()

cmatrix plask::optical::modal::ExpansionBessel::epsTss ( size_t  layer)

Definition at line 366 of file expansioncyl.cpp.

◆ epsV_k()

cmatrix plask::optical::modal::ExpansionBessel::epsV_k ( size_t  layer)

Definition at line 359 of file expansioncyl.cpp.

◆ fieldFactor()

virtual double plask::optical::modal::ExpansionBessel::fieldFactor ( size_t  i)
protectedpure virtual

◆ getEps()

Tensor3< dcomplex > plask::optical::modal::ExpansionBessel::getEps ( size_t  layer,
size_t  ri,
double  r,
double  matz,
double  lam,
double  glam 
)
protected

Definition at line 164 of file expansioncyl.cpp.

◆ getField()

LazyData< Vec< 3, dcomplex > > plask::optical::modal::ExpansionBessel::getField ( size_t  l,
const shared_ptr< const typename LevelsAdapter::Level > &  level,
const cvector E,
const cvector H 
)
overridevirtual

Compute electric or magnetic field on dst_mesh at certain level.

Parameters
llayer number
leveldestination level
E,Helectric and magnetic field coefficientscients
Returns
field distribution at dst_mesh
field distribution at dst_mesh

Implements plask::optical::modal::Expansion.

Definition at line 402 of file expansioncyl.cpp.

◆ getHzMatrix()

virtual cmatrix plask::optical::modal::ExpansionBessel::getHzMatrix ( const cmatrix Bz,
cmatrix Hz 
)
protectedpure virtual

◆ getKpts()

virtual std::vector< double > plask::optical::modal::ExpansionBessel::getKpts ( )
inlinevirtual

Return expansion wavevectors.

Definition at line 193 of file expansioncyl.hpp.

◆ getM()

unsigned plask::optical::modal::ExpansionBessel::getM ( ) const
inline

Definition at line 168 of file expansioncyl.hpp.

◆ getMaterialEps()

LazyData< Tensor3< dcomplex > > plask::optical::modal::ExpansionBessel::getMaterialEps ( size_t  lay,
const shared_ptr< const typename LevelsAdapter::Level > &  level,
InterpolationMethod  interp = INTERPOLATION_DEFAULT 
)
overridevirtual

Get epsilons index back from expansion.

Parameters
laylayer number
meshmesh to get parameters to
interpinterpolation method
Returns
computed refractive indices

Implements plask::optical::modal::Expansion.

Definition at line 461 of file expansioncyl.cpp.

◆ idxp()

size_t plask::optical::modal::ExpansionBessel::idxp ( size_t  i)
inline

Get $ X_p $ index.

Definition at line 182 of file expansioncyl.hpp.

◆ idxs()

size_t plask::optical::modal::ExpansionBessel::idxs ( size_t  i)
inline

Get $ X_s $ index.

Definition at line 179 of file expansioncyl.hpp.

◆ init1()

void plask::optical::modal::ExpansionBessel::init1 ( )

Init expansion.

Definition at line 36 of file expansioncyl.cpp.

◆ init2()

virtual void plask::optical::modal::ExpansionBessel::init2 ( )
pure virtual

Perform m-specific initialization.

Implemented in plask::optical::modal::ExpansionBesselFini, and plask::optical::modal::ExpansionBesselInfini.

◆ init3()

void plask::optical::modal::ExpansionBessel::init3 ( )

Estimate required integration order.

Definition at line 72 of file expansioncyl.cpp.

◆ integrateField()

double plask::optical::modal::ExpansionBessel::integrateField ( WhichField  field,
size_t  layer,
const cmatrix TE,
const cmatrix TH,
const std::function< std::pair< dcomplex, dcomplex >(size_t, size_t)> &  vertical 
)
overridevirtual

Compute ½ En·conj(Em) or ½ Hn·conj(Hm)

Parameters
fieldfield to integrate
layerlayer number
TEelectric field coefficients matrix
THmagnetic field coefficients matrix
[out]result½ E·conj(E) or ½ H·conj(H) matrix

Implements plask::optical::modal::Expansion.

Definition at line 501 of file expansioncyl.cpp.

◆ integrateParams()

virtual void plask::optical::modal::ExpansionBessel::integrateParams ( Integrals integrals,
const dcomplex *  datap,
const dcomplex *  datar,
const dcomplex *  dataz,
dcomplex  datap0,
dcomplex  datar0,
dcomplex  dataz0 
)
protectedpure virtual

◆ integratePoyntingVert()

double plask::optical::modal::ExpansionBessel::integratePoyntingVert ( const cvector E,
const cvector H 
)
overridevirtual

Compute vertical component of the Poynting vector for specified fields.

Parameters
Eelectric field coefficients vector
Hmagnetic field coefficients vector
Returns
integrated Poynting vector i.e. the total vertically emitted energy

Implements plask::optical::modal::Expansion.

Definition at line 491 of file expansioncyl.cpp.

◆ layerIntegrals()

void plask::optical::modal::ExpansionBessel::layerIntegrals ( size_t  layer,
double  lam,
double  glam 
)
overrideprotectedvirtual

Compute integrals for RE and RH matrices.

Parameters
layerlayer number
lamwavelength
glamwavelength for gain

Implements plask::optical::modal::Expansion.

Definition at line 218 of file expansioncyl.cpp.

◆ matrixSize()

size_t plask::optical::modal::ExpansionBessel::matrixSize ( ) const
overridevirtual

Return size of the expansion matrix (equal to the number of expansion coefficients)

Returns
size of the expansion matrix

Implements plask::optical::modal::Expansion.

Definition at line 34 of file expansioncyl.cpp.

◆ prepareField()

void plask::optical::modal::ExpansionBessel::prepareField ( )
overridevirtual

Prepare for computatiations of the fields.

Reimplemented from plask::optical::modal::Expansion.

Definition at line 396 of file expansioncyl.cpp.

◆ reset()

void plask::optical::modal::ExpansionBessel::reset ( )
virtual

Free allocated memory.

Reimplemented in plask::optical::modal::ExpansionBesselFini.

Definition at line 63 of file expansioncyl.cpp.

◆ setM()

void plask::optical::modal::ExpansionBessel::setM ( unsigned  n)
inline

Definition at line 169 of file expansioncyl.hpp.

Member Data Documentation

◆ diagonals

std::vector<bool> plask::optical::modal::ExpansionBessel::diagonals
protected

Information if the layer is diagonal.

Definition at line 114 of file expansioncyl.hpp.

◆ initialized

bool plask::optical::modal::ExpansionBessel::initialized

Expansion is initialized.

Definition at line 30 of file expansioncyl.hpp.

◆ kpts

std::vector<double> plask::optical::modal::ExpansionBessel::kpts

Argument coefficients for Bessel expansion base (zeros of Bessel function for finite domain)

Definition at line 39 of file expansioncyl.hpp.

◆ layers_integrals

std::vector<Integrals> plask::optical::modal::ExpansionBessel::layers_integrals
protected

Computed integrals.

Definition at line 146 of file expansioncyl.hpp.

◆ m

int plask::optical::modal::ExpansionBessel::m

Angular dependency index.

Definition at line 28 of file expansioncyl.hpp.

◆ m_changed

bool plask::optical::modal::ExpansionBessel::m_changed

m has changed and init2 must be called

Definition at line 32 of file expansioncyl.hpp.

◆ mesh

shared_ptr<RectangularMesh<2> > plask::optical::modal::ExpansionBessel::mesh

Mesh for getting material data.

Definition at line 42 of file expansioncyl.hpp.

◆ rbounds

OrderedAxis plask::optical::modal::ExpansionBessel::rbounds

Horizontal axis with separate integration intervals.

material functions contain discontinuities at these points

Definition at line 36 of file expansioncyl.hpp.

◆ segments

std::vector<Segment> plask::optical::modal::ExpansionBessel::segments
protected

Integration segments.

Definition at line 111 of file expansioncyl.hpp.


The documentation for this struct was generated from the following files: