HΦ  3.2.0
output.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 #include "output.h"
17 #include "FileIO.h"
18 #include "wrapperMPI.h"
19 
20 
26 
27 int output(struct BindStruct *X) {
28 
29  FILE *fp;
30  char sdt[D_FileNameMax];
31  long int i, i_max;
32  i_max = X->Check.idim_max;
33 
34  if (X->Def.iCalcType == FullDiag) {
35  switch (X->Def.iCalcModel) {
36  case Spin:
37  case Hubbard:
38  case Kondo:
39  sprintf(sdt, cFileNamePhys_FullDiag, X->Def.CDataFileHead, X->Def.Nup, X->Def.Ndown);
40  break;
41  case SpinGC:
42  case HubbardGC:
43  case KondoGC:
45  break;
46  default:
47  break;
48  }
49  if (childfopenMPI(sdt, "w", &fp) != 0) {
50  return -1;
51  }
52  fprintf(fp, " <H> <N> <Sz> <S2> <D> \n");
53  for (i = 0; i < i_max; i++) {
54  fprintf(fp, " %10lf %10lf %10lf %10lf %10lf\n", X->Phys.all_energy[i], X->Phys.all_num_up[i]+X->Phys.all_num_down[i], X->Phys.all_sz[i],
55  X->Phys.all_s2[i], X->Phys.all_doublon[i]);
56  }
57  fclose(fp);
58  }
59  else{
60  fprintf(stdoutMPI, "Error: output function is used only for FullDiag mode.");
61  return -1;
62  }
63 
64  return 0;
65 }
66 
72 
73 int outputHam(struct BindStruct *X){
74  long int i=0;
75  long int j=0;
76  long int imax = X->Check.idim_max;
77  long int ihermite=0;
78  char cHeader[256];
79  FILE *fp;
80  char sdt[D_FileNameMax];
81 
82 #pragma omp parallel for default(none) reduction(+:ihermite) firstprivate(imax) private(i, j) shared(Ham)
83  for (i=1; i<=imax; i++){
84  for (j=1; j<=i; j++){
85  if(cabs(Ham[i][j])>1.0e-13){
86  ihermite += 1;
87  }
88  }
89  }
90 
91  strcpy(cHeader, "%%%%MatrixMarket matrix coordinate complex hermitian\n");
93  if(childfopenMPI(sdt,"w",&fp)!=0){
94  return -1;
95  }
96  fprintf(fp, "%s", cHeader);
97  fprintf(fp, "%ld %ld %ld \n", imax, imax, ihermite);
98  for (i=1; i<=imax; i++){
99  for (j=1; j<=i; j++){
100  if(cabs(Ham[i][j])>1.0e-13){
101  fprintf(fp, "%ld %ld %lf %lf\n",i,j,creal(Ham[i][j]),cimag(Ham[i][j]));
102  }
103  }
104  }
105  fclose(fp);
106  return 0;
107 }
double * all_num_down
[CheckList::idim_max+1] Number of spin-down electrons for FullDiag and LOBPCG. malloc in setmem_large...
Definition: struct.h:381
double * all_s2
[CheckList::idim_max+1] for FullDiag and LOBPCG. malloc in setmem_large().
Definition: struct.h:377
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
const char * cFileNamePhys_FullDiag
Definition: global.c:74
const char * cFileNamePhys_FullDiag_GC
Definition: global.c:75
unsigned long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.h:303
int output(struct BindStruct *X)
output function for FullDiag mode
Definition: output.c:27
#define D_FileNameMax
Definition: global.h:23
int childfopenMPI(const char *_cPathChild, const char *_cmode, FILE **_fp)
Only the root process open file in output/ directory.
Definition: FileIO.c:27
unsigned int Nup
Number of spin-up electrons in this process.
Definition: struct.h:58
struct PhysList Phys
Physical quantities.
Definition: struct.h:413
unsigned int Ndown
Number of spin-down electrons in this process.
Definition: struct.h:59
Bind.
Definition: struct.h:409
double complex ** Ham
Definition: global.h:73
double * all_doublon
[CheckList::idim_max+1] Doublon for FullDiag and LOBPCG. malloc in setmem_large().
Definition: struct.h:373
int outputHam(struct BindStruct *X)
output Hamiltonian only used for FullDiag mode
Definition: output.c:73
const char * cFileNamePhys_FullDiag_Ham
Definition: global.c:78
struct EDMainCalStruct X
Definition: struct.h:432
double * all_sz
[CheckList::idim_max+1] for FullDiag and LOBPCG. malloc in setmem_large().
Definition: struct.h:375
int iCalcModel
Switch for model. 0:Hubbard, 1:Spin, 2:Kondo, 3:HubbardGC, 4:SpinGC, 5:KondoGC, 6:HubbardNConserved.
Definition: struct.h:198
struct CheckList Check
Size of the Hilbert space.
Definition: struct.h:411
char * CDataFileHead
Read from Calcmod in readdef.h. Header of output file such as Green&#39;s function.
Definition: struct.h:42
double * all_num_up
[CheckList::idim_max+1] Number of spin-up electrons for FullDiag and LOBPCG. malloc in setmem_large()...
Definition: struct.h:379
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:165
double * all_energy
[CheckList::idim_max+1] Energy for FullDiag and LOBPCG. malloc in setmem_large(). ...
Definition: struct.h:371