NESSi  v1.0.2
The NonEquilibrium Systems Simulation Library

◆ Bcast_timestep()

template<typename T >
void cntr::herm_matrix< T >::Bcast_timestep ( int  tstp,
int  root 
)

Broadcasts the herm_matrix at a given time step to all tasks.

Purpose

Broadcasts the herm_matrix at a given time step tstp to all tasks.

Parameters
tstp

Time step which should be broadcasted.

root

The task rank from which the herm_matrix should be broadcasted.

Definition at line 3473 of file cntr_herm_matrix_impl.hpp.

3473  {
3474  int numtasks, taskid;
3475  herm_matrix_timestep<T> Gtemp;
3476  MPI_Comm_size(MPI_COMM_WORLD, &numtasks);
3477  MPI_Comm_rank(MPI_COMM_WORLD, &taskid);
3478  Gtemp.resize(tstp, ntau_, size1_);
3479  if (taskid == root)
3480  this->get_timestep(tstp, Gtemp);
3481  if (sizeof(T) == sizeof(double))
3482  MPI_Bcast(Gtemp.data_, Gtemp.total_size_,
3483  MPI_DOUBLE_COMPLEX, root, MPI_COMM_WORLD);
3484  else
3485  MPI_Bcast(Gtemp.data_, Gtemp.total_size_, MPI_COMPLEX,
3486  root, MPI_COMM_WORLD);
3487  if (taskid != root)
3488  this->set_timestep(tstp, Gtemp);
3489 }
void get_timestep(int tstp, herm_matrix_timestep< T > &timestep) const
void set_timestep(int tstp, herm_matrix &g1)
Sets all components at time step tstp to the components of a given herm_matrix. ...