ARTS  2.4.0(git:4fb77825)
transmissionmatrix.h
Go to the documentation of this file.
1 /* Copyright (C) 2018
2  * Richard Larsson <ric.larsson@gmail.com>
3  *
4  * This program is free software; you can redistribute it and/or modify it
5  * under the terms of the GNU General Public License as published by the
6  * Free Software Foundation; either version 2, or (at your option) any
7  * later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program; if not, write to the Free Software
16  * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
17  * USA. */
18 
29 #ifndef transmissionmatrix_h
30 #define transmissionmatrix_h
31 
32 #include <Eigen/Dense>
33 #include <Eigen/StdVector>
34 #include "jacobian.h"
35 #include "propagationmatrix.h"
36 
39  private:
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>>>
46  T1;
47 
48  public:
49 
55  TransmissionMatrix(Index nf = 0, Index stokes = 1)
56  : stokes_dim(stokes),
57  T4(stokes_dim == 4 ? nf : 0, Eigen::Matrix4d::Identity()),
58  T3(stokes_dim == 3 ? nf : 0, Eigen::Matrix3d::Identity()),
59  T2(stokes_dim == 2 ? nf : 0, Eigen::Matrix2d::Identity()),
60  T1(stokes_dim == 1 ? nf : 0, Eigen::Matrix<double, 1, 1>::Identity()) {
61  assert(stokes_dim < 5 and stokes_dim > 0);
62  }
63 
69  : stokes_dim(std::move(tm.stokes_dim)),
70  T4(std::move(tm.T4)),
71  T3(std::move(tm.T3)),
72  T2(std::move(tm.T2)),
73  T1(std::move(tm.T1)) {}
74 
79  TransmissionMatrix(const TransmissionMatrix& tm) = default;
80 
87 
94  stokes_dim = std::move(tm.stokes_dim);
95  T4 = std::move(tm.T4);
96  T3 = std::move(tm.T3);
97  T2 = std::move(tm.T2);
98  T1 = std::move(tm.T1);
99  return *this;
100  }
101 
107  explicit TransmissionMatrix(const PropagationMatrix& pm, const Numeric& r = 1.0);
108 
109  operator Tensor3() const {
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);
121  return T;
122  }
123 
129  const Eigen::Matrix4d& Mat4(size_t i) const { return T4[i]; }
130 
136  const Eigen::Matrix3d& Mat3(size_t i) const { return T3[i]; }
137 
143  const Eigen::Matrix2d& Mat2(size_t i) const { return T2[i]; }
144 
150  const Eigen::Matrix<double, 1, 1>& Mat1(size_t i) const { return T1[i]; }
151 
157  Eigen::MatrixXd Mat(size_t i) const {
158  switch (stokes_dim) {
159  case 1:
160  return Mat1(i);
161  case 2:
162  return Mat2(i);
163  case 3:
164  return Mat3(i);
165  case 4:
166  return Mat4(i);
167  }
168  std::terminate();
169  }
170 
176  Eigen::Matrix4d& Mat4(size_t i) { return T4[i]; }
177 
183  Eigen::Matrix3d& Mat3(size_t i) { return T3[i]; }
184 
190  Eigen::Matrix2d& Mat2(size_t i) { return T2[i]; }
191 
197  Eigen::Matrix<double, 1, 1>& Mat1(size_t i) { return T1[i]; }
198 
200  void setIdentity() {
201  std::fill(T4.begin(), T4.end(), Eigen::Matrix4d::Identity());
202  std::fill(T3.begin(), T3.end(), Eigen::Matrix3d::Identity());
203  std::fill(T2.begin(), T2.end(), Eigen::Matrix2d::Identity());
204  std::fill(T1.begin(), T1.end(), Eigen::Matrix<double, 1, 1>::Identity());
205  }
206 
208  void setZero() {
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());
213  }
214 
222  void mul(const TransmissionMatrix& A, const TransmissionMatrix& B) {
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];
227  }
228 
229 
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];
242  }
243 
251  Numeric operator()(const Index i, const Index j, const Index k) const {
252  switch (stokes_dim) {
253  case 4:
254  return T4[i](j, k);
255  case 3:
256  return T3[i](j, k);
257  case 2:
258  return T2[i](j, k);
259  default:
260  return T1[i](j, k);
261  }
262  }
263 
265  Index StokesDim() const { return stokes_dim; }
266 
268  Index Frequencies() const {
269  switch (stokes_dim) {
270  case 4:
271  return Index(T4.size());
272  case 3:
273  return Index(T3.size());
274  case 2:
275  return Index(T2.size());
276  default:
277  return Index(T1.size());
278  }
279  }
280 
287  for (size_t i = 0; i < T4.size(); i++)
288  T4[i].noalias() = lstm.scale * lstm.bas.Mat4(i);
289  for (size_t i = 0; i < T3.size(); i++)
290  T3[i].noalias() = lstm.scale * lstm.bas.Mat3(i);
291  for (size_t i = 0; i < T2.size(); i++)
292  T2[i].noalias() = lstm.scale * lstm.bas.Mat2(i);
293  for (size_t i = 0; i < T1.size(); i++)
294  T1[i].noalias() = lstm.scale * lstm.bas.Mat1(i);
295  return *this;
296  }
297 
304  std::transform(T4.begin(), T4.end(), T4.begin(), [scale](auto& T){return T*scale;});
305  std::transform(T3.begin(), T3.end(), T3.begin(), [scale](auto& T){return T*scale;});
306  std::transform(T2.begin(), T2.end(), T2.begin(), [scale](auto& T){return T*scale;});
307  std::transform(T1.begin(), T1.end(), T1.begin(), [scale](auto& T){return T*scale;});
308  return *this;
309  }
310 
317  operator=(lstm.bas);
318  operator*=(lstm.scale);
319  return *this;
320  }
321 
323  friend std::ostream& operator<<(std::ostream& os,
324  const TransmissionMatrix& tm);
325 
327  friend std::istream& operator>>(std::istream& data, TransmissionMatrix& tm);
328 
330  template <size_t N>
331  struct Weights {Eigen::Matrix<Numeric, int(N), int(N)> far, close;};
332 
341  template <size_t N>
342  Weights<N> linear_in_tau_weights(size_t) const noexcept {
343  // FIXME: Make the complete implementation here by using "if constexpr" in C++17
344  static_assert (N not_eq 0 and N < 5, "Bad size N");
345  }
346 };
347 
349  const Numeric od = std::log(this -> operator()(i, 0, 0));
352 }
353 
354 template <> inline TransmissionMatrix::Weights<2> TransmissionMatrix::linear_in_tau_weights(size_t i) const noexcept {
355  const Numeric od = std::log(this -> operator()(i, 0, 0));
358 }
359 
360 template <> inline TransmissionMatrix::Weights<3> TransmissionMatrix::linear_in_tau_weights(size_t i) const noexcept {
361  const Numeric od = std::log(this -> operator()(i, 0, 0));
364 }
365 
366 template <> inline TransmissionMatrix::Weights<4> TransmissionMatrix::linear_in_tau_weights(size_t i) const noexcept {
367  const Numeric od = std::log(this -> operator()(i, 0, 0));
370 }
371 
379  const Numeric& x) {
380  return LazyScale<TransmissionMatrix>(tm, x);
381 }
382 
390  const TransmissionMatrix& tm) {
391  return LazyScale<TransmissionMatrix>(tm, x);
392 }
393 
396  private:
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>>>
403  R1;
404 
405  public:
406 
412  RadiationVector(Index nf = 0, Index stokes = 1)
413  : stokes_dim(stokes),
414  R4(stokes_dim == 4 ? nf : 0, Eigen::Vector4d::Zero()),
415  R3(stokes_dim == 3 ? nf : 0, Eigen::Vector3d::Zero()),
416  R2(stokes_dim == 2 ? nf : 0, Eigen::Vector2d::Zero()),
417  R1(stokes_dim == 1 ? nf : 0, Eigen::Matrix<double, 1, 1>::Zero()) {
418  assert(stokes_dim < 5 and stokes_dim > 0);
419  }
420 
426  : stokes_dim(std::move(rv.stokes_dim)),
427  R4(std::move(rv.R4)),
428  R3(std::move(rv.R3)),
429  R2(std::move(rv.R2)),
430  R1(std::move(rv.R1)) {}
431 
436  RadiationVector(const RadiationVector& rv) = default;
437 
444 
451  stokes_dim = std::move(rv.stokes_dim);
452  R4 = std::move(rv.R4);
453  R3 = std::move(rv.R3);
454  R2 = std::move(rv.R2);
455  R1 = std::move(rv.R1);
456  return *this;
457  }
458 
463  void leftMul(const TransmissionMatrix& T) {
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];
468  }
469 
474  void SetZero(size_t i) {
475  switch (stokes_dim) {
476  case 4:
477  R4[i].noalias() = Eigen::Vector4d::Zero();
478  break;
479  case 3:
480  R3[i].noalias() = Eigen::Vector3d::Zero();
481  break;
482  case 2:
483  R2[i].noalias() = Eigen::Vector2d::Zero();
484  break;
485  case 1:
486  R1[i][0] = 0;
487  break;
488  }
489  }
490 
496  const Eigen::Vector4d& Vec4(size_t i) const { return R4[i]; }
497 
503  const Eigen::Vector3d& Vec3(size_t i) const { return R3[i]; }
504 
510  const Eigen::Vector2d& Vec2(size_t i) const { return R2[i]; }
511 
517  const Eigen::Matrix<double, 1, 1>& Vec1(size_t i) const { return R1[i]; }
518 
524  Eigen::VectorXd Vec(size_t i) const {
525  switch (stokes_dim) {
526  case 4:
527  return Vec4(i);
528  case 3:
529  return Vec3(i);
530  case 2:
531  return Vec2(i);
532  default:
533  return Vec1(i);
534  }
535  }
536 
537 
543  Eigen::Vector4d& Vec4(size_t i) { return R4[i]; }
544 
550  Eigen::Vector3d& Vec3(size_t i) { return R3[i]; }
551 
557  Eigen::Vector2d& Vec2(size_t i) { return R2[i]; }
558 
564  Eigen::Matrix<double, 1, 1>& Vec1(size_t i) { return R1[i]; }
565 
571  void rem_avg(const RadiationVector& O1, const RadiationVector& O2) {
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]);
580  }
581 
582 
588  void add_avg(const RadiationVector& O1, const RadiationVector& O2) {
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]);
597  }
598 
599 
606  void add_weighted(const TransmissionMatrix& T, const RadiationVector& close, const RadiationVector& far) {
607  for (size_t i = 0; i < R4.size(); i++) {
608  const auto w = T.linear_in_tau_weights<4>(i);
609  R4[i].noalias() += w.far * far.R4[i] + w.close * close.R4[i];
610  }
611  for (size_t i = 0; i < R3.size(); i++) {
612  const auto w = T.linear_in_tau_weights<3>(i);
613  R3[i].noalias() += w.far * far.R3[i] + w.close * close.R3[i];
614  }
615  for (size_t i = 0; i < R2.size(); i++) {
616  const auto w = T.linear_in_tau_weights<2>(i);
617  R2[i].noalias() += w.far * far.R2[i] + w.close * close.R2[i];
618  }
619  for (size_t i = 0; i < R1.size(); i++) {
620  const auto w = T.linear_in_tau_weights<1>(i);
621  R1[i].noalias() += w.far * far.R1[i] + w.close * close.R1[i];
622  }
623  }
624 
634  const TransmissionMatrix& dT,
635  const TransmissionMatrix& T,
636  const RadiationVector& ImJ,
637  const RadiationVector& dJ) {
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] -
640  T.Mat4(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] -
643  T.Mat3(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] -
646  T.Mat2(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] -
649  T.Mat1(i) * dJ.R1[i]);
650  }
651 
652 
660  const TransmissionMatrix& dT,
661  const RadiationVector& I) {
662  for (size_t i = 0; i < R4.size(); i++)
663  R4[i].noalias() += PiT.Mat4(i) * dT.Mat4(i) * I.R4[i];
664  for (size_t i = 0; i < R3.size(); i++)
665  R3[i].noalias() += PiT.Mat3(i) * dT.Mat3(i) * I.R3[i];
666  for (size_t i = 0; i < R2.size(); i++)
667  R2[i].noalias() += PiT.Mat2(i) * dT.Mat2(i) * I.R2[i];
668  for (size_t i = 0; i < R1.size(); i++)
669  R1[i].noalias() += PiT.Mat1(i) * dT.Mat1(i) * I.R1[i];
670  }
671 
672 
683  const RadiationVector& x) {
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];
692  }
693 
702  const TransmissionMatrix& PiT,
703  const TransmissionMatrix& Z,
704  const TransmissionMatrix& dZ) {
705  for (size_t i = 0; i < R4.size(); i++)
706  R4[i] = PiT.Mat4(i) * (Z.Mat4(i) * R4[i] + dZ.Mat4(i) * I.R4[i]);
707  for (size_t i = 0; i < R3.size(); i++)
708  R3[i] = PiT.Mat3(i) * (Z.Mat3(i) * R3[i] + dZ.Mat3(i) * I.R3[i]);
709  for (size_t i = 0; i < R2.size(); i++)
710  R2[i] = PiT.Mat2(i) * (Z.Mat2(i) * R2[i] + dZ.Mat2(i) * I.R2[i]);
711  for (size_t i = 0; i < R1.size(); i++)
712  R4[i][0] = PiT.Mat1(i)[0] *
713  (Z.Mat1(i)[0] * R1[i][0] + dZ.Mat1(i)[0] * I.R1[i][0]);
714  }
715 
724  const TransmissionMatrix& Tr,
725  const TransmissionMatrix& Tf,
726  const TransmissionMatrix& Z) {
727  for (size_t i = 0; i < R4.size(); i++)
728  R4[i].noalias() = Tr.Mat4(i) * Z.Mat4(i) * Tf.Mat4(i) * I0.R4[i];
729  for (size_t i = 0; i < R3.size(); i++)
730  R3[i].noalias() = Tr.Mat3(i) * Z.Mat3(i) * Tf.Mat3(i) * I0.R3[i];
731  for (size_t i = 0; i < R2.size(); i++)
732  R2[i].noalias() = Tr.Mat2(i) * Z.Mat2(i) * Tf.Mat2(i) * I0.R2[i];
733  for (size_t i = 0; i < R1.size(); i++)
734  R1[i].noalias() = Tr.Mat1(i) * Z.Mat1(i) * Tf.Mat1(i) * I0.R1[i];
735  }
736 
745  const TransmissionMatrix& Tr,
746  const TransmissionMatrix& Tf,
747  const TransmissionMatrix& dZ) {
748  for (size_t i = 0; i < R4.size(); i++)
749  R4[i].noalias() += Tr.Mat4(i) * dZ.Mat4(i) * Tf.Mat4(i) * I0.R4[i];
750  for (size_t i = 0; i < R3.size(); i++)
751  R3[i].noalias() += Tr.Mat3(i) * dZ.Mat3(i) * Tf.Mat3(i) * I0.R3[i];
752  for (size_t i = 0; i < R2.size(); i++)
753  R2[i].noalias() += Tr.Mat2(i) * dZ.Mat2(i) * Tf.Mat2(i) * I0.R2[i];
754  for (size_t i = 0; i < R1.size(); i++)
755  R1[i].noalias() += Tr.Mat1(i) * dZ.Mat1(i) * Tf.Mat1(i) * I0.R1[i];
756  }
757 
764  assert(M.ncols() == stokes_dim and M.nrows() == Frequencies());
765  for (size_t i = 0; i < R4.size(); i++) {
766  R4[i][0] = M(i, 0);
767  R4[i][1] = M(i, 1);
768  R4[i][2] = M(i, 2);
769  R4[i][3] = M(i, 3);
770  }
771  for (size_t i = 0; i < R3.size(); i++) {
772  R3[i][0] = M(i, 0);
773  R3[i][1] = M(i, 1);
774  R3[i][2] = M(i, 2);
775  }
776  for (size_t i = 0; i < R2.size(); i++) {
777  R2[i][0] = M(i, 0);
778  R2[i][1] = M(i, 1);
779  }
780  for (size_t i = 0; i < R1.size(); i++) {
781  R1[i][0] = M(i, 0);
782  }
783  return *this;
784  }
785 
792  const Numeric& operator()(const Index i, const Index j) const {
793  switch (stokes_dim) {
794  case 4:
795  return R4[i][j];
796  case 3:
797  return R3[i][j];
798  case 2:
799  return R2[i][j];
800  default:
801  return R1[i][j];
802  }
803  }
804 
809  operator Matrix() const {
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);
818  return M;
819  }
820 
828  void setSource(const StokesVector& a,
829  const ConstVectorView& B,
830  const StokesVector& S,
831  Index i) {
832  assert(a.NumberOfAzimuthAngles() == 1);
833  assert(a.NumberOfZenithAngles() == 1);
834  assert(S.NumberOfAzimuthAngles() == 1);
835  assert(S.NumberOfZenithAngles() == 1);
836  switch (stokes_dim) {
837  case 4:
838  if (not S.IsEmpty())
839  R4[i].noalias() =
840  Eigen::Vector4d(a.Kjj()[i], a.K12()[i], a.K13()[i], a.K14()[i]) *
841  B[i] +
842  Eigen::Vector4d(S.Kjj()[i], S.K12()[i], S.K13()[i], S.K14()[i]);
843  else
844  R4[i].noalias() =
845  Eigen::Vector4d(a.Kjj()[i], a.K12()[i], a.K13()[i], a.K14()[i]) *
846  B[i];
847  break;
848  case 3:
849  if (not S.IsEmpty())
850  R3[i].noalias() =
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]);
853  else
854  R3[i].noalias() =
855  Eigen::Vector3d(a.Kjj()[i], a.K12()[i], a.K13()[i]) * B[i];
856  break;
857  case 2:
858  if (not S.IsEmpty())
859  R2[i].noalias() = Eigen::Vector2d(a.Kjj()[i], a.K12()[i]) * B[i] +
860  Eigen::Vector2d(S.Kjj()[i], S.K12()[i]);
861  else
862  R2[i].noalias() = Eigen::Vector2d(a.Kjj()[i], a.K12()[i]) * B[i];
863  break;
864  default:
865  if (not S.IsEmpty())
866  R1[i][0] = a.Kjj()[i] * B[i] + S.Kjj()[i];
867  else
868  R1[i][0] = a.Kjj()[i] * B[i];
869  }
870  }
871 
873  Index StokesDim() const { return stokes_dim; }
874 
876  Index Frequencies() const {
877  switch (stokes_dim) {
878  case 4:
879  return Index(R4.size());
880  case 3:
881  return Index(R3.size());
882  case 2:
883  return Index(R2.size());
884  default:
885  return Index(R1.size());
886  }
887  }
888 
890  friend std::ostream& operator<<(std::ostream& os, const RadiationVector& rv);
891 
893  friend std::istream& operator>>(std::istream& data, RadiationVector& rv);
894 };
895 
902 
904 std::ostream& operator<<(std::ostream& os, const TransmissionMatrix& tm);
905 
907 std::ostream& operator<<(std::ostream& os,
908  const ArrayOfTransmissionMatrix& atm);
909 
911 std::ostream& operator<<(std::ostream& os,
913 
915 std::ostream& operator<<(std::ostream& os, const RadiationVector& rv);
916 
918 std::ostream& operator<<(std::ostream& os, const ArrayOfRadiationVector& arv);
919 
921 std::ostream& operator<<(std::ostream& os,
922  const ArrayOfArrayOfRadiationVector& aarv);
923 
925 std::istream& operator>>(std::istream& is, TransmissionMatrix& tm);
926 
928 std::istream& operator>>(std::istream& is, RadiationVector& rv);
929 
931 enum class BackscatterSolver {
932  CommutativeTransmission,
934 };
935 
938  Forward,
939  Reverse,
940 };
941 
944  Emission,
945  Transmission,
947 };
948 
967  const RadiationVector& J1,
968  const RadiationVector& J2,
969  const ArrayOfRadiationVector& dJ1,
970  const ArrayOfRadiationVector& dJ2,
971  const TransmissionMatrix& T,
972  const TransmissionMatrix& PiT,
973  const ArrayOfTransmissionMatrix& dT1,
974  const ArrayOfTransmissionMatrix& dT2,
975  const RadiativeTransferSolver solver);
976 
994  const PropagationMatrix& K,
995  const StokesVector& a,
996  const StokesVector& S,
997  const ArrayOfPropagationMatrix& dK,
998  const ArrayOfStokesVector& da,
999  const ArrayOfStokesVector& dS,
1000  const ConstVectorView B,
1001  const ConstVectorView dB_dT,
1003  const bool& jacobian_do);
1004 
1022  const PropagationMatrix& K1,
1023  const PropagationMatrix& K2,
1024  const ArrayOfPropagationMatrix& dK1,
1025  const ArrayOfPropagationMatrix& dK2,
1026  const Numeric& r,
1027  const Numeric& dr_dtemp1,
1028  const Numeric& dr_dtemp2,
1029  const Index temp_deriv_pos);
1030 
1038  const ArrayOfTransmissionMatrix& T,
1039  const CumulativeTransmission type) /*[[expects: T.nelem()>0]]*/;
1040 
1058  const RadiationVector & I_incoming,
1059  const ArrayOfTransmissionMatrix& T,
1060  const ArrayOfTransmissionMatrix& PiTf,
1061  const ArrayOfTransmissionMatrix& PiTr,
1062  const ArrayOfTransmissionMatrix& Z,
1066  const BackscatterSolver solver);
1067 
1078  ConstMatrixView m);
1079 
1090  ConstTensor5View t, const ArrayOfMatrix& aom);
1091 
1092 #endif // transmissionmatrix_h
LazyScale
Class to help with hidden temporary variables for operations of type Numeric times Class.
Definition: propagationmatrix.h:46
RadiationVector::setBackscatterTransmissionDerivative
void setBackscatterTransmissionDerivative(const RadiationVector &I0, const TransmissionMatrix &Tr, const TransmissionMatrix &Tf, const TransmissionMatrix &dZ)
Set this to backscatter transmission scatter derivative.
Definition: transmissionmatrix.h:744
Matrix
The Matrix class.
Definition: matpackI.h:1193
cumulative_backscatter_derivative
ArrayOfArrayOfTransmissionMatrix cumulative_backscatter_derivative(ConstTensor5View t, const ArrayOfMatrix &aom)
Accumulated backscatter derivative (???)
Definition: transmissionmatrix.cc:1755
transform
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...
Definition: matpackI.cc:1476
CumulativeTransmission
CumulativeTransmission
Intended to hold various ways to accumulate the transmission matrix.
Definition: transmissionmatrix.h:937
TransmissionMatrix::TransmissionMatrix
TransmissionMatrix(TransmissionMatrix &&tm) noexcept
Construct a new Transmission Matrix object.
Definition: transmissionmatrix.h:68
TransmissionMatrix::T3
std::vector< Eigen::Matrix3d, Eigen::aligned_allocator< Eigen::Matrix3d > > T3
Definition: transmissionmatrix.h:42
RadiationVector::add_avg
void add_avg(const RadiationVector &O1, const RadiationVector &O2)
Add the average of two other RadiationVector to *this.
Definition: transmissionmatrix.h:588
StokesVector
Stokes vector is as Propagation matrix but only has 4 possible values.
Definition: propagationmatrix.h:1075
RadiationVector::Vec2
const Eigen::Vector2d & Vec2(size_t i) const
Return Vector at position.
Definition: transmissionmatrix.h:510
RadiationVector::Vec4
const Eigen::Vector4d & Vec4(size_t i) const
Return Vector at position.
Definition: transmissionmatrix.h:496
RadiationVector::operator()
const Numeric & operator()(const Index i, const Index j) const
Access operator.
Definition: transmissionmatrix.h:792
PropagationMatrix::IsEmpty
bool IsEmpty() const
Asks if the class is empty.
Definition: propagationmatrix.h:221
TransmissionMatrix::operator<<
friend std::ostream & operator<<(std::ostream &os, const TransmissionMatrix &tm)
Output operator.
Definition: transmissionmatrix.cc:1804
ARTS::Var::jacobian_quantities
ArrayOfRetrievalQuantity jacobian_quantities(Workspace &ws) noexcept
Definition: autoarts.h:3924
w
Complex w(Complex z) noexcept
The Faddeeva function.
Definition: linefunctions.cc:42
TransmissionMatrix::TransmissionMatrix
TransmissionMatrix(Index nf=0, Index stokes=1)
Construct a new Transmission Matrix object.
Definition: transmissionmatrix.h:55
update_radiation_vector
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.
Definition: transmissionmatrix.cc:1484
Tensor3
The Tensor3 class.
Definition: matpackIII.h:339
RadiationVector::R3
std::vector< Eigen::Vector3d, Eigen::aligned_allocator< Eigen::Vector3d > > R3
Definition: transmissionmatrix.h:399
PropagationMatrix::NumberOfAzimuthAngles
Index NumberOfAzimuthAngles() const
The number of azimuth angles of the propagation matrix.
Definition: propagationmatrix.h:210
RadiationVector::operator=
RadiationVector & operator=(const RadiationVector &rv)=default
Assign old radiation vector to this.
stepwise_transmission
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.
Definition: transmissionmatrix.cc:1353
operator>>
std::istream & operator>>(std::istream &is, TransmissionMatrix &tm)
Input operator.
Definition: transmissionmatrix.cc:1844
ARTS::Group::Tensor3
Tensor3 Tensor3
Definition: autoarts.h:104
TransmissionMatrix::Weights
Struct of far and close weights.
Definition: transmissionmatrix.h:331
oem::Matrix
invlib::Matrix< ArtsMatrix > Matrix
invlib wrapper type for ARTS matrices.
Definition: oem.h:34
TransmissionMatrix::operator()
Numeric operator()(const Index i, const Index j, const Index k) const
Access value in matrix.
Definition: transmissionmatrix.h:251
TransmissionMatrix::linear_in_tau_weights
Weights< N > linear_in_tau_weights(size_t) const noexcept
Definition: transmissionmatrix.h:342
BackscatterSolver::CommutativeTransmission
@ CommutativeTransmission
RadiationVector::SetZero
void SetZero(size_t i)
Set Radiation Vector to Zero at position.
Definition: transmissionmatrix.h:474
RadiativeTransferSolver::Emission
@ Emission
RadiationVector::RadiationVector
RadiationVector(const RadiationVector &rv)=default
Construct a new Radiation Vector object.
TransmissionMatrix::Frequencies
Index Frequencies() const
Number of frequencies.
Definition: transmissionmatrix.h:268
RadiationVector::Vec2
Eigen::Vector2d & Vec2(size_t i)
Return Vector at position.
Definition: transmissionmatrix.h:557
RadiationVector::operator=
RadiationVector & operator=(const ConstMatrixView &M)
Set *this from matrix.
Definition: transmissionmatrix.h:763
ARTS::Var::jacobian_do
Index jacobian_do(Workspace &ws) noexcept
Definition: autoarts.h:3912
PropagationMatrix
Definition: propagationmatrix.h:87
operator<<
std::ostream & operator<<(std::ostream &os, const TransmissionMatrix &tm)
Output operator.
Definition: transmissionmatrix.cc:1804
data
G0 G2 FVC Y DV Numeric Numeric Numeric Zeeman LowerQuantumNumbers void * data
Definition: arts_api_classes.cc:232
TransmissionMatrix::Mat2
Eigen::Matrix2d & Mat2(size_t i)
Get Matrix at position.
Definition: transmissionmatrix.h:190
CumulativeTransmission::Forward
@ Forward
TransmissionMatrix::Mat4
Eigen::Matrix4d & Mat4(size_t i)
Get Matrix at position.
Definition: transmissionmatrix.h:176
TransmissionMatrix::Weights::close
Eigen::Matrix< Numeric, int(N), int(N)> close
Definition: transmissionmatrix.h:331
PropagationMatrix::K13
VectorView K13(const Index iz=0, const Index ia=0)
Vector view to K(0, 2) elements.
Definition: propagationmatrix.h:819
stepwise_source
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.
Definition: transmissionmatrix.cc:1371
RadiationVector::addDerivTransmission
void addDerivTransmission(const TransmissionMatrix &PiT, const TransmissionMatrix &dT, const RadiationVector &I)
Add the transmission derivative to this.
Definition: transmissionmatrix.h:659
RadiationVector::addDerivEmission
void addDerivEmission(const TransmissionMatrix &PiT, const TransmissionMatrix &dT, const TransmissionMatrix &T, const RadiationVector &ImJ, const RadiationVector &dJ)
Add the emission derivative to this.
Definition: transmissionmatrix.h:633
RadiationVector::setDerivReflection
void setDerivReflection(const RadiationVector &I, const TransmissionMatrix &PiT, const TransmissionMatrix &Z, const TransmissionMatrix &dZ)
Sets *this to the reflection derivative.
Definition: transmissionmatrix.h:701
RadiationVector::stokes_dim
Index stokes_dim
Definition: transmissionmatrix.h:397
Array
This can be used to make arrays out of anything.
Definition: array.h:108
RadiationVector::operator>>
friend std::istream & operator>>(std::istream &data, RadiationVector &rv)
Input operator.
Definition: transmissionmatrix.cc:1860
RadiationVector::setSource
void setSource(const StokesVector &a, const ConstVectorView &B, const StokesVector &S, Index i)
Set this to source vector at position.
Definition: transmissionmatrix.h:828
RadiativeTransferSolver
RadiativeTransferSolver
Intended to hold various forward solvers.
Definition: transmissionmatrix.h:943
RadiationVector::Vec3
Eigen::Vector3d & Vec3(size_t i)
Return Vector at position.
Definition: transmissionmatrix.h:550
BackscatterSolver
BackscatterSolver
Intended to hold various backscatter solvers.
Definition: transmissionmatrix.h:931
TransmissionMatrix::mul
void mul(const TransmissionMatrix &A, const TransmissionMatrix &B)
Set this to a multiple of A by B.
Definition: transmissionmatrix.h:222
PropagationMatrix::NumberOfZenithAngles
Index NumberOfZenithAngles() const
The number of zenith angles of the propagation matrix.
Definition: propagationmatrix.h:207
TransmissionMatrix::operator=
TransmissionMatrix & operator=(const LazyScale< TransmissionMatrix > &lstm)
Assign lazily.
Definition: transmissionmatrix.h:316
TransmissionMatrix::StokesDim
Index StokesDim() const
Stokes dimensionaility.
Definition: transmissionmatrix.h:265
ArrayOfArrayOfRadiationVector
Array< ArrayOfRadiationVector > ArrayOfArrayOfRadiationVector
Definition: transmissionmatrix.h:900
PropagationMatrix::K12
VectorView K12(const Index iz=0, const Index ia=0)
Vector view to K(0, 1) elements.
Definition: propagationmatrix.h:809
cumulative_transmission
ArrayOfTransmissionMatrix cumulative_transmission(const ArrayOfTransmissionMatrix &T, const CumulativeTransmission type)
Accumulate the transmission matrix over all layers.
Definition: transmissionmatrix.cc:1527
ArrayOfArrayOfTransmissionMatrix
Array< ArrayOfTransmissionMatrix > ArrayOfArrayOfTransmissionMatrix
Definition: transmissionmatrix.h:897
TransmissionMatrix::Mat3
Eigen::Matrix3d & Mat3(size_t i)
Get Matrix at position.
Definition: transmissionmatrix.h:183
ArrayOfArrayOfArrayOfTransmissionMatrix
Array< ArrayOfArrayOfTransmissionMatrix > ArrayOfArrayOfArrayOfTransmissionMatrix
Definition: transmissionmatrix.h:898
TransmissionMatrix::operator*=
TransmissionMatrix & operator*=(const Numeric &scale)
Scale self.
Definition: transmissionmatrix.h:303
set_backscatter_radiation_vector
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.
Definition: transmissionmatrix.cc:1548
RadiationVector::RadiationVector
RadiationVector(Index nf=0, Index stokes=1)
Construct a new Radiation Vector object.
Definition: transmissionmatrix.h:412
LazyScale::bas
const base & bas
A const reference to a value.
Definition: propagationmatrix.h:57
TransmissionMatrix::setZero
void setZero()
Set to zero matrix.
Definition: transmissionmatrix.h:208
ArrayOfRadiationVector
Array< RadiationVector > ArrayOfRadiationVector
Definition: transmissionmatrix.h:899
RadiationVector::Vec1
Eigen::Matrix< double, 1, 1 > & Vec1(size_t i)
Return Vector at position.
Definition: transmissionmatrix.h:564
RadiationVector::Vec1
const Eigen::Matrix< double, 1, 1 > & Vec1(size_t i) const
Return Vector at position.
Definition: transmissionmatrix.h:517
RadiationVector::R4
std::vector< Eigen::Vector4d, Eigen::aligned_allocator< Eigen::Vector4d > > R4
Definition: transmissionmatrix.h:398
jacobian.h
Routines for setting up the jacobian.
Numeric
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: matpack.h:33
TransmissionMatrix::Mat3
const Eigen::Matrix3d & Mat3(size_t i) const
Get Matrix at position.
Definition: transmissionmatrix.h:136
RadiationVector::Vec3
const Eigen::Vector3d & Vec3(size_t i) const
Return Vector at position.
Definition: transmissionmatrix.h:503
TransmissionMatrix::T1
std::vector< Eigen::Matrix< double, 1, 1 >, Eigen::aligned_allocator< Eigen::Matrix< double, 1, 1 > > > T1
Definition: transmissionmatrix.h:46
RadiationVector::add_weighted
void add_weighted(const TransmissionMatrix &T, const RadiationVector &close, const RadiationVector &far)
Add the weighted source of two RadiationVector to *this.
Definition: transmissionmatrix.h:606
RadiationVector::RadiationVector
RadiationVector(RadiationVector &&rv) noexcept
Construct a new Radiation Vector object.
Definition: transmissionmatrix.h:425
TransmissionMatrix::Mat1
Eigen::Matrix< double, 1, 1 > & Mat1(size_t i)
Get Matrix at position.
Definition: transmissionmatrix.h:197
RadiationVector::leftMul
void leftMul(const TransmissionMatrix &T)
Multiply radiation vector from the left.
Definition: transmissionmatrix.h:463
RadiationVector::Vec
Eigen::VectorXd Vec(size_t i) const
Return Vector at position by copy.
Definition: transmissionmatrix.h:524
propagationmatrix.h
Stuff related to the propagation matrix.
TransmissionMatrix::setIdentity
void setIdentity()
Set to identity matrix.
Definition: transmissionmatrix.h:200
PropagationMatrix::K14
VectorView K14(const Index iz=0, const Index ia=0)
Vector view to K(0, 3) elements.
Definition: propagationmatrix.h:829
RadiationVector::addMultiplied
void addMultiplied(const TransmissionMatrix &A, const RadiationVector &x)
Add multiply.
Definition: transmissionmatrix.h:682
ConstMatrixView
A constant view of a Matrix.
Definition: matpackI.h:982
cumulative_backscatter
ArrayOfTransmissionMatrix cumulative_backscatter(ConstTensor5View t, ConstMatrixView m)
Accumulated backscatter (???)
Definition: transmissionmatrix.cc:1714
oem::Identity
invlib::MatrixIdentity< Matrix > Identity
Definition: oem.h:39
ArrayOfTransmissionMatrix
Array< TransmissionMatrix > ArrayOfTransmissionMatrix
Definition: transmissionmatrix.h:896
N
#define N
Definition: rng.cc:164
RadiationVector
Radiation Vector for Stokes dimension 1-4.
Definition: transmissionmatrix.h:395
LazyScale::scale
const Numeric & scale
A const reference to a Numeric.
Definition: propagationmatrix.h:60
TransmissionMatrix::T4
std::vector< Eigen::Matrix4d, Eigen::aligned_allocator< Eigen::Matrix4d > > T4
Definition: transmissionmatrix.h:41
RadiationVector::StokesDim
Index StokesDim() const
Get Stokes dimension.
Definition: transmissionmatrix.h:873
RadiationVector::R2
std::vector< Eigen::Vector2d, Eigen::aligned_allocator< Eigen::Vector2d > > R2
Definition: transmissionmatrix.h:400
TransmissionMatrix::Mat1
const Eigen::Matrix< double, 1, 1 > & Mat1(size_t i) const
Get Matrix at position.
Definition: transmissionmatrix.h:150
PropagationMatrix::Kjj
VectorView Kjj(const Index iz=0, const Index ia=0)
Vector view to diagonal elements.
Definition: propagationmatrix.h:799
TransmissionMatrix::T2
std::vector< Eigen::Matrix2d, Eigen::aligned_allocator< Eigen::Matrix2d > > T2
Definition: transmissionmatrix.h:43
RadiationVector::setBackscatterTransmission
void setBackscatterTransmission(const RadiationVector &I0, const TransmissionMatrix &Tr, const TransmissionMatrix &Tf, const TransmissionMatrix &Z)
Set this to backscatter transmission.
Definition: transmissionmatrix.h:723
TransmissionMatrix::operator+=
TransmissionMatrix & operator+=(const LazyScale< TransmissionMatrix > &lstm)
Assign to *this lazily.
Definition: transmissionmatrix.h:286
TransmissionMatrix::TransmissionMatrix
TransmissionMatrix(const TransmissionMatrix &tm)=default
Construct a new Transmission Matrix object.
TransmissionMatrix::Weights::far
Eigen::Matrix< Numeric, int(N), int(N)> far
Definition: transmissionmatrix.h:331
ArrayOfArrayOfArrayOfRadiationVector
Array< ArrayOfArrayOfRadiationVector > ArrayOfArrayOfArrayOfRadiationVector
Definition: transmissionmatrix.h:901
RadiationVector::R1
std::vector< Eigen::Matrix< double, 1, 1 >, Eigen::aligned_allocator< Eigen::Matrix< double, 1, 1 > > > R1
Definition: transmissionmatrix.h:403
ARTS::Var::x
Vector x(Workspace &ws) noexcept
Definition: autoarts.h:7346
Index
INDEX Index
The type to use for all integer numbers and indices.
Definition: matpack.h:39
RadiationVector::rem_avg
void rem_avg(const RadiationVector &O1, const RadiationVector &O2)
Remove the average of two other RadiationVector from *this.
Definition: transmissionmatrix.h:571
TransmissionMatrix::Mat2
const Eigen::Matrix2d & Mat2(size_t i) const
Get Matrix at position.
Definition: transmissionmatrix.h:143
TransmissionMatrix::operator>>
friend std::istream & operator>>(std::istream &data, TransmissionMatrix &tm)
Input operator.
Definition: transmissionmatrix.cc:1844
RadiationVector::Frequencies
Index Frequencies() const
Get frequency count.
Definition: transmissionmatrix.h:876
TransmissionMatrix::Mat
Eigen::MatrixXd Mat(size_t i) const
Get Matrix at position by copy.
Definition: transmissionmatrix.h:157
M
#define M
Definition: rng.cc:165
TransmissionMatrix::operator=
TransmissionMatrix & operator=(TransmissionMatrix &&tm) noexcept
Move operator.
Definition: transmissionmatrix.h:93
RadiationVector::operator=
RadiationVector & operator=(RadiationVector &&rv) noexcept
Assign old radiation vector to this.
Definition: transmissionmatrix.h:450
TransmissionMatrix
Class to keep track of Transmission Matrices for Stokes Dim 1-4.
Definition: transmissionmatrix.h:38
TransmissionMatrix::stokes_dim
Index stokes_dim
Definition: transmissionmatrix.h:40
RadiationVector::Vec4
Eigen::Vector4d & Vec4(size_t i)
Return Vector at position.
Definition: transmissionmatrix.h:543
RadiationVector::operator<<
friend std::ostream & operator<<(std::ostream &os, const RadiationVector &rv)
Output operator.
Definition: transmissionmatrix.cc:1824
ConstVectorView
A constant view of a Vector.
Definition: matpackI.h:476
TransmissionMatrix::operator=
TransmissionMatrix & operator=(const TransmissionMatrix &tm)=default
Assignment operator.
TransmissionMatrix::mul_aliased
void mul_aliased(const TransmissionMatrix &A, const TransmissionMatrix &B)
Set this to a multiple of A by B.
Definition: transmissionmatrix.h:237
ConstTensor5View
A constant view of a Tensor5.
Definition: matpackV.h:143
operator*
LazyScale< TransmissionMatrix > operator*(const TransmissionMatrix &tm, const Numeric &x)
Lazy scale of Transmission Matrix.
Definition: transmissionmatrix.h:378
TransmissionMatrix::Mat4
const Eigen::Matrix4d & Mat4(size_t i) const
Get Matrix at position.
Definition: transmissionmatrix.h:129