HPhi++  3.1.0
wrapperMPI.cpp File Reference

MPI wrapper for init, finalize, bcast, etc. More...

#include <mpi.h>
#include <cstdio>
#include <cstdlib>
#include <cstring>
#include "wrapperMPI.hpp"
#include <cmath>
#include <complex>
#include "splash.hpp"
#include "global.hpp"
#include "common/setmemory.hpp"

Go to the source code of this file.

Functions

void InitializeMPI (int argc, char *argv[])
 MPI initialization wrapper Process ID (myrank), Number of processes (nproc), Number of threads (nthreads), and pointer to the standard output (stdoutMPI) are specified here. More...
 
void FinalizeMPI ()
 MPI Finitialization wrapper. More...
 
void exitMPI (int errorcode)
 MPI Abortation wrapper. More...
 
FILE * fopenMPI (const char *FileName, const char *mode)
 MPI file I/O (open) wrapper. Only the root node (myrank = 0) should be open/read/write (small) parameter files. More...
 
char * fgetsMPI (char *InputString, int maxcount, FILE *fp)
 MPI file I/O (get a line, fgets) wrapper. Only the root node (myrank = 0) reads and broadcast string. More...
 
void BarrierMPI ()
 MPI barrier wrapper. More...
 
long int MaxMPI_li (long int idim)
 MPI wrapper function to obtain maximum unsigned long integer across processes. More...
 
double MaxMPI_d (double dvalue)
 MPI wrapper function to obtain maximum Double across processes. More...
 
std::complex< double > SumMPI_dc (std::complex< double > norm)
 MPI wrapper function to obtain sum of Double complex across processes. More...
 
double SumMPI_d (double norm)
 MPI wrapper function to obtain sum of Double across processes. More...
 
void SumMPI_dv (int nnorm, double *norm)
 MPI wrapper function to obtain sum of Double array across processes. More...
 
void SumMPI_cv (int nnorm, std::complex< double > *norm)
 MPI wrapper function to obtain sum of Double array across processes. More...
 
long int SumMPI_li (long int idim)
 MPI wrapper function to obtain sum of unsigned long integer across processes. More...
 
int SumMPI_i (int idim)
 MPI wrapper function to obtain sum of integer across processes. More...
 
long int BcastMPI_li (int root, long int idim)
 MPI wrapper function to broadcast long integer across processes. More...
 
double NormMPI_dc (long int idim, std::complex< double > *_v1)
 Compute norm of process-distributed vector \(|{\bf v}_1|^2\). More...
 
void NormMPI_dv (long int ndim, int nstate, std::complex< double > **_v1, double *dnorm)
 Compute norm of process-distributed vector \(|{\bf v}_1|^2\). More...
 
std::complex< double > VecProdMPI (long int ndim, std::complex< double > *v1, std::complex< double > *v2)
 Compute conjugate scaler product of process-distributed vector \({\bf v}_1^* \cdot {\bf v}_2\). More...
 
void MultiVecProdMPI (long int ndim, int nstate, std::complex< double > **v1, std::complex< double > **v2, std::complex< double > *prod)
 Compute conjugate scaler product of process-distributed vector \({\bf v}_1^* \cdot {\bf v}_2\). More...
 
void SendRecv_cv (int origin, long int nMsgS, long int nMsgR, std::complex< double > *vecs, std::complex< double > *vecr)
 Wrapper of MPI_Sendrecv for std::complex<double> number. When we pass a message longer than 2^31-1 (max of int: 2147483647), we need to divide it. More...
 
void SendRecv_iv (int origin, long int nMsgS, long int nMsgR, long int *vecs, long int *vecr)
 Wrapper of MPI_Sendrecv for long integer number. When we pass a message longer than 2^31-1 (max of int: 2147483647), we need to divide it. More...
 
long int SendRecv_i (int origin, long int isend)
 Wrapper of MPI_Sendrecv for long integer number. More...
 

Detailed Description

MPI wrapper for init, finalize, bcast, etc.

Definition in file wrapperMPI.cpp.

Function Documentation

◆ BarrierMPI()

void BarrierMPI ( )

MPI barrier wrapper.

Author
Mitsuaki Kawamura (The University of Tokyo)

Definition at line 160 of file wrapperMPI.cpp.

160  {
161 #ifdef __MPI
162  MPI_Barrier(MPI_COMM_WORLD);
163 #endif
164 }/*void BarrierMPI()*/

◆ BcastMPI_li()

long int BcastMPI_li ( int  root,
long int  idim 
)

MPI wrapper function to broadcast long integer across processes.

Returns
Broadcasted value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]rootThe source process of the broadcast
[in]idimValue to be broadcasted

Definition at line 305 of file wrapperMPI.cpp.

Referenced by Initialize_wave().

308  {
309  long int idim0;
310  idim0 = idim;
311 #ifdef __MPI
312  MPI_Bcast(&idim0, 1, MPI_LONG, root, MPI_COMM_WORLD);
313 #endif
314  return(idim0);
315 }/*long int BcastMPI_li*/

◆ exitMPI()

void exitMPI ( int  errorcode)

MPI Abortation wrapper.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]errorcodeError-code to be returned as that of this program

Definition at line 86 of file wrapperMPI.cpp.

Referenced by CalcByLOBPCG(), CalcByTEM(), CalcByTPQ(), CalcSpectrum(), CalcSpectrumByBiCG(), CheckMPI(), CheckMPI_Summary(), dsfmt_chk_init_by_array(), dsfmt_chk_init_gen_rand(), Initialize_wave(), InitializeMPI(), main(), MakeExcitedList(), MaxMPI_d(), MaxMPI_li(), Output_restart(), phys(), Read_sz(), ReadDefFileIdxPara(), SeedSwitch(), SendRecv_cv(), SendRecv_i(), SendRecv_iv(), SumMPI_cv(), SumMPI_d(), SumMPI_dc(), SumMPI_dv(), SumMPI_i(), SumMPI_li(), sz(), and X_GC_child_CisAjt_Hubbard_MPI().

89 {
90  fflush(stdout);
91 #ifdef __MPI
92  fprintf(stdout,"\n\n ####### [HPhi] You DO NOT have to WORRY about the following MPI-ERROR MESSAGE. #######\n\n");
93  int ierr;
94  ierr = MPI_Abort(MPI_COMM_WORLD, errorcode);
95  ierr = MPI_Finalize();
96  if (ierr != 0) fprintf(stderr, "\n MPI_Finalize() = %d\n\n", ierr);
97 #endif
98  exit(errorcode);
99 }/*void exitMPI*/

◆ fgetsMPI()

char* fgetsMPI ( char *  InputString,
int  maxcount,
FILE *  fp 
)

MPI file I/O (get a line, fgets) wrapper. Only the root node (myrank = 0) reads and broadcast string.

Returns
The same as that of fgets
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[out]InputStringread line.
[in]maxcountLength of string
[in]fpfile pointer

Definition at line 122 of file wrapperMPI.cpp.

References myrank.

Referenced by CalcSpectrumByBiCG(), GetFileName(), inputHam(), Read_sz(), ReadcalcmodFile(), ReadDefFileIdxPara(), ReadDefFileNInt(), and SetOmega().

126  {
127  int inull;
128  char *ctmp;
129 
130  ctmp = InputString;
131  inull = 0;
132  if (myrank == 0) {
133  ctmp = fgets(InputString, maxcount, fp);
134  if (ctmp == NULL){
135  inull = 1;
136  }
137 
138  while(*InputString == '\n' || strncmp(InputString, "#", 1)==0){
139  ctmp = fgets(InputString, maxcount, fp);
140  if (ctmp == NULL){
141  inull=1;
142  break;
143  }
144  }
145  }
146 #ifdef __MPI
147  MPI_Bcast(InputString, maxcount, MPI_CHAR, 0, MPI_COMM_WORLD);
148  MPI_Bcast(&inull, 1, MPI_INT, 0, MPI_COMM_WORLD);
149 #endif
150  if (myrank != 0 && inull == 1) {
151  ctmp = NULL;
152  }
153 
154  return ctmp;
155 }/*char* fgetsMPI*/
int myrank
Process ID, defined in InitializeMPI()
Definition: global.cpp:73

◆ FinalizeMPI()

void FinalizeMPI ( )

MPI Finitialization wrapper.

Author
Mitsuaki Kawamura (The University of Tokyo)

Definition at line 74 of file wrapperMPI.cpp.

References myrank, and stdoutMPI.

Referenced by main(), and MakeExcitedList().

74  {
75 #ifdef __MPI
76  int ierr;
77  ierr = MPI_Finalize();
78  if (ierr != 0) fprintf(stderr, "\n MPI_Finalize() = %d\n\n", ierr);
79 #endif
80  if (myrank != 0) fclose(stdoutMPI);
81 }
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.cpp:75
int myrank
Process ID, defined in InitializeMPI()
Definition: global.cpp:73

◆ fopenMPI()

FILE* fopenMPI ( const char *  FileName,
const char *  mode 
)

MPI file I/O (open) wrapper. Only the root node (myrank = 0) should be open/read/write (small) parameter files.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]FileNameInput/output file
[in]mode"w", "r", etc.

Definition at line 105 of file wrapperMPI.cpp.

References myrank.

Referenced by childfopenMPI(), GetFileName(), ReadcalcmodFile(), ReadDefFileIdxPara(), and ReadDefFileNInt().

108  {
109  FILE* fp;
110 
111  if (myrank == 0) fp = fopen(FileName, mode);
112  else fp = fopen("/dev/null", "w");
113 
114  return fp;
115 }/*FILE* fopenMPI*/
int myrank
Process ID, defined in InitializeMPI()
Definition: global.cpp:73

◆ InitializeMPI()

void InitializeMPI ( int  argc,
char *  argv[] 
)

MPI initialization wrapper Process ID (myrank), Number of processes (nproc), Number of threads (nthreads), and pointer to the standard output (stdoutMPI) are specified here.

Author
Mitsuaki Kawamura (The University of Tokyo)

Definition at line 44 of file wrapperMPI.cpp.

References exitMPI(), myrank, nproc, nthreads, splash(), and stdoutMPI.

Referenced by main().

44  {
45 #ifdef __MPI
46  int ierr;
47  ierr = MPI_Init(&argc, &argv);
48  ierr = MPI_Comm_size(MPI_COMM_WORLD, &nproc);
49  ierr = MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
50  if(ierr != 0) exitMPI(ierr);
51 #else
52  nproc = 1;
53  myrank = 0;
54 #endif
55  if (myrank == 0) stdoutMPI = stdout;
56  else stdoutMPI = fopen("/dev/null", "w");
57  splash();
58 
59 #pragma omp parallel default(none) shared(nthreads)
60 #pragma omp master
61 #ifdef _OPENMP
62  nthreads = omp_get_num_threads();
63 #else
64  nthreads=1;
65 #endif
66  fprintf(stdoutMPI, "\n\n##### Parallelization Info. #####\n\n");
67  fprintf(stdoutMPI, " OpenMP threads : %d\n", nthreads);
68  fprintf(stdoutMPI, " MPI PEs : %d \n\n", nproc);
69 }/*void InitializeMPI(int argc, char *argv[])*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86
int nproc
Number of processors, defined in InitializeMPI()
Definition: global.cpp:72
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.cpp:75
void splash()
Print logo mark and version number.
Definition: splash.cpp:25
int nthreads
Number of Threads, defined in InitializeMPI()
Definition: global.cpp:74
int myrank
Process ID, defined in InitializeMPI()
Definition: global.cpp:73

◆ MaxMPI_d()

double MaxMPI_d ( double  dvalue)

MPI wrapper function to obtain maximum Double across processes.

Returns
Maximum value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]dvalueValue to be maximized

Definition at line 188 of file wrapperMPI.cpp.

References exitMPI().

Referenced by check().

190  {
191 #ifdef __MPI
192  int ierr;
193  ierr = MPI_Allreduce(MPI_IN_PLACE, &dvalue, 1,
194  MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
195  if(ierr != 0) exitMPI(-1);
196 #endif
197  return(dvalue);
198 }/*double MaxMPI_d*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ MaxMPI_li()

long int MaxMPI_li ( long int  idim)

MPI wrapper function to obtain maximum unsigned long integer across processes.

Returns
Maximum value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]idimValue to be maximized

Definition at line 171 of file wrapperMPI.cpp.

References exitMPI().

Referenced by check(), MakeExcitedList(), and setmem_large().

173  {
174 #ifdef __MPI
175  int ierr;
176  ierr = MPI_Allreduce(MPI_IN_PLACE, &idim, 1,
177  MPI_LONG, MPI_MAX, MPI_COMM_WORLD);
178  if(ierr != 0) exitMPI(-1);
179 #endif
180  return(idim);
181 }/*long int MaxMPI_li*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ MultiVecProdMPI()

void MultiVecProdMPI ( long int  ndim,
int  nstate,
std::complex< double > **  v1,
std::complex< double > **  v2,
std::complex< double > *  prod 
)

Compute conjugate scaler product of process-distributed vector \({\bf v}_1^* \cdot {\bf v}_2\).

Parameters
[in]ndimLocal dimension of vector
[in]v1[ndim] vector to be producted
[in]v2[ndim] vector to be producted

Definition at line 401 of file wrapperMPI.cpp.

References SumMPI_cv().

Referenced by expec_cisajs_Hubbard(), expec_cisajs_HubbardGC(), expec_cisajs_SpinGCGeneral(), expec_cisajs_SpinGCHalf(), expec_cisajs_SpinGeneral(), expec_cisajs_SpinHalf(), expec_cisajscktalt_Hubbard(), expec_cisajscktalt_HubbardGC(), expec_cisajscktalt_SpinGCGeneral(), expec_cisajscktalt_SpinGCHalf(), expec_cisajscktalt_SpinGeneral(), and expec_cisajscktalt_SpinHalf().

407  {
408  long int idim;
409  int istate;
410 
411  for (istate = 0; istate < nstate; istate++) prod[istate] = 0.0;
412  for (idim = 1; idim <= ndim; idim++) {
413  for (istate = 0; istate < nstate; istate++) {
414  prod[istate] += conj(v1[idim][istate])*v2[idim][istate];
415  }
416  }
417  SumMPI_cv(nstate, prod);
418 }/*void MultiVecProdMPI*/
std::complex< double > ** v1
Definition: global.cpp:21
void SumMPI_cv(int nnorm, std::complex< double > *norm)
MPI wrapper function to obtain sum of Double array across processes.
Definition: wrapperMPI.cpp:254

◆ NormMPI_dc()

double NormMPI_dc ( long int  idim,
std::complex< double > *  _v1 
)

Compute norm of process-distributed vector \(|{\bf v}_1|^2\).

Returns
Norm \(|{\bf v}_1|^2\)
Parameters
[in]idimLocal dimension of vector
[in]_v1[idim] vector to be producted

Definition at line 321 of file wrapperMPI.cpp.

References SumMPI_d().

Referenced by CalcSpectrumByBiCG().

324  {
325  double dnorm =0;
326  long int i;
327 
328  dnorm = 0.0;
329 #pragma omp parallel for default(none) private(i) \
330 shared(_v1, idim) reduction(+:dnorm)
331  for (i = 1; i <= idim; i++)
332  dnorm += real(conj(_v1[i])*_v1[i]);
333 
334 #ifdef __MPI
335  dnorm = SumMPI_d(dnorm);
336 #endif
337  return dnorm;
338 }/*double NormMPI_dc*/
double SumMPI_d(double norm)
MPI wrapper function to obtain sum of Double across processes.
Definition: wrapperMPI.cpp:222

◆ NormMPI_dv()

void NormMPI_dv ( long int  ndim,
int  nstate,
std::complex< double > **  _v1,
double *  dnorm 
)

Compute norm of process-distributed vector \(|{\bf v}_1|^2\).

Returns
Norm \(|{\bf v}_1|^2\)
Parameters
[in]ndimLocal dimension of vector
[in]_v1[idim] vector to be producted

Definition at line 344 of file wrapperMPI.cpp.

References SumMPI_dv().

Referenced by Initialize_wave(), LOBPCG_Main(), and Multiply().

349  {
350  long int idim;
351  int istate;
352 
353  for (istate = 0; istate < nstate; istate++) dnorm[istate] = 0.0;
354  for (idim = 1; idim <= ndim; idim++) {
355  for (istate = 0; istate < nstate; istate++) {
356  dnorm[istate] += real(conj(_v1[idim][istate])*_v1[idim][istate]);
357  }
358  }
359  SumMPI_dv(nstate, dnorm);
360  for (istate = 0; istate < nstate; istate++) dnorm[istate] = sqrt(dnorm[istate]);
361 }/*double NormMPI_cv*/
void SumMPI_dv(int nnorm, double *norm)
MPI wrapper function to obtain sum of Double array across processes.
Definition: wrapperMPI.cpp:238

◆ SendRecv_cv()

void SendRecv_cv ( int  origin,
long int  nMsgS,
long int  nMsgR,
std::complex< double > *  vecs,
std::complex< double > *  vecr 
)

Wrapper of MPI_Sendrecv for std::complex<double> number. When we pass a message longer than 2^31-1 (max of int: 2147483647), we need to divide it.

Definition at line 424 of file wrapperMPI.cpp.

References exitMPI().

Referenced by X_Ajt_MPI(), X_child_CisAit_GeneralSpin_MPIdouble(), X_child_CisAit_spin_MPIdouble(), X_child_CisAitCjuAjv_GeneralSpin_MPIdouble(), X_child_CisAitCjuAjv_GeneralSpin_MPIsingle(), X_child_CisAjt_MPIdouble(), X_child_CisAjt_MPIsingle(), X_child_CisAjtCkuAku_Hubbard_MPI(), X_child_CisAjtCkuAlv_Hubbard_MPI(), X_child_general_hopp_MPIdouble(), X_child_general_hopp_MPIsingle(), X_child_general_int_spin_MPIdouble(), X_child_general_int_spin_MPIsingle(), X_child_general_int_spin_TotalS_MPIdouble(), X_Cis_MPI(), X_GC_Ajt_MPI(), X_GC_child_CisAisCjuAjv_GeneralSpin_MPIdouble(), X_GC_child_CisAisCjuAjv_GeneralSpin_MPIsingle(), X_GC_child_CisAisCjuAjv_spin_MPIdouble(), X_GC_child_CisAisCjuAjv_spin_MPIsingle(), X_GC_child_CisAit_GeneralSpin_MPIdouble(), X_GC_child_CisAit_spin_MPIdouble(), X_GC_child_CisAitCiuAiv_spin_MPIdouble(), X_GC_child_CisAitCiuAiv_spin_MPIsingle(), X_GC_child_CisAitCjuAju_GeneralSpin_MPIdouble(), X_GC_child_CisAitCjuAju_spin_MPIdouble(), X_GC_child_CisAitCjuAjv_GeneralSpin_MPIdouble(), X_GC_child_CisAitCjuAjv_GeneralSpin_MPIsingle(), X_GC_child_CisAjtCkuAku_Hubbard_MPI(), X_GC_child_CisAjtCkuAlv_Hubbard_MPI(), X_GC_child_general_hopp_MPIdouble(), X_GC_child_general_hopp_MPIsingle(), and X_GC_Cis_MPI().

430  {
431 #ifdef __MPI
432  int ierr, two31m1 = 2147483647, modMsg, nMsgS2, nMsgR2;
433  long int nMsg, nnMsg, iMsg, sMsgR, sMsgS;
434  MPI_Status statusMPI;
435 
436  if (nMsgS > nMsgR) nMsg = nMsgS;
437  else nMsg = nMsgR;
438  nnMsg = nMsg / two31m1;
439  modMsg = nMsg % two31m1;
440  if (modMsg != 0) nnMsg += 1;
441 
442  sMsgS = 0;
443  sMsgR = 0;
444  for (iMsg = 0; iMsg < nnMsg; iMsg++) {
445  nMsgS2 = nMsgS / nnMsg;
446  nMsgR2 = nMsgR / nnMsg;
447  if (iMsg < nMsgS % nnMsg) nMsgS2 += 1;
448  if (iMsg < nMsgR % nnMsg) nMsgR2 += 1;
449 
450  ierr = MPI_Sendrecv(&vecs[sMsgS], nMsgS2, MPI_DOUBLE_COMPLEX, origin, 0,
451  &vecr[sMsgR], nMsgR2, MPI_DOUBLE_COMPLEX, origin, 0,
452  MPI_COMM_WORLD, &statusMPI);
453  if (ierr != 0) exitMPI(-1);
454 
455  sMsgS += nMsgS2;
456  sMsgR += nMsgR2;
457  }
458 #endif
459 }/*void SendRecv_cv*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SendRecv_i()

long int SendRecv_i ( int  origin,
long int  isend 
)

Wrapper of MPI_Sendrecv for long integer number.

Definition at line 504 of file wrapperMPI.cpp.

References exitMPI().

Referenced by X_Ajt_MPI(), X_child_CisAit_GeneralSpin_MPIdouble(), X_child_CisAit_spin_MPIdouble(), X_child_CisAitCjuAjv_GeneralSpin_MPIdouble(), X_child_CisAitCjuAjv_GeneralSpin_MPIsingle(), X_child_CisAjt_MPIdouble(), X_child_CisAjt_MPIsingle(), X_child_CisAjtCkuAku_Hubbard_MPI(), X_child_CisAjtCkuAlv_Hubbard_MPI(), X_child_general_hopp_MPIdouble(), X_child_general_hopp_MPIsingle(), X_child_general_int_spin_MPIdouble(), X_child_general_int_spin_MPIsingle(), X_child_general_int_spin_TotalS_MPIdouble(), X_Cis_MPI(), X_GC_Ajt_MPI(), X_GC_child_CisAisCjuAjv_spin_MPIdouble(), X_GC_child_CisAisCjuAjv_spin_MPIsingle(), X_GC_child_CisAit_spin_MPIdouble(), X_GC_child_CisAitCiuAiv_spin_MPIdouble(), X_GC_child_CisAitCiuAiv_spin_MPIsingle(), X_GC_child_CisAitCjuAju_spin_MPIdouble(), X_GC_child_CisAjtCkuAku_Hubbard_MPI(), X_GC_child_CisAjtCkuAlv_Hubbard_MPI(), X_GC_child_general_hopp_MPIdouble(), X_GC_child_general_hopp_MPIsingle(), and X_GC_Cis_MPI().

507  {
508 #ifdef __MPI
509  int ierr;
510  MPI_Status statusMPI;
511  long int ircv;
512  ierr = MPI_Sendrecv(&isend, 1, MPI_LONG, origin, 0,
513  &ircv, 1, MPI_LONG, origin, 0,
514  MPI_COMM_WORLD, &statusMPI);
515  if (ierr != 0) exitMPI(ierr);
516  return ircv;
517 #else
518  return isend;
519 #endif
520 }/*void SendRecv_i*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SendRecv_iv()

void SendRecv_iv ( int  origin,
long int  nMsgS,
long int  nMsgR,
long int *  vecs,
long int *  vecr 
)

Wrapper of MPI_Sendrecv for long integer number. When we pass a message longer than 2^31-1 (max of int: 2147483647), we need to divide it.

Definition at line 465 of file wrapperMPI.cpp.

References exitMPI().

Referenced by X_Ajt_MPI(), X_child_CisAit_GeneralSpin_MPIdouble(), X_child_CisAit_spin_MPIdouble(), X_child_CisAitCjuAjv_GeneralSpin_MPIdouble(), X_child_CisAitCjuAjv_GeneralSpin_MPIsingle(), X_child_CisAjt_MPIdouble(), X_child_CisAjt_MPIsingle(), X_child_CisAjtCkuAku_Hubbard_MPI(), X_child_CisAjtCkuAlv_Hubbard_MPI(), X_child_general_hopp_MPIdouble(), X_child_general_hopp_MPIsingle(), X_child_general_int_spin_MPIdouble(), X_child_general_int_spin_MPIsingle(), X_child_general_int_spin_TotalS_MPIdouble(), and X_Cis_MPI().

471  {
472 #ifdef __MPI
473  int ierr, two31m1 = 2147483647, modMsg, nMsgS2, nMsgR2;
474  long int nMsg, nnMsg, iMsg, sMsgR, sMsgS;
475  MPI_Status statusMPI;
476 
477  if (nMsgS > nMsgR) nMsg = nMsgS;
478  else nMsg = nMsgR;
479  nnMsg = nMsg / two31m1;
480  modMsg = nMsg % two31m1;
481  if (modMsg != 0) nnMsg += 1;
482 
483  sMsgS = 0;
484  sMsgR = 0;
485  for (iMsg = 0; iMsg < nnMsg; iMsg++) {
486  nMsgS2 = nMsgS / nnMsg;
487  nMsgR2 = nMsgR / nnMsg;
488  if (iMsg < nMsgS % nnMsg) nMsgS2 += 1;
489  if (iMsg < nMsgR % nnMsg) nMsgR2 += 1;
490 
491  ierr = MPI_Sendrecv(&vecs[sMsgS], nMsgS2, MPI_LONG, origin, 0,
492  &vecr[sMsgR], nMsgR2, MPI_LONG, origin, 0,
493  MPI_COMM_WORLD, &statusMPI);
494  if (ierr != 0) exitMPI(-1);
495 
496  sMsgS += nMsgS2;
497  sMsgR += nMsgR2;
498  }
499 #endif
500 }/*void SendRecv_iv*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SumMPI_cv()

void SumMPI_cv ( int  nnorm,
std::complex< double > *  norm 
)

MPI wrapper function to obtain sum of Double array across processes.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]normValue to be summed

Definition at line 254 of file wrapperMPI.cpp.

References exitMPI().

Referenced by LOBPCG_Main(), and MultiVecProdMPI().

257  {
258 #ifdef __MPI
259  int ierr;
260  ierr = MPI_Allreduce(MPI_IN_PLACE, norm, nnorm,
261  MPI_DOUBLE_COMPLEX, MPI_SUM, MPI_COMM_WORLD);
262  if (ierr != 0) exitMPI(-1);
263 #endif
264 }/*void SumMPI_cv*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SumMPI_d()

double SumMPI_d ( double  norm)

MPI wrapper function to obtain sum of Double across processes.

Returns
Sumed value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]normValue to be summed

Definition at line 222 of file wrapperMPI.cpp.

References exitMPI().

Referenced by MultiplyForTEM(), and NormMPI_dc().

224  {
225 #ifdef __MPI
226  int ierr;
227  ierr = MPI_Allreduce(MPI_IN_PLACE, &norm, 1,
228  MPI_DOUBLE_PRECISION, MPI_SUM, MPI_COMM_WORLD);
229  if(ierr != 0) exitMPI(-1);
230 #endif
231  return(norm);
232 }/*double SumMPI_d*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SumMPI_dc()

std::complex<double> SumMPI_dc ( std::complex< double >  norm)

MPI wrapper function to obtain sum of Double complex across processes.

Returns
Sumed value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]normValue to be summed

Definition at line 205 of file wrapperMPI.cpp.

References exitMPI().

Referenced by VecProdMPI().

207  {
208 #ifdef __MPI
209  int ierr;
210  ierr = MPI_Allreduce(MPI_IN_PLACE, &norm, 1,
211  MPI_DOUBLE_COMPLEX, MPI_SUM, MPI_COMM_WORLD);
212  if(ierr != 0) exitMPI(-1);
213 #endif
214  return(norm);
215 }/*std::complex<double> SumMPI_dc*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SumMPI_dv()

void SumMPI_dv ( int  nnorm,
double *  norm 
)

MPI wrapper function to obtain sum of Double array across processes.

Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]normValue to be summed

Definition at line 238 of file wrapperMPI.cpp.

References exitMPI().

Referenced by expec_energy_flct(), expec_energy_flct_GeneralSpin(), expec_energy_flct_GeneralSpinGC(), expec_energy_flct_HalfSpin(), expec_energy_flct_HalfSpinGC(), expec_energy_flct_Hubbard(), expec_energy_flct_HubbardGC(), LOBPCG_Main(), NormMPI_dv(), totalspin_Hubbard(), totalspin_HubbardGC(), totalspin_Spin(), and totalspin_SpinGC().

241  {
242 #ifdef __MPI
243  int ierr;
244  ierr = MPI_Allreduce(MPI_IN_PLACE, norm, nnorm,
245  MPI_DOUBLE_PRECISION, MPI_SUM, MPI_COMM_WORLD);
246  if (ierr != 0) exitMPI(-1);
247 #endif
248 }/*void SumMPI_dv*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SumMPI_i()

int SumMPI_i ( int  idim)

MPI wrapper function to obtain sum of integer across processes.

Returns
Sumed value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]idimValue to be summed

Definition at line 288 of file wrapperMPI.cpp.

References exitMPI().

Referenced by CheckMPI_Summary().

290  {
291 #ifdef __MPI
292  int ierr;
293  ierr = MPI_Allreduce(MPI_IN_PLACE, &idim, 1,
294  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
295  if(ierr != 0) exitMPI(-1);
296 #endif
297  return(idim);
298 }/*int SumMPI_i*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ SumMPI_li()

long int SumMPI_li ( long int  idim)

MPI wrapper function to obtain sum of unsigned long integer across processes.

Returns
Sumed value across processes.
Author
Mitsuaki Kawamura (The University of Tokyo)
Parameters
[in]idimValue to be summed

Definition at line 271 of file wrapperMPI.cpp.

References exitMPI().

Referenced by CheckMPI_Summary(), and Initialize_wave().

273  {
274 #ifdef __MPI
275  int ierr;
276  ierr = MPI_Allreduce(MPI_IN_PLACE, &idim, 1,
277  MPI_LONG, MPI_SUM, MPI_COMM_WORLD);
278  if(ierr != 0) exitMPI(-1);
279 #endif
280  return(idim);
281 }/*long int SumMPI_li*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.cpp:86

◆ VecProdMPI()

std::complex<double> VecProdMPI ( long int  ndim,
std::complex< double > *  v1,
std::complex< double > *  v2 
)

Compute conjugate scaler product of process-distributed vector \({\bf v}_1^* \cdot {\bf v}_2\).

Returns
Conjugate scaler product \({\bf v}_1^* \cdot {\bf v}_2\)
Parameters
[in]ndimLocal dimension of vector
[in]v1[ndim] vector to be producted
[in]v2[ndim] vector to be producted

Definition at line 367 of file wrapperMPI.cpp.

References nthreads, and SumMPI_dc().

Referenced by CalcSpectrumByBiCG().

371  {
372  long int idim;
373  std::complex<double> prod, *prod_thr;
374  int mythread;
375 
376  prod_thr = cd_1d_allocate(nthreads);
377 #pragma omp parallel default(none) shared(v1,v2,ndim,prod,prod_thr) private(idim,mythread)
378  {
379 #ifdef _OPENMP
380  mythread = omp_get_thread_num();
381 #else
382  mythread = 0;
383 #endif
384 #pragma omp for
385  for (idim = 1; idim <= ndim; idim++)
386  prod_thr[mythread] += conj(v1[idim]) * v2[idim];
387  }
388  prod = 0.0;
389  for (mythread = 0; mythread < nthreads; mythread++)
390  prod += prod_thr[mythread];
391  free_cd_1d_allocate(prod_thr);
392 
393  prod = SumMPI_dc(prod);
394 
395  return(prod);
396 }/*std::complex<double> VecProdMPI*/
std::complex< double > ** v1
Definition: global.cpp:21
int nthreads
Number of Threads, defined in InitializeMPI()
Definition: global.cpp:74
std::complex< double > SumMPI_dc(std::complex< double > norm)
MPI wrapper function to obtain sum of Double complex across processes.
Definition: wrapperMPI.cpp:205