23 inline void get_mat(
const int m, std::complex<T> &G_mat, herm_matrix_timestep_view<T> &G,
24 herm_matrix_timestep_view<T> &Gcc){
25 assert(m <= G.ntau());
26 assert(G.ntau() == Gcc.ntau());
27 assert(G.size1() == Gcc.size1());
28 assert(G.size2() == Gcc.size2());
29 assert(G.sig() == Gcc.sig());
56 void get_mat(
const int m, cdmatrix &G_mat, herm_matrix_timestep_view<T> &G,
57 herm_matrix_timestep_view<T> &Gcc){
58 assert(m <= G.ntau());
59 assert(G.ntau() == Gcc.ntau());
60 assert(G.size1() == Gcc.size1());
61 assert(G.size2() == Gcc.size2());
62 assert(G.sig() == Gcc.sig());
63 int size1=G.size1(), size2=G.size2();
67 map_ptr2matrix(size1, size2, mat, G_mat);
92 inline void get_mat(
const int m, std::complex<T> &G_mat, herm_matrix_timestep_view<T> &G){
93 assert(m <= G.ntau());
117 template <
typename T>
118 void get_mat(
const int m, cdmatrix &G_mat, herm_matrix_timestep_view<T> &G){
119 assert(m <= G.ntau());
120 std::complex<T> *mat;
121 int size1=G.size1(), size2=G.size2();
124 map_ptr2matrix(size1, size2, mat, G_mat);
153 template <
typename T>
154 inline void get_les(
const int i,
const int j, std::complex<T> &G_les, herm_matrix_timestep_view<T> &G,
155 herm_matrix_timestep_view<T> &Gcc){
156 assert(i <= G.tstp() && j <= G.tstp());
157 assert(i == G.tstp() || j == G.tstp());
158 assert(i <= Gcc.tstp() && j <= Gcc.tstp());
159 assert(i == Gcc.tstp() || j == Gcc.tstp());
160 assert(G.tstp() == Gcc.tstp());
161 assert(G.ntau() == Gcc.ntau());
162 assert(G.size1() == Gcc.size1());
163 assert(G.size2() == Gcc.size2());
164 assert(G.sig() == Gcc.sig());
167 G_les = *G.lesptr(i);
168 }
else if (G.tstp() == i) {
169 G_les = *Gcc.lesptr(j);
170 G_les = -std::conj(G_les);
198 template <
typename T>
199 void get_les(
const int i,
const int j, cdmatrix &G_les, herm_matrix_timestep_view<T> &G,
200 herm_matrix_timestep_view<T> &Gcc){
201 assert(i <= G.tstp() && j <= G.tstp());
202 assert(i == G.tstp() || j == G.tstp());
203 assert(i <= Gcc.tstp() && j <= Gcc.tstp());
204 assert(i == Gcc.tstp() || j == Gcc.tstp());
205 assert(G.tstp() == Gcc.tstp());
206 assert(G.ntau() == Gcc.ntau());
207 assert(G.size1() == Gcc.size1());
208 assert(G.size2() == Gcc.size2());
209 assert(G.sig() == Gcc.sig());
210 std::complex<T> *les;
211 int size1=G.size1(), size2=G.size2();
215 map_ptr2matrix(size1, size2, les, G_les);
216 }
else if (G.tstp() == i) {
218 map_ptr2matrix(size1, size2, les, G_les);
219 G_les.adjointInPlace();
251 template <
typename T>
252 inline void get_ret(
const int i,
const int j, std::complex<T> &G_ret, herm_matrix_timestep_view<T> &G,
253 herm_matrix_timestep_view<T> &Gcc){
254 assert(i <= G.tstp() && j <= G.tstp());
255 assert(i == G.tstp() || j == G.tstp());
256 assert(i <= Gcc.tstp() && j <= Gcc.tstp());
257 assert(i == Gcc.tstp() || j == Gcc.tstp());
258 assert(G.tstp() == Gcc.tstp());
259 assert(G.ntau() == Gcc.ntau());
260 assert(G.size1() == Gcc.size1());
261 assert(G.size2() == Gcc.size2());
262 assert(G.sig() == Gcc.sig());
265 G_ret = *G.retptr(j);
266 }
else if (G.tstp() == j) {
267 G_ret = *Gcc.retptr(i);
268 G_ret = -std::conj(G_ret);
297 template <
typename T>
298 void get_ret(
const int i,
const int j, cdmatrix &G_ret, herm_matrix_timestep_view<T> &G,
299 herm_matrix_timestep_view<T> &Gcc){
300 assert(i <= G.tstp() && j <= G.tstp());
301 assert(i == G.tstp() || j == G.tstp());
302 assert(i <= Gcc.tstp() && j <= Gcc.tstp());
303 assert(i == Gcc.tstp() || j == Gcc.tstp());
304 assert(G.tstp() == Gcc.tstp());
305 assert(G.ntau() == Gcc.ntau());
306 assert(G.size1() == Gcc.size1());
307 assert(G.size2() == Gcc.size2());
308 assert(G.sig() == Gcc.sig());
309 std::complex<T> *ret;
310 int size1=G.size1(), size2=G.size2();
314 map_ptr2matrix(size1, size2, ret, G_ret);
315 }
else if (G.tstp() == j) {
317 map_ptr2matrix(size1, size2, ret, G_ret);
318 G_ret.adjointInPlace();
350 template <
typename T>
351 inline void get_tv(
const int i,
const int m, std::complex<T> &G_tv, herm_matrix_timestep_view<T> &G,
352 herm_matrix_timestep_view<T> &Gcc){
353 assert(i == G.tstp());
354 assert(m <= G.ntau());
355 assert(G.tstp() == Gcc.tstp());
356 assert(G.ntau() == Gcc.ntau());
357 assert(G.size1() == Gcc.size1());
358 assert(G.size2() == Gcc.size2());
359 assert(G.sig() == Gcc.sig());
388 template <
typename T>
389 void get_tv(
const int i,
const int m, cdmatrix &G_tv, herm_matrix_timestep_view<T> &G, herm_matrix_timestep_view<T> &Gcc){
390 assert(i == G.tstp());
391 assert(m <= G.ntau());
392 assert(G.tstp() == Gcc.tstp());
393 assert(G.ntau() == Gcc.ntau());
394 assert(G.size1() == Gcc.size1());
395 assert(G.size2() == Gcc.size2());
396 assert(G.sig() == Gcc.sig());
397 int size1 = G.size1(), size2 = G.size2();
401 map_ptr2matrix(size1, size2, tv, G_tv);
434 template <
typename T>
435 inline void get_vt(
const int m,
const int i, std::complex<T> &G_vt, herm_matrix_timestep_view<T> &G,
436 herm_matrix_timestep_view<T> &Gcc){
437 assert(i == G.tstp());
438 assert(m <= G.ntau());
439 assert(G.tstp() == Gcc.tstp());
440 assert(G.ntau() == Gcc.ntau());
441 assert(G.size1() == Gcc.size1());
442 assert(G.size2() == Gcc.size2());
443 assert(G.sig() == Gcc.sig());
445 G_vt = *Gcc.tvptr(G.ntau() - m);
446 G_vt = std::complex<T>(-G.sig(),0.0) * std::conj(G_vt);
473 template <
typename T>
474 void get_vt(
const int m,
const int i, cdmatrix &G_vt, herm_matrix_timestep_view<T> &G,
475 herm_matrix_timestep_view<T> &Gcc){
476 assert(i == G.tstp());
477 assert(m <= G.ntau());
478 assert(G.tstp() == Gcc.tstp());
479 assert(G.ntau() == Gcc.ntau());
480 assert(G.size1() == Gcc.size1());
481 assert(G.size2() == Gcc.size2());
482 assert(G.sig() == Gcc.sig());
483 int size1 = G.size1(), size2 = G.size2();
486 vt = Gcc.tvptr(G.ntau() - m);
487 map_ptr2matrix(size1, size2, vt, G_vt);
488 G_vt.adjointInPlace();
489 G_vt = std::complex<T>(-G.sig(),0.0) * G_vt;
519 template <
typename T>
520 inline void get_gtr(
const int i,
const int j, std::complex<T> &G_gtr, herm_matrix_timestep_view<T> &G,
521 herm_matrix_timestep_view<T> &Gcc){
522 std::complex<T> G_les, G_ret;
526 G_gtr = G_les + G_ret;
553 template <
typename T>
554 void get_gtr(
const int i,
const int j, cdmatrix &G_gtr, herm_matrix_timestep_view<T> &G,
555 herm_matrix_timestep_view<T> &Gcc){
556 cdmatrix G_les, G_ret;
560 G_gtr.noalias() = G_les + G_ret;
591 template <
typename T>
592 inline void get_les(
const int i,
const int j, std::complex<T> &G_les, herm_matrix_timestep_view<T> &G){
618 template <
typename T>
619 void get_les(
const int i,
const int j, cdmatrix &G_les, herm_matrix_timestep_view<T> &G){
648 template <
typename T>
649 inline void get_ret(
const int i,
const int j, std::complex<T> &G_ret, herm_matrix_timestep_view<T> &G){
676 template <
typename T>
677 void get_ret(
const int i,
const int j, cdmatrix &G_ret, herm_matrix_timestep_view<T> &G){
706 template <
typename T>
707 inline void get_tv(
const int i,
const int m, std::complex<T> &G_tv, herm_matrix_timestep_view<T> &G){
733 template <
typename T>
734 void get_tv(
const int i,
const int m, cdmatrix &G_tv, herm_matrix_timestep_view<T> &G){
766 template <
typename T>
767 inline void get_vt(
const int m,
const int i, std::complex<T> &G_vt, herm_matrix_timestep_view<T> &G){
793 template <
typename T>
794 void get_vt(
const int m,
const int i, cdmatrix &G_vt, herm_matrix_timestep_view<T> &G){
823 template <
typename T>
824 inline void get_gtr(
const int i,
const int j, std::complex<T> &G_gtr, herm_matrix_timestep_view<T> &G){
850 template <
typename T>
851 void get_gtr(
const int i,
const int j, cdmatrix &G_gtr, herm_matrix_timestep_view<T> &G){
void get_tv(const int i, const int m, std::complex< T > &G_tv, herm_matrix_timestep_view< T > &G, herm_matrix_timestep_view< T > &Gcc)
Returns the left-mixing component of a general contour function at given times. ...
void get_vt(const int m, const int i, std::complex< T > &G_vt, herm_matrix_timestep_view< T > &G, herm_matrix_timestep_view< T > &Gcc)
Returns the right-mixing component of a general contour function at given times. ...
void get_mat(const int m, std::complex< T > &G_mat, herm_matrix_timestep_view< T > &G, herm_matrix_timestep_view< T > &Gcc)
Returns the Matsubara component of a general contour function at given times.
void get_gtr(const int i, const int j, std::complex< T > &G_gtr, herm_matrix_timestep_view< T > &G, herm_matrix_timestep_view< T > &Gcc)
Returns the greater component of a general contour function at given times.
void get_les(const int i, const int j, std::complex< T > &G_les, herm_matrix_timestep_view< T > &G, herm_matrix_timestep_view< T > &Gcc)
Returns the lesser component of a general contour function at given times.
void get_ret(const int i, const int j, std::complex< T > &G_ret, herm_matrix_timestep_view< T > &G, herm_matrix_timestep_view< T > &Gcc)
Returns the retarded component of a general contour function at given times.