ARTS  2.3.1285(git:92a29ea9-dirty)
RadiationVector Class Reference

Radiation Vector for Stokes dimension 1-4. More...

#include <transmissionmatrix.h>

Public Member Functions

 RadiationVector (Index nf=0, Index stokes=1)
 Construct a new Radiation Vector object. More...
 
 RadiationVector (RadiationVector &&rv) noexcept
 Construct a new Radiation Vector object. More...
 
 RadiationVector (const RadiationVector &rv)=default
 Construct a new Radiation Vector object. More...
 
RadiationVectoroperator= (const RadiationVector &rv)=default
 Assign old radiation vector to this. More...
 
RadiationVectoroperator= (RadiationVector &&rv) noexcept
 Assign old radiation vector to this. More...
 
void leftMul (const TransmissionMatrix &T)
 Multiply radiation vector from the left. More...
 
void SetZero (size_t i)
 Set Radiation Vector to Zero at position. More...
 
const Eigen::Vector4d & Vec4 (size_t i) const
 Return Vector at position. More...
 
const Eigen::Vector3d & Vec3 (size_t i) const
 Return Vector at position. More...
 
const Eigen::Vector2d & Vec2 (size_t i) const
 Return Vector at position. More...
 
const Eigen::Matrix< double, 1, 1 > & Vec1 (size_t i) const
 Return Vector at position. More...
 
Eigen::VectorXd Vec (size_t i) const
 Return Vector at position by copy. More...
 
Eigen::Vector4d & Vec4 (size_t i)
 Return Vector at position. More...
 
Eigen::Vector3d & Vec3 (size_t i)
 Return Vector at position. More...
 
Eigen::Vector2d & Vec2 (size_t i)
 Return Vector at position. More...
 
Eigen::Matrix< double, 1, 1 > & Vec1 (size_t i)
 Return Vector at position. More...
 
void rem_avg (const RadiationVector &O1, const RadiationVector &O2)
 Remove the average of two other RadiationVector from *this. More...
 
void add_avg (const RadiationVector &O1, const RadiationVector &O2)
 Add the average of two other RadiationVector to *this. More...
 
void add_weighted (const TransmissionMatrix &T, const RadiationVector &close, const RadiationVector &far)
 Add the weighted source of two RadiationVector to *this. More...
 
void addDerivEmission (const TransmissionMatrix &PiT, const TransmissionMatrix &dT, const TransmissionMatrix &T, const RadiationVector &ImJ, const RadiationVector &dJ)
 Add the emission derivative to this. More...
 
void addDerivTransmission (const TransmissionMatrix &PiT, const TransmissionMatrix &dT, const RadiationVector &I)
 Add the transmission derivative to this. More...
 
void addMultiplied (const TransmissionMatrix &A, const RadiationVector &x)
 Add multiply. More...
 
void setDerivReflection (const RadiationVector &I, const TransmissionMatrix &PiT, const TransmissionMatrix &Z, const TransmissionMatrix &dZ)
 Sets *this to the reflection derivative. More...
 
void setBackscatterTransmission (const RadiationVector &I0, const TransmissionMatrix &Tr, const TransmissionMatrix &Tf, const TransmissionMatrix &Z)
 Set this to backscatter transmission. More...
 
void setBackscatterTransmissionDerivative (const RadiationVector &I0, const TransmissionMatrix &Tr, const TransmissionMatrix &Tf, const TransmissionMatrix &dZ)
 Set this to backscatter transmission scatter derivative. More...
 
RadiationVectoroperator= (const ConstMatrixView &M)
 Set *this from matrix. More...
 
const Numericoperator() (const Index i, const Index j) const
 Access operator. More...
 
 operator Matrix () const
 Convert *this to Matrix class. More...
 
void setSource (const StokesVector &a, const ConstVectorView &B, const StokesVector &S, Index i)
 Set this to source vector at position. More...
 
Index StokesDim () const
 Get Stokes dimension. More...
 
Index Frequencies () const
 Get frequency count. More...
 

Private Attributes

Index stokes_dim
 
std::vector< Eigen::Vector4d, Eigen::aligned_allocator< Eigen::Vector4d > > R4
 
std::vector< Eigen::Vector3d, Eigen::aligned_allocator< Eigen::Vector3d > > R3
 
std::vector< Eigen::Vector2d, Eigen::aligned_allocator< Eigen::Vector2d > > R2
 
std::vector< Eigen::Matrix< double, 1, 1 >, Eigen::aligned_allocator< Eigen::Matrix< double, 1, 1 > > > R1
 

Friends

std::ostream & operator<< (std::ostream &os, const RadiationVector &rv)
 Output operator. More...
 
std::istream & operator>> (std::istream &data, RadiationVector &rv)
 Input operator. More...
 

Detailed Description

Radiation Vector for Stokes dimension 1-4.

Definition at line 395 of file transmissionmatrix.h.

Constructor & Destructor Documentation

◆ RadiationVector() [1/3]

RadiationVector::RadiationVector ( Index  nf = 0,
Index  stokes = 1 
)
inline

Construct a new Radiation Vector object.

Parameters
[in]nfNumber of frequencies
[in]stokesStokes dimension

Definition at line 412 of file transmissionmatrix.h.

◆ RadiationVector() [2/3]

RadiationVector::RadiationVector ( RadiationVector &&  rv)
inlinenoexcept

Construct a new Radiation Vector object.

Parameters
[in]rvOld Radiation Vector to move from

Definition at line 425 of file transmissionmatrix.h.

References TransmissionMatrix::operator=().

◆ RadiationVector() [3/3]

RadiationVector::RadiationVector ( const RadiationVector rv)
default

Construct a new Radiation Vector object.

Parameters
[in]rvOld Vector to copy from

Member Function Documentation

◆ add_avg()

void RadiationVector::add_avg ( const RadiationVector O1,
const RadiationVector O2 
)
inline

Add the average of two other RadiationVector to *this.

Parameters
[in]O1Input 1
[in]O2Input 2

Definition at line 588 of file transmissionmatrix.h.

References i, R1, R2, R3, and R4.

Referenced by update_radiation_vector().

◆ add_weighted()

void RadiationVector::add_weighted ( const TransmissionMatrix T,
const RadiationVector close,
const RadiationVector far 
)
inline

Add the weighted source of two RadiationVector to *this.

Parameters
[in]TThe transmission matrix
[in]closeInput 1
[in]farInput 2

Definition at line 606 of file transmissionmatrix.h.

References i, TransmissionMatrix::linear_in_tau_weights(), R1, R2, R3, R4, and w().

Referenced by update_radiation_vector().

◆ addDerivEmission()

void RadiationVector::addDerivEmission ( const TransmissionMatrix PiT,
const TransmissionMatrix dT,
const TransmissionMatrix T,
const RadiationVector ImJ,
const RadiationVector dJ 
)
inline

Add the emission derivative to this.

Parameters
[in]PiTAccumulated transmission to space
[in]dTDerivative of transmission matrix
[in]TTransmission matrix
[in]ImJIntensity minus the emission vector
[in]dJDerivative of the emission vector

Definition at line 633 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), TransmissionMatrix::Mat4(), R1, R2, R3, and R4.

◆ addDerivTransmission()

void RadiationVector::addDerivTransmission ( const TransmissionMatrix PiT,
const TransmissionMatrix dT,
const RadiationVector I 
)
inline

Add the transmission derivative to this.

Parameters
[in]PiTAccumulated transmission to space
[in]dTDerivative of transmission matrix
[in]IIntensity vector

Definition at line 659 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), TransmissionMatrix::Mat4(), R1, R2, R3, and R4.

◆ addMultiplied()

void RadiationVector::addMultiplied ( const TransmissionMatrix A,
const RadiationVector x 
)
inline

Add multiply.

Performs essentially this += A * x

Assumes this and x are not aliases

Parameters
[in]ADerivative of transmission matrix
[in]xIntensity vector

Definition at line 682 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), TransmissionMatrix::Mat4(), R1, R2, R3, and R4.

◆ Frequencies()

Index RadiationVector::Frequencies ( ) const
inline

Get frequency count.

Definition at line 876 of file transmissionmatrix.h.

References data, TransmissionMatrix::operator<<, and TransmissionMatrix::operator>>.

Referenced by xml_write_to_stream().

◆ leftMul()

void RadiationVector::leftMul ( const TransmissionMatrix T)
inline

Multiply radiation vector from the left.

Parameters
[in]TTranmission Vector

Definition at line 463 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), and TransmissionMatrix::Mat4().

Referenced by update_radiation_vector().

◆ operator Matrix()

RadiationVector::operator Matrix ( ) const
inline

Convert *this to Matrix class.

Returns
Matrix

Definition at line 809 of file transmissionmatrix.h.

References TransmissionMatrix::Frequencies(), i, and M.

◆ operator()()

const Numeric& RadiationVector::operator() ( const Index  i,
const Index  j 
) const
inline

Access operator.

Parameters
[in]iPosition in outer vector
[in]jPosition in inner vector
Returns
const Numeric&

Definition at line 792 of file transmissionmatrix.h.

References i.

◆ operator=() [1/3]

RadiationVector& RadiationVector::operator= ( const RadiationVector rv)
default

Assign old radiation vector to this.

Parameters
[in]rvold vetor to copy
Returns
RadiationVector& *this

◆ operator=() [2/3]

RadiationVector& RadiationVector::operator= ( RadiationVector &&  rv)
inlinenoexcept

Assign old radiation vector to this.

Parameters
[in]rvold vetor to move from
Returns
RadiationVector& *this

Definition at line 450 of file transmissionmatrix.h.

References R1, R2, R3, R4, and stokes_dim.

◆ operator=() [3/3]

RadiationVector& RadiationVector::operator= ( const ConstMatrixView M)
inline

Set *this from matrix.

Parameters
[in]MMatrix
Returns
RadiationVector& *this

Definition at line 763 of file transmissionmatrix.h.

References TransmissionMatrix::Frequencies(), i, M, ConstMatrixView::ncols(), and ConstMatrixView::nrows().

◆ rem_avg()

void RadiationVector::rem_avg ( const RadiationVector O1,
const RadiationVector O2 
)
inline

Remove the average of two other RadiationVector from *this.

Parameters
[in]O1Input 1
[in]O2Input 2

Definition at line 571 of file transmissionmatrix.h.

References i, R1, R2, R3, and R4.

Referenced by update_radiation_vector().

◆ setBackscatterTransmission()

void RadiationVector::setBackscatterTransmission ( const RadiationVector I0,
const TransmissionMatrix Tr,
const TransmissionMatrix Tf,
const TransmissionMatrix Z 
)
inline

Set this to backscatter transmission.

Parameters
[in]I0Incoming intensity vector
[in]TrReflection transmission
[in]TfForward transmission
[in]ZReflection matrix

Definition at line 723 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), TransmissionMatrix::Mat4(), R1, R2, R3, and R4.

◆ setBackscatterTransmissionDerivative()

void RadiationVector::setBackscatterTransmissionDerivative ( const RadiationVector I0,
const TransmissionMatrix Tr,
const TransmissionMatrix Tf,
const TransmissionMatrix dZ 
)
inline

Set this to backscatter transmission scatter derivative.

Parameters
[in]I0Incoming intensity vector
[in]TrReflection transmission
[in]TfForward transmission
[in]dZReflection matrix derivative

Definition at line 744 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), TransmissionMatrix::Mat4(), R1, R2, R3, and R4.

◆ setDerivReflection()

void RadiationVector::setDerivReflection ( const RadiationVector I,
const TransmissionMatrix PiT,
const TransmissionMatrix Z,
const TransmissionMatrix dZ 
)
inline

Sets *this to the reflection derivative.

Parameters
[in]IIntensity vector
[in]PiTAccumulated transmission to space
[in]ZReflection matrix
[in]dZDerivative of reflection matrix

Definition at line 701 of file transmissionmatrix.h.

References i, TransmissionMatrix::Mat1(), TransmissionMatrix::Mat2(), TransmissionMatrix::Mat3(), TransmissionMatrix::Mat4(), R1, R2, R3, and R4.

◆ setSource()

void RadiationVector::setSource ( const StokesVector a,
const ConstVectorView B,
const StokesVector S,
Index  i 
)
inline

Set this to source vector at position.

Parameters
[in]aAbsorption vector
[in]BPlanck vector
[in]SScattering source vector
[in]iPosition

Definition at line 828 of file transmissionmatrix.h.

References i, PropagationMatrix::IsEmpty(), PropagationMatrix::K12(), PropagationMatrix::K13(), PropagationMatrix::K14(), PropagationMatrix::Kjj(), PropagationMatrix::NumberOfAzimuthAngles(), and PropagationMatrix::NumberOfZenithAngles().

Referenced by stepwise_source().

◆ SetZero()

void RadiationVector::SetZero ( size_t  i)
inline

Set Radiation Vector to Zero at position.

Parameters
[in]iposition

Definition at line 474 of file transmissionmatrix.h.

References i.

Referenced by stepwise_source().

◆ StokesDim()

Index RadiationVector::StokesDim ( ) const
inline

Get Stokes dimension.

Definition at line 873 of file transmissionmatrix.h.

References TransmissionMatrix::stokes_dim.

Referenced by stepwise_source(), and xml_write_to_stream().

◆ Vec()

Eigen::VectorXd RadiationVector::Vec ( size_t  i) const
inline

Return Vector at position by copy.

Parameters
[in]iposition
Returns
const Eigen::Matrix<double, 1, 1>& Vector

Definition at line 524 of file transmissionmatrix.h.

◆ Vec1() [1/2]

const Eigen::Matrix<double, 1, 1>& RadiationVector::Vec1 ( size_t  i) const
inline

Return Vector at position.

Parameters
[in]iposition
Returns
const Eigen::Matrix<double, 1, 1>& Vector

Definition at line 517 of file transmissionmatrix.h.

References i.

Referenced by stepwise_source().

◆ Vec1() [2/2]

Eigen::Matrix<double, 1, 1>& RadiationVector::Vec1 ( size_t  i)
inline

Return Vector at position.

Parameters
[in]iposition
Returns
Eigen::Matrix<double, 1, 1>& Vector

Definition at line 564 of file transmissionmatrix.h.

References i.

◆ Vec2() [1/2]

const Eigen::Vector2d& RadiationVector::Vec2 ( size_t  i) const
inline

Return Vector at position.

Parameters
[in]iposition
Returns
const Eigen::Vector2d& Vector

Definition at line 510 of file transmissionmatrix.h.

References i.

Referenced by stepwise_source().

◆ Vec2() [2/2]

Eigen::Vector2d& RadiationVector::Vec2 ( size_t  i)
inline

Return Vector at position.

Parameters
[in]iposition
Returns
Eigen::Vector2d& Vector

Definition at line 557 of file transmissionmatrix.h.

References i.

◆ Vec3() [1/2]

const Eigen::Vector3d& RadiationVector::Vec3 ( size_t  i) const
inline

Return Vector at position.

Parameters
[in]iposition
Returns
const Eigen::Vector3d& Vector

Definition at line 503 of file transmissionmatrix.h.

References i.

Referenced by stepwise_source(), and test_transmissionmatrix().

◆ Vec3() [2/2]

Eigen::Vector3d& RadiationVector::Vec3 ( size_t  i)
inline

Return Vector at position.

Parameters
[in]iposition
Returns
Eigen::Vector3d& Vector

Definition at line 550 of file transmissionmatrix.h.

References i.

◆ Vec4() [1/2]

const Eigen::Vector4d& RadiationVector::Vec4 ( size_t  i) const
inline

Return Vector at position.

Parameters
[in]iposition
Returns
const Eigen::Vector4d& Vector

Definition at line 496 of file transmissionmatrix.h.

References i.

Referenced by stepwise_source().

◆ Vec4() [2/2]

Eigen::Vector4d& RadiationVector::Vec4 ( size_t  i)
inline

Return Vector at position.

Parameters
[in]iposition
Returns
Eigen::Vector4d& Vector

Definition at line 543 of file transmissionmatrix.h.

References i.

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const RadiationVector rv 
)
friend

Output operator.

Definition at line 1824 of file transmissionmatrix.cc.

◆ operator>>

std::istream& operator>> ( std::istream &  data,
RadiationVector rv 
)
friend

Input operator.

Definition at line 1860 of file transmissionmatrix.cc.

Member Data Documentation

◆ R1

std::vector<Eigen::Matrix<double, 1, 1>, Eigen::aligned_allocator<Eigen::Matrix<double, 1, 1> > > RadiationVector::R1
private

◆ R2

std::vector<Eigen::Vector2d, Eigen::aligned_allocator<Eigen::Vector2d> > RadiationVector::R2
private

◆ R3

std::vector<Eigen::Vector3d, Eigen::aligned_allocator<Eigen::Vector3d> > RadiationVector::R3
private

◆ R4

std::vector<Eigen::Vector4d, Eigen::aligned_allocator<Eigen::Vector4d> > RadiationVector::R4
private

◆ stokes_dim

Index RadiationVector::stokes_dim
private

Definition at line 397 of file transmissionmatrix.h.

Referenced by operator=().


The documentation for this class was generated from the following file: