HΦ  3.1.0
HPhiMain.c
Go to the documentation of this file.
1 /* HPhi - Quantum Lattice Model Simulator */
2 /* Copyright (C) 2015 The University of Tokyo */
3 
4 /* This program is free software: you can redistribute it and/or modify */
5 /* it under the terms of the GNU General Public License as published by */
6 /* the Free Software Foundation, either version 3 of the License, or */
7 /* (at your option) any later version. */
8 
9 /* This program is distributed in the hope that it will be useful, */
10 /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
11 /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
12 /* GNU General Public License for more details. */
13 
14 /* You should have received a copy of the GNU General Public License */
15 /* along with this program. If not, see <http://www.gnu.org/licenses/>. */
16 /*-------------------------------------------------------------*/
17 
18 #include <sz.h>
19 #include <HPhiTrans.h>
20 #include <output_list.h>
21 #include <diagonalcalc.h>
22 #include <CalcByLanczos.h>
23 #include <CalcByLOBPCG.h>
24 #include <CalcByFullDiag.h>
25 #include <CalcByTPQ.h>
26 #include <CalcSpectrum.h>
27 #include <check.h>
28 #include "CalcByTEM.h"
29 #include "readdef.h"
30 #include "StdFace_main.h"
31 #include "wrapperMPI.h"
32 #include "splash.h"
33 #include "CalcTime.h"
34 
177 int main(int argc, char* argv[]){
178 
179  int mode=0;
180  char cFileListName[D_FileNameMax];
181 
182  stdoutMPI = stdout;
183  if(JudgeDefType(argc, argv, &mode)!=0){
184  FinalizeMPI();
185  return 0;
186  }
187 
188  if (mode == STANDARD_DRY_MODE) {
189  myrank = 0;
190  nproc = 1;
191  stdoutMPI = stdout;
192  splash();
193  }
194  else InitializeMPI(argc, argv);
195 
196  //Timer
197  InitTimer();
198  if (mode != STANDARD_DRY_MODE) StartTimer(0);
199 
200  //MakeDirectory for output
201  struct stat tmpst;
202  if (myrank == 0) {
203  if (stat(cParentOutputFolder, &tmpst) != 0) {
204  if (mkdir(cParentOutputFolder, 0777) != 0) {
205  fprintf(stdoutMPI, "%s", cErrOutput);
206  exitMPI(-1);
207  }
208  }
209  }/*if (myrank == 0)*/
210 
211  strcpy(cFileListName, argv[2]);
212 
213  if(mode==STANDARD_MODE || mode == STANDARD_DRY_MODE){
214  if (myrank == 0) StdFace_main(argv[2]);
215  strcpy(cFileListName, "namelist.def");
216  if (mode == STANDARD_DRY_MODE){
217  fprintf(stdout, "Dry run is Finished. \n\n");
218  return 0;
219  }
220  }
221 
222  setmem_HEAD(&X.Bind);
223  if(ReadDefFileNInt(cFileListName, &(X.Bind.Def), &(X.Bind.Boost))!=0){
224  fprintf(stdoutMPI, "%s", cErrDefFile);
225  FinalizeMPI();
226  return -1;
227  }
228 
229  if (X.Bind.Def.nvec < X.Bind.Def.k_exct){
230  fprintf(stdoutMPI, "%s", cErrnvec);
232  FinalizeMPI();
233  return -1;
234  }
235  fprintf(stdoutMPI, "%s", cProFinishDefFiles);
236 
237  /*ALLOCATE-------------------------------------------*/
238  setmem_def(&X.Bind, &X.Bind.Boost);
239  /*-----------------------------------------------------*/
240 
241  /*Read Def files.*/
243  if(ReadDefFileIdxPara(&(X.Bind.Def), &(X.Bind.Boost))!=0){
244  fprintf(stdoutMPI, "%s", cErrIndices);
245  FinalizeMPI();
246  return -1;
247  }
249  fprintf(stdoutMPI, "%s", cProFinishDefCheck);
250 
251  /*Set convergence Factor*/
253 
254  /*---------------------------*/
255  if(HPhiTrans(&(X.Bind))!=0) {
256  exitMPI(-1);
257  }
258 
259  //Start Calculation
260  if(X.Bind.Def.iFlgCalcSpec == CALCSPEC_NOT) {
261 
262  if(check(&(X.Bind))==MPIFALSE){
263  FinalizeMPI();
264  return -1;
265  }
266 
267  /*LARGE VECTORS ARE ALLOCATED*/
268  if (setmem_large(&X.Bind) != 0) {
270  exitMPI(-1);
271  }
272 
273  StartTimer(1000);
274  if(sz(&(X.Bind), list_1, list_2_1, list_2_2)!=0){
275  exitMPI(-1);
276  }
277 
278  StopTimer(1000);
279  if(X.Bind.Def.WRITE==1){
280  output_list(&(X.Bind));
281  FinalizeMPI();
282  return 0;
283  }
284  StartTimer(2000);
285  diagonalcalc(&(X.Bind));
286  StopTimer(2000);
287 
288  switch (X.Bind.Def.iCalcType) {
289  case Lanczos:
290  StartTimer(4000);
291  if (CalcByLanczos(&X) != TRUE) {
292  FinalizeMPI();
293  StopTimer(4000);
294  return 0;
295  }
296  StopTimer(4000);
297  break;
298 
299  case CG:
300  if (CalcByLOBPCG(&X) != TRUE) {
301  FinalizeMPI();
302  return 0;
303  }
304  break;
305 
306  case FullDiag:
307  StartTimer(5000);
308  if (X.Bind.Def.iFlgScaLAPACK ==0 && nproc != 1) {
309  fprintf(stdoutMPI, "Error: Full Diagonalization by LAPACK is only allowed for one process.\n");
310  FinalizeMPI();
311  }
312  if (CalcByFullDiag(&X) != TRUE) {
313  FinalizeMPI();
314  }
315  StopTimer(5000);
316  break;
317 
318  case TPQCalc:
319  StartTimer(3000);
321  FinalizeMPI();
322  StopTimer(3000);
323  return 0;
324  }
325  StopTimer(3000);
326  break;
327 
328  case TimeEvolution:
329  if(CalcByTEM(X.Bind.Def.Param.ExpecInterval, &X)!=0){
330  FinalizeMPI();
331  return 0;
332  }
333 
334  default:
335  FinalizeMPI();
336  StopTimer(0);
337  return 0;
338  }
339  }
340  else{
341  StartTimer(6000);
342  if (CalcSpectrum(&X) != TRUE) {
343  FinalizeMPI();
344  StopTimer(6000);
345  return 0;
346  }
347  StopTimer(6000);
348  }
349 
350  StopTimer(0);
351  OutputTimer(&(X.Bind));
352  FinalizeMPI();
353  return 0;
354 }
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:409
void StartTimer(int n)
function for initializing elapse time [start]
Definition: time.c:71
int JudgeDefType(const int argc, char *argv[], int *mode)
function of judging a type of define files.
Definition: readdef.c:2312
const char * cReadDefStart
Definition: LogMessage.c:29
void SetConvergenceFactor(struct DefineList *X)
function to set convergence factors
Definition: readdef.c:2442
int sz(struct BindStruct *X, long unsigned int *list_1_, long unsigned int *list_2_1_, long unsigned int *list_2_2_)
generating Hilbert space
Definition: sz.c:63
void StopTimer(int n)
function for calculating elapse time [elapse time=StartTimer-StopTimer]
Definition: time.c:83
char * cErrnvec
Definition: ErrorMessage.c:26
int CalcByTEM(const int ExpecInterval, struct EDMainCalStruct *X)
main function of time evolution calculation
Definition: CalcByTEM.c:53
int NumAve
Definition: global.h:65
#define D_FileNameMax
Definition: global.h:23
char * cErrDefFile
Definition: ErrorMessage.c:24
int ReadDefFileIdxPara(struct DefineList *X, struct BoostList *xBoost)
function of reading def files to get keyword index
Definition: readdef.c:976
void setmem_def(struct BindStruct *X, struct BoostList *xBoost)
Set size of memories for Def and Phys in BindStruct.
Definition: xsetmem.c:56
int CalcByLanczos(struct EDMainCalStruct *X)
A main function to calculate eigenvalues and eigenvectors by Lanczos method.
Definition: CalcByLanczos.c:54
int iFlgScaLAPACK
ScaLAPACK mode ( only for FullDiag )
Definition: struct.h:234
unsigned int nvec
Read from Calcmod in readdef.h.
Definition: struct.h:46
int ReadDefFileNInt(char *xNameListFile, struct DefineList *X, struct BoostList *xBoost)
Function of reading information about "ModPara" file and total number of parameters from other def fi...
Definition: readdef.c:456
int check(struct BindStruct *X)
A program to check size of dimension for Hilbert-space.
Definition: check.c:51
int CalcByLOBPCG(struct EDMainCalStruct *X)
Driver routine for LOB(P)CG method.
Definition: CalcByLOBPCG.c:634
int setmem_large(struct BindStruct *X)
Set size of memories for Hamiltonian (Ham, L_vec), vectors(vg, v0, v1, v2, vec, alpha, beta), lists (list_1, list_2_1, list_2_2, list_Diagonal) and Phys(BindStruct.PhysList) struct in the case of Full Diag mode.
Definition: xsetmem.c:224
#define TRUE
Definition: global.h:26
const char * cReadDefFinish
Definition: LogMessage.c:30
void InitTimer()
function for initializing Timer[]
Definition: time.c:53
int iErrCodeMem
Error Message in HPhiMain.c.
Definition: ErrorMessage.c:20
int ExpecInterval
Definition: struct.h:36
void splash()
Print logo mark and version number.
Definition: splash.c:25
int nproc
Number of processors, defined in InitializeMPI()
Definition: global.h:161
int CalcByFullDiag(struct EDMainCalStruct *X)
Parent function for FullDiag mode.
int diagonalcalc(struct BindStruct *X)
Calculate diagonal components and obtain the list, list_diagonal.
Definition: diagonalcalc.c:85
int CalcByTPQ(const int NumAve, const int ExpecInterval, struct EDMainCalStruct *X)
A main function to calculate physical quqntities by TPQ method.
Definition: CalcByTPQ.c:52
struct BoostList Boost
For Boost.
Definition: struct.h:413
#define MPIFALSE
Definition: global.h:24
char * cErrnvecShow
Definition: ErrorMessage.c:27
long unsigned int * list_2_1
Definition: global.h:49
void setmem_HEAD(struct BindStruct *X)
Set size of memories headers of output files.
Definition: xsetmem.c:42
void OutputTimer(struct BindStruct *X)
function for outputting elapse time for each function
Definition: time.c:95
void FinalizeMPI()
MPI Finitialization wrapper.
Definition: wrapperMPI.c:74
int CalcSpectrum(struct EDMainCalStruct *X)
A main function to calculate spectrum.
Definition: CalcSpectrum.c:90
char * cErrLargeMem
Definition: ErrorMessage.c:35
char * cErrOutput
Definition: ErrorMessage.c:23
struct ParamList Param
Definition: struct.h:240
long unsigned int * list_1
Definition: global.h:47
int main(int argc, char *argv[])
Main program for HPhi.
Definition: HPhiMain.c:177
const char * cProFinishDefCheck
const char * cFileNameTimeKeep
Definition: global.c:23
long unsigned int * list_2_2
Definition: global.h:50
struct EDMainCalStruct X
Definition: struct.h:431
const char * cProFinishDefFiles
int WRITE
It is ALWAYS 0 ???
Definition: struct.h:54
const char * cParentOutputFolder
Definition: global.c:20
char * cErrIndices
Definition: ErrorMessage.c:28
int myrank
Process ID, defined in InitializeMPI()
Definition: global.h:162
int iFlgCalcSpec
Input parameter CalcSpec in teh CalcMod file.
Definition: struct.h:215
int HPhiTrans(struct BindStruct *X)
Function of checking transfers not to count the same type of operators. .
Definition: HPhiTrans.c:45
int output_list(struct BindStruct *X)
Output list_1 for canonical ensembles.
Definition: output_list.c:40
void StdFace_main(char *fname)
Main routine for the standard mode.
void InitializeMPI(int argc, char *argv[])
MPI initialization wrapper Process ID (myrank), Number of processes (nproc), Number of threads (nthre...
Definition: wrapperMPI.c:43
struct BindStruct Bind
Binded struct.
Definition: struct.h:419
int TimeKeeper(struct BindStruct *X, const char *cFileName, const char *cTimeKeeper_Message, const char *cWriteType)
Functions for writing a time log.
Definition: log.c:42
int iCalcType
Switch for calculation type. 0:Lanczos, 1:TPQCalc, 2:FullDiag.
Definition: struct.h:192
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:164
unsigned int k_exct
Read from Calcmod in readdef.h.
Definition: struct.h:47