HΦ  3.2.0
output.c File Reference
#include "output.h"
#include "FileIO.h"
#include "wrapperMPI.h"
+ Include dependency graph for output.c:

Go to the source code of this file.

Functions

int output (struct BindStruct *X)
 output function for FullDiag mode More...
 
int outputHam (struct BindStruct *X)
 output Hamiltonian only used for FullDiag mode More...
 

Function Documentation

◆ output()

int output ( struct BindStruct X)

output function for FullDiag mode

Parameters
[in]XStruct to get information about file header names, dimension of hirbert space, calc type, physical quantities.
Return values
0normally finished.
-1abnormally finished.

Definition at line 27 of file output.c.

References PhysList::all_doublon, PhysList::all_energy, PhysList::all_num_down, PhysList::all_num_up, PhysList::all_s2, PhysList::all_sz, DefineList::CDataFileHead, cFileNamePhys_FullDiag, cFileNamePhys_FullDiag_GC, BindStruct::Check, childfopenMPI(), D_FileNameMax, BindStruct::Def, DefineList::iCalcModel, DefineList::iCalcType, CheckList::idim_max, DefineList::Ndown, DefineList::Nup, BindStruct::Phys, and stdoutMPI.

Referenced by CalcByFullDiag().

27  {
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 }
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
#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
double * all_doublon
[CheckList::idim_max+1] Doublon for FullDiag and LOBPCG. malloc in setmem_large().
Definition: struct.h:373
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
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ outputHam()

int outputHam ( struct BindStruct X)

output Hamiltonian only used for FullDiag mode

Note
global: [in] Ham
Parameters
[in]XStruct to get information about file header names, dimension of hirbert space.
Return values
0normally finished.
-1abnormally finished.

Definition at line 73 of file output.c.

References DefineList::CDataFileHead, cFileNamePhys_FullDiag_Ham, BindStruct::Check, childfopenMPI(), D_FileNameMax, BindStruct::Def, Ham, and CheckList::idim_max.

Referenced by CalcByFullDiag().

73  {
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 }
struct DefineList Def
Definision of system (Hamiltonian) etc.
Definition: struct.h:410
unsigned long int idim_max
The dimension of the Hilbert space of this process.
Definition: struct.h:303
#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
double complex ** Ham
Definition: global.h:73
const char * cFileNamePhys_FullDiag_Ham
Definition: global.c:78
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
+ Here is the call graph for this function:
+ Here is the caller graph for this function: