Go to the documentation of this file.
29 #ifndef transmissionmatrix_h
30 #define transmissionmatrix_h
32 #include <Eigen/Dense>
33 #include <Eigen/StdVector>
41 std::vector<Eigen::Matrix4d, Eigen::aligned_allocator<Eigen::Matrix4d>>
T4;
42 std::vector<Eigen::Matrix3d, Eigen::aligned_allocator<Eigen::Matrix3d>>
T3;
43 std::vector<Eigen::Matrix2d, Eigen::aligned_allocator<Eigen::Matrix2d>>
T2;
44 std::vector<Eigen::Matrix<double, 1, 1>,
45 Eigen::aligned_allocator<Eigen::Matrix<double, 1, 1>>>
61 assert(stokes_dim < 5 and stokes_dim > 0);
73 T1(std::move(tm.T1)) {}
95 T4 = std::move(tm.T4);
96 T3 = std::move(tm.T3);
97 T2 = std::move(tm.T2);
98 T1 = std::move(tm.T1);
111 for (
size_t i = 0; i <
T4.size(); i++)
112 for (
size_t j = 0; j < 4; j++)
113 for (
size_t k = 0; k < 4; k++) T(i, j, k) =
T4[i](j, k);
114 for (
size_t i = 0; i <
T3.size(); i++)
115 for (
size_t j = 0; j < 3; j++)
116 for (
size_t k = 0; k < 3; k++) T(i, j, k) =
T3[i](j, k);
117 for (
size_t i = 0; i <
T2.size(); i++)
118 for (
size_t j = 0; j < 2; j++)
119 for (
size_t k = 0; k < 2; k++) T(i, j, k) =
T2[i](j, k);
120 for (
size_t i = 0; i <
T1.size(); i++) T(i, 0, 0) =
T1[i](0, 0);
129 const Eigen::Matrix4d&
Mat4(
size_t i)
const {
return T4[i]; }
136 const Eigen::Matrix3d&
Mat3(
size_t i)
const {
return T3[i]; }
143 const Eigen::Matrix2d&
Mat2(
size_t i)
const {
return T2[i]; }
150 const Eigen::Matrix<double, 1, 1>&
Mat1(
size_t i)
const {
return T1[i]; }
157 Eigen::MatrixXd
Mat(
size_t i)
const {
176 Eigen::Matrix4d&
Mat4(
size_t i) {
return T4[i]; }
183 Eigen::Matrix3d&
Mat3(
size_t i) {
return T3[i]; }
190 Eigen::Matrix2d&
Mat2(
size_t i) {
return T2[i]; }
197 Eigen::Matrix<double, 1, 1>&
Mat1(
size_t i) {
return T1[i]; }
209 std::fill(
T4.begin(),
T4.end(), Eigen::Matrix4d::Zero());
210 std::fill(
T3.begin(),
T3.end(), Eigen::Matrix3d::Zero());
211 std::fill(
T2.begin(),
T2.end(), Eigen::Matrix2d::Zero());
212 std::fill(
T1.begin(),
T1.end(), Eigen::Matrix<double, 1, 1>::Zero());
223 for (
size_t i = 0; i <
T4.size(); i++)
T4[i].noalias() = A.
T4[i] * B.
T4[i];
224 for (
size_t i = 0; i <
T3.size(); i++)
T3[i].noalias() = A.
T3[i] * B.
T3[i];
225 for (
size_t i = 0; i <
T2.size(); i++)
T2[i].noalias() = A.
T2[i] * B.
T2[i];
226 for (
size_t i = 0; i <
T1.size(); i++)
T1[i].noalias() = A.
T1[i] * B.
T1[i];
238 for (
size_t i = 0; i <
T4.size(); i++)
T4[i] = A.
T4[i] * B.
T4[i];
239 for (
size_t i = 0; i <
T3.size(); i++)
T3[i] = A.
T3[i] * B.
T3[i];
240 for (
size_t i = 0; i <
T2.size(); i++)
T2[i] = A.
T2[i] * B.
T2[i];
241 for (
size_t i = 0; i <
T1.size(); i++)
T1[i] = A.
T1[i] * B.
T1[i];
287 for (
size_t i = 0; i <
T4.size(); i++)
289 for (
size_t i = 0; i <
T3.size(); i++)
291 for (
size_t i = 0; i <
T2.size(); i++)
293 for (
size_t i = 0; i <
T1.size(); i++)
323 friend std::ostream&
operator<<(std::ostream& os,
344 static_assert (
N not_eq 0 and
N < 5,
"Bad size N");
349 const Numeric od = std::log(
this ->
operator()(i, 0, 0));
355 const Numeric od = std::log(
this ->
operator()(i, 0, 0));
361 const Numeric od = std::log(
this ->
operator()(i, 0, 0));
367 const Numeric od = std::log(
this ->
operator()(i, 0, 0));
398 std::vector<Eigen::Vector4d, Eigen::aligned_allocator<Eigen::Vector4d>>
R4;
399 std::vector<Eigen::Vector3d, Eigen::aligned_allocator<Eigen::Vector3d>>
R3;
400 std::vector<Eigen::Vector2d, Eigen::aligned_allocator<Eigen::Vector2d>>
R2;
401 std::vector<Eigen::Matrix<double, 1, 1>,
402 Eigen::aligned_allocator<Eigen::Matrix<double, 1, 1>>>
418 assert(stokes_dim < 5 and stokes_dim > 0);
427 R4(std::move(rv.R4)),
428 R3(std::move(rv.R3)),
429 R2(std::move(rv.R2)),
430 R1(std::move(rv.R1)) {}
452 R4 = std::move(rv.R4);
453 R3 = std::move(rv.R3);
454 R2 = std::move(rv.R2);
455 R1 = std::move(rv.R1);
464 for (
size_t i = 0; i <
R4.size(); i++)
R4[i] = T.
Mat4(i) *
R4[i];
465 for (
size_t i = 0; i <
R3.size(); i++)
R3[i] = T.
Mat3(i) *
R3[i];
466 for (
size_t i = 0; i <
R2.size(); i++)
R2[i] = T.
Mat2(i) *
R2[i];
467 for (
size_t i = 0; i <
R1.size(); i++)
R1[i] = T.
Mat1(i) *
R1[i];
477 R4[i].noalias() = Eigen::Vector4d::Zero();
480 R3[i].noalias() = Eigen::Vector3d::Zero();
483 R2[i].noalias() = Eigen::Vector2d::Zero();
496 const Eigen::Vector4d&
Vec4(
size_t i)
const {
return R4[i]; }
503 const Eigen::Vector3d&
Vec3(
size_t i)
const {
return R3[i]; }
510 const Eigen::Vector2d&
Vec2(
size_t i)
const {
return R2[i]; }
517 const Eigen::Matrix<double, 1, 1>&
Vec1(
size_t i)
const {
return R1[i]; }
524 Eigen::VectorXd
Vec(
size_t i)
const {
543 Eigen::Vector4d&
Vec4(
size_t i) {
return R4[i]; }
550 Eigen::Vector3d&
Vec3(
size_t i) {
return R3[i]; }
557 Eigen::Vector2d&
Vec2(
size_t i) {
return R2[i]; }
564 Eigen::Matrix<double, 1, 1>&
Vec1(
size_t i) {
return R1[i]; }
572 for (
size_t i = 0; i <
R4.size(); i++)
573 R4[i].noalias() -= 0.5 * (O1.
R4[i] + O2.
R4[i]);
574 for (
size_t i = 0; i <
R3.size(); i++)
575 R3[i].noalias() -= 0.5 * (O1.
R3[i] + O2.
R3[i]);
576 for (
size_t i = 0; i <
R2.size(); i++)
577 R2[i].noalias() -= 0.5 * (O1.
R2[i] + O2.
R2[i]);
578 for (
size_t i = 0; i <
R1.size(); i++)
579 R1[i].noalias() -= 0.5 * (O1.
R1[i] + O2.
R1[i]);
589 for (
size_t i = 0; i <
R4.size(); i++)
590 R4[i].noalias() += 0.5 * (O1.
R4[i] + O2.
R4[i]);
591 for (
size_t i = 0; i <
R3.size(); i++)
592 R3[i].noalias() += 0.5 * (O1.
R3[i] + O2.
R3[i]);
593 for (
size_t i = 0; i <
R2.size(); i++)
594 R2[i].noalias() += 0.5 * (O1.
R2[i] + O2.
R2[i]);
595 for (
size_t i = 0; i <
R1.size(); i++)
596 R1[i].noalias() += 0.5 * (O1.
R1[i] + O2.
R1[i]);
607 for (
size_t i = 0; i <
R4.size(); i++) {
609 R4[i].noalias() +=
w.far * far.
R4[i] +
w.close * close.
R4[i];
611 for (
size_t i = 0; i <
R3.size(); i++) {
613 R3[i].noalias() +=
w.far * far.
R3[i] +
w.close * close.
R3[i];
615 for (
size_t i = 0; i <
R2.size(); i++) {
617 R2[i].noalias() +=
w.far * far.
R2[i] +
w.close * close.
R2[i];
619 for (
size_t i = 0; i <
R1.size(); i++) {
621 R1[i].noalias() +=
w.far * far.
R1[i] +
w.close * close.
R1[i];
638 for (
size_t i = 0; i <
R4.size(); i++)
639 R4[i].noalias() += PiT.
Mat4(i) * (dT.
Mat4(i) * ImJ.
R4[i] + dJ.
R4[i] -
641 for (
size_t i = 0; i <
R3.size(); i++)
642 R3[i].noalias() += PiT.
Mat3(i) * (dT.
Mat3(i) * ImJ.
R3[i] + dJ.
R3[i] -
644 for (
size_t i = 0; i <
R2.size(); i++)
645 R2[i].noalias() += PiT.
Mat2(i) * (dT.
Mat2(i) * ImJ.
R2[i] + dJ.
R2[i] -
647 for (
size_t i = 0; i <
R1.size(); i++)
648 R1[i].noalias() += PiT.
Mat1(i) * (dT.
Mat1(i) * ImJ.
R1[i] + dJ.
R1[i] -
662 for (
size_t i = 0; i <
R4.size(); i++)
664 for (
size_t i = 0; i <
R3.size(); i++)
666 for (
size_t i = 0; i <
R2.size(); i++)
668 for (
size_t i = 0; i <
R1.size(); i++)
684 for (
size_t i = 0; i <
R4.size(); i++)
685 R4[i].noalias() += A.
Mat4(i) *
x.R4[i];
686 for (
size_t i = 0; i <
R3.size(); i++)
687 R3[i].noalias() += A.
Mat3(i) *
x.R3[i];
688 for (
size_t i = 0; i <
R2.size(); i++)
689 R2[i].noalias() += A.
Mat2(i) *
x.R2[i];
690 for (
size_t i = 0; i <
R1.size(); i++)
691 R1[i].noalias() += A.
Mat1(i) *
x.R1[i];
705 for (
size_t i = 0; i <
R4.size(); i++)
707 for (
size_t i = 0; i <
R3.size(); i++)
709 for (
size_t i = 0; i <
R2.size(); i++)
711 for (
size_t i = 0; i <
R1.size(); i++)
712 R4[i][0] = PiT.
Mat1(i)[0] *
727 for (
size_t i = 0; i <
R4.size(); i++)
729 for (
size_t i = 0; i <
R3.size(); i++)
731 for (
size_t i = 0; i <
R2.size(); i++)
733 for (
size_t i = 0; i <
R1.size(); i++)
748 for (
size_t i = 0; i <
R4.size(); i++)
750 for (
size_t i = 0; i <
R3.size(); i++)
752 for (
size_t i = 0; i <
R2.size(); i++)
754 for (
size_t i = 0; i <
R1.size(); i++)
765 for (
size_t i = 0; i <
R4.size(); i++) {
771 for (
size_t i = 0; i <
R3.size(); i++) {
776 for (
size_t i = 0; i <
R2.size(); i++) {
780 for (
size_t i = 0; i <
R1.size(); i++) {
811 for (
size_t i = 0; i <
R4.size(); i++)
812 for (
size_t j = 0; j < 4; j++)
M(i, j) =
R4[i](j);
813 for (
size_t i = 0; i <
R3.size(); i++)
814 for (
size_t j = 0; j < 3; j++)
M(i, j) =
R3[i](j);
815 for (
size_t i = 0; i <
R2.size(); i++)
816 for (
size_t j = 0; j < 2; j++)
M(i, j) =
R2[i](j);
817 for (
size_t i = 0; i <
R1.size(); i++)
M(i, 0) =
R1[i](0);
840 Eigen::Vector4d(a.
Kjj()[i], a.
K12()[i], a.
K13()[i], a.
K14()[i]) *
842 Eigen::Vector4d(S.
Kjj()[i], S.
K12()[i], S.
K13()[i], S.
K14()[i]);
845 Eigen::Vector4d(a.
Kjj()[i], a.
K12()[i], a.
K13()[i], a.
K14()[i]) *
851 Eigen::Vector3d(a.
Kjj()[i], a.
K12()[i], a.
K13()[i]) * B[i] +
852 Eigen::Vector3d(S.
Kjj()[i], S.
K12()[i], S.
K13()[i]);
855 Eigen::Vector3d(a.
Kjj()[i], a.
K12()[i], a.
K13()[i]) * B[i];
859 R2[i].noalias() = Eigen::Vector2d(a.
Kjj()[i], a.
K12()[i]) * B[i] +
860 Eigen::Vector2d(S.
Kjj()[i], S.
K12()[i]);
862 R2[i].noalias() = Eigen::Vector2d(a.
Kjj()[i], a.
K12()[i]) * B[i];
866 R1[i][0] = a.
Kjj()[i] * B[i] + S.
Kjj()[i];
868 R1[i][0] = a.
Kjj()[i] * B[i];
932 CommutativeTransmission,
1029 const Index temp_deriv_pos);
1092 #endif // transmissionmatrix_h
Class to help with hidden temporary variables for operations of type Numeric times Class.
void setBackscatterTransmissionDerivative(const RadiationVector &I0, const TransmissionMatrix &Tr, const TransmissionMatrix &Tf, const TransmissionMatrix &dZ)
Set this to backscatter transmission scatter derivative.
ArrayOfArrayOfTransmissionMatrix cumulative_backscatter_derivative(ConstTensor5View t, const ArrayOfMatrix &aom)
Accumulated backscatter derivative (???)
void transform(VectorView y, double(&my_func)(double), ConstVectorView x)
A generic transform function for vectors, which can be used to implement mathematical functions opera...
CumulativeTransmission
Intended to hold various ways to accumulate the transmission matrix.
TransmissionMatrix(TransmissionMatrix &&tm) noexcept
Construct a new Transmission Matrix object.
std::vector< Eigen::Matrix3d, Eigen::aligned_allocator< Eigen::Matrix3d > > T3
void add_avg(const RadiationVector &O1, const RadiationVector &O2)
Add the average of two other RadiationVector to *this.
Stokes vector is as Propagation matrix but only has 4 possible values.
const Eigen::Vector2d & Vec2(size_t i) const
Return Vector at position.
const Eigen::Vector4d & Vec4(size_t i) const
Return Vector at position.
const Numeric & operator()(const Index i, const Index j) const
Access operator.
bool IsEmpty() const
Asks if the class is empty.
friend std::ostream & operator<<(std::ostream &os, const TransmissionMatrix &tm)
Output operator.
ArrayOfRetrievalQuantity jacobian_quantities(Workspace &ws) noexcept
Complex w(Complex z) noexcept
The Faddeeva function.
TransmissionMatrix(Index nf=0, Index stokes=1)
Construct a new Transmission Matrix object.
void update_radiation_vector(RadiationVector &I, ArrayOfRadiationVector &dI1, ArrayOfRadiationVector &dI2, const RadiationVector &J1, const RadiationVector &J2, const ArrayOfRadiationVector &dJ1, const ArrayOfRadiationVector &dJ2, const TransmissionMatrix &T, const TransmissionMatrix &PiT, const ArrayOfTransmissionMatrix &dT1, const ArrayOfTransmissionMatrix &dT2, const RadiativeTransferSolver solver)
Update the Radiation Vector.
std::vector< Eigen::Vector3d, Eigen::aligned_allocator< Eigen::Vector3d > > R3
Index NumberOfAzimuthAngles() const
The number of azimuth angles of the propagation matrix.
RadiationVector & operator=(const RadiationVector &rv)=default
Assign old radiation vector to this.
void stepwise_transmission(TransmissionMatrix &T, ArrayOfTransmissionMatrix &dT1, ArrayOfTransmissionMatrix &dT2, const PropagationMatrix &K1, const PropagationMatrix &K2, const ArrayOfPropagationMatrix &dK1, const ArrayOfPropagationMatrix &dK2, const Numeric &r, const Numeric &dr_dtemp1, const Numeric &dr_dtemp2, const Index temp_deriv_pos)
Set the stepwise transmission matrix.
std::istream & operator>>(std::istream &is, TransmissionMatrix &tm)
Input operator.
Struct of far and close weights.
invlib::Matrix< ArtsMatrix > Matrix
invlib wrapper type for ARTS matrices.
Numeric operator()(const Index i, const Index j, const Index k) const
Access value in matrix.
Weights< N > linear_in_tau_weights(size_t) const noexcept
@ CommutativeTransmission
void SetZero(size_t i)
Set Radiation Vector to Zero at position.
RadiationVector(const RadiationVector &rv)=default
Construct a new Radiation Vector object.
Index Frequencies() const
Number of frequencies.
Eigen::Vector2d & Vec2(size_t i)
Return Vector at position.
RadiationVector & operator=(const ConstMatrixView &M)
Set *this from matrix.
Index jacobian_do(Workspace &ws) noexcept
std::ostream & operator<<(std::ostream &os, const TransmissionMatrix &tm)
Output operator.
G0 G2 FVC Y DV Numeric Numeric Numeric Zeeman LowerQuantumNumbers void * data
Eigen::Matrix2d & Mat2(size_t i)
Get Matrix at position.
Eigen::Matrix4d & Mat4(size_t i)
Get Matrix at position.
Eigen::Matrix< Numeric, int(N), int(N)> close
VectorView K13(const Index iz=0, const Index ia=0)
Vector view to K(0, 2) elements.
void stepwise_source(RadiationVector &J, ArrayOfRadiationVector &dJ, const PropagationMatrix &K, const StokesVector &a, const StokesVector &S, const ArrayOfPropagationMatrix &dK, const ArrayOfStokesVector &da, const ArrayOfStokesVector &dS, const ConstVectorView B, const ConstVectorView dB_dT, const ArrayOfRetrievalQuantity &jacobian_quantities, const bool &jacobian_do)
Set the stepwise source.
void addDerivTransmission(const TransmissionMatrix &PiT, const TransmissionMatrix &dT, const RadiationVector &I)
Add the transmission derivative to this.
void addDerivEmission(const TransmissionMatrix &PiT, const TransmissionMatrix &dT, const TransmissionMatrix &T, const RadiationVector &ImJ, const RadiationVector &dJ)
Add the emission derivative to this.
void setDerivReflection(const RadiationVector &I, const TransmissionMatrix &PiT, const TransmissionMatrix &Z, const TransmissionMatrix &dZ)
Sets *this to the reflection derivative.
This can be used to make arrays out of anything.
friend std::istream & operator>>(std::istream &data, RadiationVector &rv)
Input operator.
void setSource(const StokesVector &a, const ConstVectorView &B, const StokesVector &S, Index i)
Set this to source vector at position.
RadiativeTransferSolver
Intended to hold various forward solvers.
Eigen::Vector3d & Vec3(size_t i)
Return Vector at position.
BackscatterSolver
Intended to hold various backscatter solvers.
void mul(const TransmissionMatrix &A, const TransmissionMatrix &B)
Set this to a multiple of A by B.
Index NumberOfZenithAngles() const
The number of zenith angles of the propagation matrix.
TransmissionMatrix & operator=(const LazyScale< TransmissionMatrix > &lstm)
Assign lazily.
Index StokesDim() const
Stokes dimensionaility.
Array< ArrayOfRadiationVector > ArrayOfArrayOfRadiationVector
VectorView K12(const Index iz=0, const Index ia=0)
Vector view to K(0, 1) elements.
ArrayOfTransmissionMatrix cumulative_transmission(const ArrayOfTransmissionMatrix &T, const CumulativeTransmission type)
Accumulate the transmission matrix over all layers.
Array< ArrayOfTransmissionMatrix > ArrayOfArrayOfTransmissionMatrix
Eigen::Matrix3d & Mat3(size_t i)
Get Matrix at position.
Array< ArrayOfArrayOfTransmissionMatrix > ArrayOfArrayOfArrayOfTransmissionMatrix
TransmissionMatrix & operator*=(const Numeric &scale)
Scale self.
void set_backscatter_radiation_vector(ArrayOfRadiationVector &I, ArrayOfArrayOfArrayOfRadiationVector &dI, const RadiationVector &I_incoming, const ArrayOfTransmissionMatrix &T, const ArrayOfTransmissionMatrix &PiTf, const ArrayOfTransmissionMatrix &PiTr, const ArrayOfTransmissionMatrix &Z, const ArrayOfArrayOfTransmissionMatrix &dT1, const ArrayOfArrayOfTransmissionMatrix &dT2, const ArrayOfArrayOfTransmissionMatrix &dZ, const BackscatterSolver solver)
Set the backscatter radiation vector.
RadiationVector(Index nf=0, Index stokes=1)
Construct a new Radiation Vector object.
const base & bas
A const reference to a value.
void setZero()
Set to zero matrix.
Array< RadiationVector > ArrayOfRadiationVector
Eigen::Matrix< double, 1, 1 > & Vec1(size_t i)
Return Vector at position.
const Eigen::Matrix< double, 1, 1 > & Vec1(size_t i) const
Return Vector at position.
std::vector< Eigen::Vector4d, Eigen::aligned_allocator< Eigen::Vector4d > > R4
Routines for setting up the jacobian.
NUMERIC Numeric
The type to use for all floating point numbers.
const Eigen::Matrix3d & Mat3(size_t i) const
Get Matrix at position.
const Eigen::Vector3d & Vec3(size_t i) const
Return Vector at position.
std::vector< Eigen::Matrix< double, 1, 1 >, Eigen::aligned_allocator< Eigen::Matrix< double, 1, 1 > > > T1
void add_weighted(const TransmissionMatrix &T, const RadiationVector &close, const RadiationVector &far)
Add the weighted source of two RadiationVector to *this.
RadiationVector(RadiationVector &&rv) noexcept
Construct a new Radiation Vector object.
Eigen::Matrix< double, 1, 1 > & Mat1(size_t i)
Get Matrix at position.
void leftMul(const TransmissionMatrix &T)
Multiply radiation vector from the left.
Eigen::VectorXd Vec(size_t i) const
Return Vector at position by copy.
Stuff related to the propagation matrix.
void setIdentity()
Set to identity matrix.
VectorView K14(const Index iz=0, const Index ia=0)
Vector view to K(0, 3) elements.
void addMultiplied(const TransmissionMatrix &A, const RadiationVector &x)
Add multiply.
A constant view of a Matrix.
ArrayOfTransmissionMatrix cumulative_backscatter(ConstTensor5View t, ConstMatrixView m)
Accumulated backscatter (???)
invlib::MatrixIdentity< Matrix > Identity
Array< TransmissionMatrix > ArrayOfTransmissionMatrix
Radiation Vector for Stokes dimension 1-4.
const Numeric & scale
A const reference to a Numeric.
std::vector< Eigen::Matrix4d, Eigen::aligned_allocator< Eigen::Matrix4d > > T4
Index StokesDim() const
Get Stokes dimension.
std::vector< Eigen::Vector2d, Eigen::aligned_allocator< Eigen::Vector2d > > R2
const Eigen::Matrix< double, 1, 1 > & Mat1(size_t i) const
Get Matrix at position.
VectorView Kjj(const Index iz=0, const Index ia=0)
Vector view to diagonal elements.
std::vector< Eigen::Matrix2d, Eigen::aligned_allocator< Eigen::Matrix2d > > T2
void setBackscatterTransmission(const RadiationVector &I0, const TransmissionMatrix &Tr, const TransmissionMatrix &Tf, const TransmissionMatrix &Z)
Set this to backscatter transmission.
TransmissionMatrix & operator+=(const LazyScale< TransmissionMatrix > &lstm)
Assign to *this lazily.
TransmissionMatrix(const TransmissionMatrix &tm)=default
Construct a new Transmission Matrix object.
Eigen::Matrix< Numeric, int(N), int(N)> far
Array< ArrayOfArrayOfRadiationVector > ArrayOfArrayOfArrayOfRadiationVector
std::vector< Eigen::Matrix< double, 1, 1 >, Eigen::aligned_allocator< Eigen::Matrix< double, 1, 1 > > > R1
Vector x(Workspace &ws) noexcept
INDEX Index
The type to use for all integer numbers and indices.
void rem_avg(const RadiationVector &O1, const RadiationVector &O2)
Remove the average of two other RadiationVector from *this.
const Eigen::Matrix2d & Mat2(size_t i) const
Get Matrix at position.
friend std::istream & operator>>(std::istream &data, TransmissionMatrix &tm)
Input operator.
Index Frequencies() const
Get frequency count.
Eigen::MatrixXd Mat(size_t i) const
Get Matrix at position by copy.
TransmissionMatrix & operator=(TransmissionMatrix &&tm) noexcept
Move operator.
RadiationVector & operator=(RadiationVector &&rv) noexcept
Assign old radiation vector to this.
Class to keep track of Transmission Matrices for Stokes Dim 1-4.
Eigen::Vector4d & Vec4(size_t i)
Return Vector at position.
friend std::ostream & operator<<(std::ostream &os, const RadiationVector &rv)
Output operator.
A constant view of a Vector.
TransmissionMatrix & operator=(const TransmissionMatrix &tm)=default
Assignment operator.
void mul_aliased(const TransmissionMatrix &A, const TransmissionMatrix &B)
Set this to a multiple of A by B.
A constant view of a Tensor5.
LazyScale< TransmissionMatrix > operator*(const TransmissionMatrix &tm, const Numeric &x)
Lazy scale of Transmission Matrix.
const Eigen::Matrix4d & Mat4(size_t i) const
Get Matrix at position.