HΦ  3.2.0
wrapperMPI.c File Reference

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

#include <mpi.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "wrapperMPI.h"
#include <math.h>
#include <complex.h>
#include "splash.h"
#include "global.h"
+ Include dependency graph for wrapperMPI.c:

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...
 
unsigned long int MaxMPI_li (unsigned 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...
 
double complex SumMPI_dc (double complex 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...
 
unsigned long int SumMPI_li (unsigned 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...
 
unsigned long int BcastMPI_li (int root, unsigned long int idim)
 MPI wrapper function to broadcast unsigned long integer across processes. More...
 
double NormMPI_dc (unsigned long int idim, double complex *_v1)
 Compute norm of process-distributed vector \(|{\bf v}_1|^2\). More...
 
double complex VecProdMPI (long unsigned int ndim, double complex *v1, double complex *v2)
 Compute conjugate scaler product of process-distributed vector \({\bf v}_1^* \cdot {\bf v}_2\). More...
 

Detailed Description

MPI wrapper for init, finalize, bcast, etc.

Definition in file wrapperMPI.c.

Function Documentation

◆ BarrierMPI()

void BarrierMPI ( )

MPI barrier wrapper.

Author
Mitsuaki Kawamura (The University of Tokyo)

Definition at line 160 of file wrapperMPI.c.

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

◆ BcastMPI_li()

unsigned long int BcastMPI_li ( int  root,
unsigned long int  idim 
)

MPI wrapper function to broadcast unsigned 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 273 of file wrapperMPI.c.

Referenced by Initialize_wave(), Lanczos_EigenVector(), and SetInitialVector().

276  {
277  unsigned long int idim0;
278  idim0 = idim;
279 #ifdef MPI
280  MPI_Bcast(&idim0, 1, MPI_UNSIGNED_LONG, root, MPI_COMM_WORLD);
281 #endif
282  return(idim0);
283 }/*unsigned long int BcastMPI_li*/
+ Here is the caller graph for this function:

◆ 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.c.

Referenced by CalcByLanczos(), CalcByLOBPCG(), CalcByTEM(), CalcByTPQ(), CalcSpectrum(), CalcSpectrumByBiCG(), CalcSpectrumByTPQ(), 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(), SumMPI_d(), SumMPI_dc(), SumMPI_i(), SumMPI_li(), sz(), 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_CisAjt_Hubbard_MPI(), 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().

89 {
90  int ierr;
91  fflush(stdout);
92 #ifdef MPI
93  fprintf(stdout,"\n\n ####### [HPhi] You DO NOT have to WORRY about the following MPI-ERROR MESSAGE. #######\n\n");
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.c.

References myrank.

Referenced by GetFileName(), inputHam(), Read_sz(), ReadcalcmodFile(), ReadDefFileIdxPara(), ReadDefFileNInt(), ReadTMComponents(), ReadTMComponents_BiCG(), 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.h:163
+ Here is the caller graph for this function:

◆ FinalizeMPI()

void FinalizeMPI ( )

MPI Finitialization wrapper.

Author
Mitsuaki Kawamura (The University of Tokyo)

Definition at line 74 of file wrapperMPI.c.

References myrank, and stdoutMPI.

Referenced by main(), and MakeExcitedList().

74  {
75  int ierr;
76 #ifdef MPI
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 }
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:163
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:165
+ Here is the caller graph for this function:

◆ 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.c.

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.h:163
+ Here is the caller graph for this function:

◆ 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 43 of file wrapperMPI.c.

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

Referenced by main().

43  {
44  int ierr;
45 
46 #ifdef MPI
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.c:86
int nthreads
Number of Threads, defined in InitializeMPI()
Definition: global.h:164
void splash()
Print logo mark and version number.
Definition: splash.c:25
int nproc
Number of processors, defined in InitializeMPI()
Definition: global.h:162
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:163
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:165
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ 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.c.

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.c:86
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ MaxMPI_li()

unsigned long int MaxMPI_li ( unsigned 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.c.

References exitMPI().

Referenced by check(), GetPairExcitedStateGeneralSpin(), GetPairExcitedStateHalfSpin(), GetPairExcitedStateHubbard(), GetSingleExcitedStateHubbard(), GetSingleExcitedStateHubbardGC(), and setmem_large().

173  {
174 #ifdef MPI
175  int ierr;
176  ierr = MPI_Allreduce(MPI_IN_PLACE, &idim, 1,
177  MPI_UNSIGNED_LONG, MPI_MAX, MPI_COMM_WORLD);
178  if(ierr != 0) exitMPI(-1);
179 #endif
180  return(idim);
181 }/*unsigned long int MaxMPI_li*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ NormMPI_dc()

double NormMPI_dc ( unsigned long int  idim,
double complex *  _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 289 of file wrapperMPI.c.

References SumMPI_dc().

Referenced by CalcSpectrum().

292  {
293  double complex cdnorm=0;
294  double dnorm =0;
295  unsigned long int i;
296  //DEBUG
297 #pragma omp parallel for default(none) private(i) firstprivate(myrank) shared(_v1, idim) reduction(+: cdnorm)
298  for(i=1;i<=idim;i++){
299  cdnorm += conj(_v1[i])*_v1[i];
300  }
301 #ifdef MPI
302  cdnorm = SumMPI_dc(cdnorm);
303 #endif
304  dnorm=creal(cdnorm);
305  dnorm=sqrt(dnorm);
306 
307  return dnorm;
308 }/*double NormMPI_dc*/
double complex SumMPI_dc(double complex norm)
MPI wrapper function to obtain sum of Double complex across processes.
Definition: wrapperMPI.c:205
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ 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.c.

References exitMPI().

Referenced by CG_EigenVector(), expec_energy_flct_GeneralSpin(), expec_energy_flct_GeneralSpinGC(), expec_energy_flct_HalfSpin(), expec_energy_flct_HalfSpinGC(), expec_energy_flct_Hubbard(), expec_energy_flct_HubbardGC(), Lanczos_EigenVector(), and PowerLanczos().

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.c:86
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SumMPI_dc()

double complex SumMPI_dc ( double complex  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.c.

References exitMPI().

Referenced by CG_EigenVector(), 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(), expec_cisajscktalt_SpinHalf(), expec_energy_flct(), Lanczos_EigenValue(), Lanczos_EigenVector(), Lanczos_GetTridiagonalMatrixComponents(), mltply(), Multiply(), MultiplyForTEM(), NormMPI_dc(), PowerLanczos(), SetDiagonalTEChemi(), SetDiagonalTEInterAll(), SetDiagonalTETransfer(), SetInitialVector(), totalspin_Spin(), totalspin_SpinGC(), totalSz_HubbardGC(), totalSz_SpinGC(), and 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 }/*double complex SumMPI_dc*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ 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 256 of file wrapperMPI.c.

References exitMPI().

Referenced by CheckMPI_Summary().

258  {
259 #ifdef MPI
260  int ierr;
261  ierr = MPI_Allreduce(MPI_IN_PLACE, &idim, 1,
262  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
263  if(ierr != 0) exitMPI(-1);
264 #endif
265  return(idim);
266 }/*int SumMPI_i*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SumMPI_li()

unsigned long int SumMPI_li ( unsigned 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 239 of file wrapperMPI.c.

References exitMPI().

Referenced by CheckMPI_Summary(), Initialize_wave(), Lanczos_EigenValue(), Lanczos_EigenVector(), Lanczos_GetTridiagonalMatrixComponents(), and SetInitialVector().

241  {
242 #ifdef MPI
243  int ierr;
244  ierr = MPI_Allreduce(MPI_IN_PLACE, &idim, 1,
245  MPI_UNSIGNED_LONG, MPI_SUM, MPI_COMM_WORLD);
246  if(ierr != 0) exitMPI(-1);
247 #endif
248  return(idim);
249 }/*unsigned long int SumMPI_li*/
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ VecProdMPI()

double complex VecProdMPI ( long unsigned int  ndim,
double complex *  v1,
double complex *  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 314 of file wrapperMPI.c.

References SumMPI_dc().

Referenced by CalcSpectrumByBiCG(), Initialize_wave(), InitShadowRes(), and LOBPCG_Main().

318  {
319  long unsigned int idim;
320  double complex prod;
321 
322  prod = 0.0;
323 #pragma omp parallel for default(none) shared(v1,v2,ndim) private(idim) reduction(+: prod)
324  for (idim = 1; idim <= ndim; idim++) prod += conj(v1[idim]) * v2[idim];
325  prod = SumMPI_dc(prod);
326 
327  return(prod);
328 }/*double complex VecProdMPI*/
double complex SumMPI_dc(double complex norm)
MPI wrapper function to obtain sum of Double complex across processes.
Definition: wrapperMPI.c:205
double complex * v1
Definition: global.h:35
double complex * v2
Definition: global.h:36
+ Here is the call graph for this function:
+ Here is the caller graph for this function: