Accu
Functions
header.h File Reference

Go to the source code of this file.

Functions

void Send_ProcD (const int to, const int nin, const double xin[], const MPI::Intracomm &icomm)
void Recv_ProcD (const int from, int &nout, double xout[], const int maxbuf, const MPI::Intracomm &icomm)
void ExchangeD (const int yourid, const int nin, const double xin[], int &nout, double xout[], const int maxbuf, const MPI::Intracomm &icomm)
int ReadIn ()
void DebugD (const int n, const double x[])
void PivotD (const int n, double x[])
void vdcopy (const int n, double xp[], const int ix, const double yp[], const int iy)
void vdplus (const int n, double x[], const int ix, const double y[], const int iy, const double z[], const int iz)
void vddiv (const int n, double x[], const int ix, const double y[], const int iy, const double z[], const int iz)
void vdmult (const int n, double x[], const int ix, const double y[], const int iy, const double z[], const int iz)
double dscapr (const int n, const double x[], const double y[])
void vdaxpy (const int n, double x[], const double y[], const double &a, const double z[])
void vicopy (const int n, int xp[], const int ix, const int yp[], const int iy)
void IniGeom (const int myid, const int procx, const int procy, int neigh[], int &color)
void IniCoord (const int myid, const int procx, const int procy, double &xl, double &xr, double &yb, double &yt)
void AddBound (const int ib, const int nx, const int ny, double w[], const double s[])
void GetBound (const int ib, const int nx, const int ny, const double w[], double s[])
void VecAccu (const int nx, const int ny, double w[], const int neigh[], const int color, const MPI::Intracomm &icomm)
void FreeVecAccu ()

Function Documentation

void AddBound ( const int  ib,
const int  nx,
const int  ny,
double  w[],
const double  s[] 
)

Computes w := w + s at the interface/boundary nodes on the boundary ib . South (ib==1), East (ib==2), North (ib==3), West (ib==4)

Parameters:
[in]ibmy local boundary
[in]nxnumber of discretization intervals in x-direction
[in]nynumber of discretization intervals in y-direction
[in,out]wvector for all nodes of local discretization
[in]sshort vector with values on boundary ib
void DebugD ( const int  n,
const double  x[] 
)

Reads a integer number with the rank of a process from the terminal and writes the elements of vector x of this process onto the terminal. This will be repeated while rank number is valid, i.e., -1 stops the loop.

Parameters:
[in]nnumber of elements of vector x (of my process)
[in]xlocal vector (double) on my process
double dscapr ( const int  n,
const double  x[],
const double  y[] 
)
void ExchangeD ( const int  yourid,
const int  nin,
const double  xin[],
int &  nout,
double  xout[],
const int  maxbuf,
const MPI::Intracomm &  icomm 
)

Exchanges vectors between my process and the process with rank .

Parameters:
[in]youridrank of process in communicator icomm to send the message
[in]ninnumber of elements in vector xin
[in]xinmy vector
[out]noutnumber of received elements stored in vector xout
[out]xoutvector
[in]maxbufmaximal number of elements of xout
[in]icommcommunicator to use
void FreeVecAccu ( )

Frees the dynamic memory allocated in the first call in VecAccu.

See also:
VecAccu
void GetBound ( const int  ib,
const int  nx,
const int  ny,
const double  w[],
double  s[] 
)

Copies the values of w corresponding to boundary ib onto vector s. South (ib==1), East (ib==2), North (ib==3), West (ib==4). The vector s has to be long enough!!

Parameters:
[in]ibmy local boundary
[in]nxnumber of discretization intervals in x-direction
[in]nynumber of discretization intervals in y-direction
[in]wvector for all nodes of local discretization
[out]sshort vector with values on boundary ib
void IniCoord ( const int  myid,
const int  procx,
const int  procy,
double &  xl,
double &  xr,
double &  yb,
double &  yt 
)

The quadratic domain [0,1] x [0,1] is divided into procx * procy subdomains numbered rowise. According to my process with rank number myid, the coordinates of the lower left corner (xl, yb) and of the upper right corner (xr, yt) are generated.

Parameters:
[in]myidmy rank
[in]procxnumber of processes in x-direction
[in]procynumber of processes in y-direction
[out]xlx-coordinate of left boundary
[out]xrx-coordinate of right boundary
[out]yby-coordinate of lower boundary
[out]yty-coordinate of upper boundary
void IniGeom ( const int  myid,
const int  procx,
const int  procy,
int  neigh[],
int &  color 
)

The quadratic domain [0,1] x [0,1] is divided into procx * procy subdomains numbered rowise. According to my process with rank number myid, the vector neigh contains the ranks of the neighbours t the South, East, North, West.

Parameters:
[in]myidmy rank
[in]procxnumber of processes in x-direction
[in]procynumber of processes in y-direction
[out]neighvector of 4 elements containing the ranks of neighbouring processes, -1 indicates a boundary in this direction
[out]colorred/black coloring of subdomains (0/1)
void PivotD ( const int  n,
double  x[] 
)

Exchanges global minimum and maximum of the global vector x .

Parameters:
[in]nnumber of elements of vector x (of my process)
[in,out]xlocal vector (double) on my process
int ReadIn ( )

Reads one integer from the terminal and broadcasts it to all processes.

Returns:
a number read by root process 0 which is broadcasted to all processes
void Recv_ProcD ( const int  from,
int &  nout,
double  xout[],
const int  maxbuf,
const MPI::Intracomm &  icomm 
)

Receives data from process with rank from into a double-vector xout .

Parameters:
[in]fromrank of process in communicator icomm to receive the message from
[out]noutnumber of received elements stored in vector xout
[out]xoutvector
[in]maxbufmaximal number of elements of xout
[in]icommcommunicator to use
void Send_ProcD ( const int  to,
const int  nin,
const double  xin[],
const MPI::Intracomm &  icomm 
)

Sends data of double-vector xin to the process with rank to .

Parameters:
[in]torank of process in communicator icomm to send the message
[in]ninnumber of elements in vector xin
[in]xinvector
[in]icommcommunicator to use
void vdaxpy ( const int  n,
double  x[],
const double  y[],
const double &  a,
const double  z[] 
)
void vdcopy ( const int  n,
double  xp[],
const int  ix,
const double  yp[],
const int  iy 
)
void vddiv ( const int  n,
double  x[],
const int  ix,
const double  y[],
const int  iy,
const double  z[],
const int  iz 
)
void vdmult ( const int  n,
double  x[],
const int  ix,
const double  y[],
const int  iy,
const double  z[],
const int  iz 
)
void vdplus ( const int  n,
double  x[],
const int  ix,
const double  y[],
const int  iy,
const double  z[],
const int  iz 
)
void VecAccu ( const int  nx,
const int  ny,
double  w[],
const int  neigh[],
const int  color,
const MPI::Intracomm &  icomm 
)

Accumulates vector w, i.e., a distributed vector is converted into an accumulated vector. Dynamic memory allocation in the first call of this function.

Parameters:
[in]nxnumber of discretization intervals in x-direction
[in]nynumber of discretization intervals in y-direction
[in,out]wvector for all nodes of local discretization
[in]neighvector of 4 elements containing the ranks of neighbouring processes, -1 indicates a boundary in this direction
[in]colorred/black coloring of subdomains (0/1)
[in]icommcommunicator
void vicopy ( const int  n,
int  xp[],
const int  ix,
const int  yp[],
const int  iy 
)
 All Files Functions Variables