HΦ  3.1.0
phys.c File Reference

File for giving a parent function to calculate physical quantities by full diagonalization method. More...

#include "phys.h"
#include "expec_energy_flct.h"
#include "expec_totalspin.h"
#include "expec_cisajs.h"
#include "expec_cisajscktaltdc.h"
#include "wrapperMPI.h"

Go to the source code of this file.

Functions

void phys (struct BindStruct *X, unsigned long int neig)
 A main function to calculate physical quantities by full diagonalization method. More...
 

Detailed Description

File for giving a parent function to calculate physical quantities by full diagonalization method.

Version
0.1, 0.2
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)

Definition in file phys.c.

Function Documentation

◆ phys()

void phys ( struct BindStruct X,
unsigned long int  neig 
)

A main function to calculate physical quantities by full diagonalization method.

Parameters
[in,out]XCalcStruct list for getting and pushing calculation information
neignumber of eigenvalues
Version
0.2

add output process of calculation results for general spin

Version
0.1
Author
Takahiro Misawa (The University of Tokyo)
Kazuyoshi Yoshimi (The University of Tokyo)
Parameters
[in,out]X
[in]neig

Definition at line 50 of file phys.c.

References exitMPI(), expec_cisajs(), expec_cisajscktaltdc(), expec_energy_flct(), expec_totalspin(), L_vec, stdoutMPI, v0, v1, and X.

Referenced by CalcByFullDiag(), and CalcByLOBPCG().

52  {
53 
54  long unsigned int i, j, i_max;
55  double tmp_N;
56 #ifdef _SCALAPACK
57  double complex *vec_tmp;
58  int ictxt, ierr, rank;
59 #endif
60 
61  i_max = X->Check.idim_max;
62  for (i = 0; i < neig; i++) {
63 #ifdef _SCALAPACK
64  if(use_scalapack){
65  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
66  vec_tmp = malloc(i_max*sizeof(double complex));
67  GetEigenVector(i, i_max, Z_vec, descZ_vec, vec_tmp);
68  if(rank == 0) {
69  for (j = 0; j < i_max; j++) {
70  v0[j + 1] = vec_tmp[j];
71  }
72  }
73  free(vec_tmp);
74  } else for (j = 0; j < i_max; j++) {
75  v0[j + 1] = L_vec[i][j];
76  }
77 #else
78  for (j = 0; j < i_max; j++) {
79  v0[j + 1] = L_vec[i][j];
80  }
81 #endif
82  X->Phys.eigen_num = i;
83  if (expec_energy_flct(X) != 0) {
84  fprintf(stderr, "Error: calc expec_energy.\n");
85  exitMPI(-1);
86  }
87  if (expec_cisajs(X, v1) != 0) {
88  fprintf(stderr, "Error: calc OneBodyG.\n");
89  exitMPI(-1);
90  }
91  if (expec_cisajscktaltdc(X, v1) != 0) {
92  fprintf(stderr, "Error: calc TwoBodyG.\n");
93  exitMPI(-1);
94  }
95 #ifdef _SCALAPACK
96  if(use_scalapack){
97  if (X->Def.iCalcType == FullDiag) {
98  fprintf(stderr, "In scalapack fulldiag, total spin is not calculated !\n");
99  X->Phys.s2=0.0;
100  X->Phys.Sz=0.0;
101  }
102  }else{
103  if (X->Def.iCalcType == FullDiag) {
104  if (expec_totalspin(X, v1) != 0) {
105  fprintf(stderr, "Error: calc TotalSpin.\n");
106  exitMPI(-1);
107  }
108  }
109  }
110 #else
111  if (X->Def.iCalcType == FullDiag) {
112  if (expec_totalspin(X, v1) != 0) {
113  fprintf(stderr, "Error: calc TotalSpin.\n");
114  exitMPI(-1);
115  }
116  }
117 #endif
118 
119  if (X->Def.iCalcModel == Spin || X->Def.iCalcModel == SpinGC) {
120  tmp_N = X->Def.NsiteMPI;
121  } else {
122  tmp_N = X->Phys.num_up + X->Phys.num_down;
123  }
124 
125  if (X->Def.iCalcType == FullDiag){
126 #ifdef _SCALAPACK
127  if (use_scalapack){
128  fprintf(stdoutMPI, "i=%5ld Energy=%10lf N=%10lf Sz=%10lf Doublon=%10lf \n", i, X->Phys.energy, tmp_N,
129  X->Phys.Sz, X->Phys.doublon);
130  }
131  else{
132  fprintf(stdoutMPI, "i=%5ld Energy=%10lf N=%10lf Sz=%10lf S2=%10lf Doublon=%10lf \n", i, X->Phys.energy, tmp_N,
133  X->Phys.Sz, X->Phys.s2, X->Phys.doublon);
134  }
135 #else
136  fprintf(stdoutMPI, "i=%5ld Energy=%10lf N=%10lf Sz=%10lf S2=%10lf Doublon=%10lf \n", i, X->Phys.energy, tmp_N,
137  X->Phys.Sz, X->Phys.s2, X->Phys.doublon);
138 
139 #endif
140  }
141  else if (X->Def.iCalcType == CG)
142  fprintf(stdoutMPI, "i=%5ld Energy=%10lf N=%10lf Sz=%10lf Doublon=%10lf \n", i, X->Phys.energy, tmp_N,
143  X->Phys.Sz, X->Phys.doublon);
144  X->Phys.all_energy[i] = X->Phys.energy;
145  X->Phys.all_doublon[i] = X->Phys.doublon;
146  X->Phys.all_sz[i] = X->Phys.Sz;
147  X->Phys.all_s2[i] = X->Phys.s2;
148  X->Phys.all_num_up[i] = X->Phys.num_up;
149  X->Phys.all_num_down[i] = X->Phys.num_down;
150  }
151 }
void exitMPI(int errorcode)
MPI Abortation wrapper.
Definition: wrapperMPI.c:86
int expec_cisajscktaltdc(struct BindStruct *X, double complex *vec)
Parent function to calculate two-body green&#39;s functions.
int expec_totalspin(struct BindStruct *X, double complex *vec)
Parent function of calculation of total spin.
double complex * v1
Definition: global.h:35
int expec_cisajs(struct BindStruct *X, double complex *vec)
function of calculation for one body green&#39;s function
Definition: expec_cisajs.c:69
double complex ** L_vec
Definition: global.h:74
struct EDMainCalStruct X
Definition: struct.h:431
int expec_energy_flct(struct BindStruct *X)
Parent function to calculate expected values of energy and physical quantities.
double complex * v0
Definition: global.h:34
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:164