HΦ  3.1.0
xsetmem.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 #include "Common.h"
18 #include "mfmemory.h"
19 #include "xsetmem.h"
20 #include "wrapperMPI.h"
21 
34 static unsigned long int mfint[7];/*for malloc*/
35 
41 void setmem_HEAD
42 (
43  struct BindStruct *X
44  )
45 {
46  X->Def.CDataFileHead = (char*)malloc(D_FileNameMax*sizeof(char));
47  X->Def.CParaFileHead = (char*)malloc(D_FileNameMax*sizeof(char));
48 }
49 
55 void setmem_def
56 (
57  struct BindStruct *X,
58  struct BoostList *xBoost
59  )
60 {
61  unsigned long int i=0;
62  unsigned long int j=0;
63  unsigned long int k=0;
64  lui_malloc1(X->Def.Tpow, 2*X->Def.Nsite+2);
65  lui_malloc1(X->Def.OrgTpow, 2*X->Def.Nsite+2);
66  for(i=0; i<2*X->Def.Nsite+2; i++){
67  X->Def.Tpow[i]=0;
68  X->Def.OrgTpow[i]=0;
69  }
70  li_malloc1(X->Def.SiteToBit, X->Def.Nsite+1);
71  for(i=0; i<X->Def.Nsite+1; i++){
72  X->Def.SiteToBit[i]=0;
73  }
74 
75  i_malloc1(X->Def.LocSpn, X->Def.Nsite);
76  d_malloc1(X->Phys.spin_real_cor, X->Def.Nsite*X->Def.Nsite);
77  d_malloc1(X->Phys.charge_real_cor, X->Def.Nsite*X->Def.Nsite);
78  d_malloc1(X->Phys.loc_spin_z, X->Def.Nsite*X->Def.Nsite);
79 
80  i_malloc1(X->Def.EDChemi, X->Def.EDNChemi+X->Def.NInterAll+X->Def.NTransfer);
81  i_malloc1(X->Def.EDSpinChemi, X->Def.EDNChemi+X->Def.NInterAll+X->Def.NTransfer);
82  d_malloc1(X->Def.EDParaChemi, X->Def.EDNChemi+X->Def.NInterAll+X->Def.NTransfer);
83 
84  i_malloc2(X->Def.GeneralTransfer, X->Def.NTransfer, 4);
85  c_malloc1(X->Def.ParaGeneralTransfer, X->Def.NTransfer);
86 
87  if(X->Def.iCalcType == TimeEvolution){
88  i_malloc2(X->Def.EDGeneralTransfer, X->Def.NTransfer+X->Def.NTETransferMax, 4);
89  c_malloc1(X->Def.EDParaGeneralTransfer, X->Def.NTransfer+X->Def.NTETransferMax);
90 
91  }
92  else {
93  i_malloc2(X->Def.EDGeneralTransfer, X->Def.NTransfer, 4);
94  c_malloc1(X->Def.EDParaGeneralTransfer, X->Def.NTransfer);
95  }
96  i_malloc2(X->Def.CoulombIntra, X->Def.NCoulombIntra, 1);
97  d_malloc1(X->Def.ParaCoulombIntra, X->Def.NCoulombIntra);
98  i_malloc2(X->Def.CoulombInter, X->Def.NCoulombInter+X->Def.NIsingCoupling, 2);
99  d_malloc1(X->Def.ParaCoulombInter, X->Def.NCoulombInter+X->Def.NIsingCoupling);
100  i_malloc2(X->Def.HundCoupling, X->Def.NHundCoupling+X->Def.NIsingCoupling, 2);
101  d_malloc1(X->Def.ParaHundCoupling, X->Def.NHundCoupling+X->Def.NIsingCoupling);
102  i_malloc2(X->Def.PairHopping, X->Def.NPairHopping, 2);
103  d_malloc1(X->Def.ParaPairHopping, X->Def.NPairHopping);
104  i_malloc2(X->Def.ExchangeCoupling, X->Def.NExchangeCoupling, 2);
105  d_malloc1(X->Def.ParaExchangeCoupling, X->Def.NExchangeCoupling);
106  i_malloc2(X->Def.PairLiftCoupling, X->Def.NPairLiftCoupling, 2);
107  d_malloc1(X->Def.ParaPairLiftCoupling, X->Def.NPairLiftCoupling);
108 
109  i_malloc2(X->Def.InterAll, X->Def.NInterAll, 8);
110  c_malloc1(X->Def.ParaInterAll, X->Def.NInterAll);
111 
112  i_malloc2(X->Def.CisAjt, X->Def.NCisAjt, 4);
113  i_malloc2(X->Def.CisAjtCkuAlvDC, X->Def.NCisAjtCkuAlvDC, 8);
114 
115  i_malloc2(X->Def.SingleExcitationOperator, X->Def.NSingleExcitationOperator, 3);
116  c_malloc1(X->Def.ParaSingleExcitationOperator, X->Def.NSingleExcitationOperator);
117  i_malloc2(X->Def.PairExcitationOperator, X->Def.NPairExcitationOperator, 5);
118  c_malloc1(X->Def.ParaPairExcitationOperator, X->Def.NPairExcitationOperator);
119 
120  d_malloc1(X->Def.ParaLaser, X->Def.NLaser);
121 
122  unsigned int ipivot,iarrayJ,ispin;
123  xBoost->list_6spin_star = (int **)malloc(sizeof(int*) * xBoost->R0 * xBoost->num_pivot);
124  for (ipivot = 0; ipivot < xBoost->R0 * xBoost->num_pivot; ipivot++) {
125  xBoost->list_6spin_star[ipivot] = (int *)malloc(sizeof(int) * 7);
126  }
127 
128  xBoost->list_6spin_pair = (int ***)malloc(sizeof(int**) * xBoost->R0 * xBoost->num_pivot);
129  for (ipivot = 0; ipivot < xBoost->R0 * xBoost->num_pivot; ipivot++) {
130  xBoost->list_6spin_pair[ipivot] = (int **)malloc(sizeof(int*) * 7);
131  for (ispin = 0; ispin < 7; ispin++) {
132  xBoost->list_6spin_pair[ipivot][ispin] = (int *)malloc(sizeof(int) * 15);
133  }
134  }
135 
136  xBoost->arrayJ = (double complex ***)malloc(sizeof(double complex**) * xBoost->NumarrayJ);
137 for (iarrayJ = 0; iarrayJ < xBoost->NumarrayJ; iarrayJ++) {
138  xBoost->arrayJ[iarrayJ] = (double complex **)malloc(sizeof(double complex*) * 3);
139  for (i = 0; i < 3; i++) {
140  xBoost->arrayJ[iarrayJ][i] = (double complex *)malloc(sizeof(double complex) * 3);
141  }
142 }
143 
144  int NInterAllSet;
145  NInterAllSet= (X->Def.iCalcType==TimeEvolution) ? X->Def.NInterAll+X->Def.NTEInterAllMax: X->Def.NInterAll;
146  i_malloc2(X->Def.InterAll_OffDiagonal, NInterAllSet, 8);
147  c_malloc1(X->Def.ParaInterAll_OffDiagonal, NInterAllSet);
148  i_malloc2(X->Def.InterAll_Diagonal, NInterAllSet, 4);
149  d_malloc1(X->Def.ParaInterAll_Diagonal, NInterAllSet);
150 
151  if (X->Def.iCalcType == TimeEvolution){
152  d_malloc1(X->Def.TETime, X->Def.NTETimeSteps);
153  //Time-dependent Transfer
154  ui_malloc1(X->Def.NTETransfer, X->Def.NTETimeSteps);
155  ui_malloc1(X->Def.NTETransferDiagonal, X->Def.NTETimeSteps);
156  i_malloc3(X->Def.TETransfer, X->Def.NTETimeSteps, X->Def.NTETransferMax, 4);
157  i_malloc3(X->Def.TETransferDiagonal, X->Def.NTETimeSteps, X->Def.NTETransferMax, 2);
158  c_malloc2(X->Def.ParaTETransfer, X->Def.NTETimeSteps, X->Def.NTETransferMax);
159  d_malloc2(X->Def.ParaTETransferDiagonal, X->Def.NTETimeSteps,X->Def.NTETransferMax);
160  //Time-dependent InterAll
161  ui_malloc1(X->Def.NTEInterAll, X->Def.NTETimeSteps);
162  ui_malloc1(X->Def.NTEInterAllDiagonal, X->Def.NTETimeSteps);
163  i_malloc3(X->Def.TEInterAll, X->Def.NTETimeSteps, X->Def.NTEInterAllMax, 8);
164  i_malloc3(X->Def.TEInterAllDiagonal, X->Def.NTETimeSteps, X->Def.NTEInterAllMax, 4);
165  c_malloc2(X->Def.ParaTEInterAll, X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
166  d_malloc2(X->Def.ParaTEInterAllDiagonal, X->Def.NTETimeSteps,X->Def.NTEInterAllMax);
167  ui_malloc1(X->Def.NTEInterAllOffDiagonal, X->Def.NTETimeSteps);
168  i_malloc3(X->Def.TEInterAllOffDiagonal, X->Def.NTETimeSteps, X->Def.NTEInterAllMax, 8);
169  c_malloc2(X->Def.ParaTEInterAllOffDiagonal, X->Def.NTETimeSteps,X->Def.NTEInterAllMax);
170 
171  //Time-dependent Chemi generated by InterAll diagonal components
172  ui_malloc1(X->Def.NTEChemi, X->Def.NTETimeSteps);
173  i_malloc2(X->Def.TEChemi, X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
174  i_malloc2(X->Def.SpinTEChemi, X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
175  d_malloc2(X->Def.ParaTEChemi, X->Def.NTETimeSteps, X->Def.NTEInterAllMax);
176 
177  for(i = 0; i < X->Def.NTETimeSteps; i++){
178  X->Def.TETime[i]=0;
179  X->Def.NTETransfer[i]=0;
180  X->Def.NTETransferDiagonal[i]=0;
181  X->Def.NTEChemi[i]=0;
182 
183  X->Def.NTEInterAll[i] = 0;
184  X->Def.NTEInterAllDiagonal[i] = 0;
185  X->Def.NTEInterAllOffDiagonal[i] = 0;
186  for(j = 0; j < X->Def.NTETransferMax; j++) {
187  X->Def.ParaTETransfer[i][j]=0;
188  X->Def.ParaTETransferDiagonal[i][j]=0;
189  for(k = 0; k< 4; k++){
190  X->Def.TETransfer[i][j][k]=0;
191  }
192  for(k = 0; k< 2; k++){
193  X->Def.TETransferDiagonal[i][j][k]=0;
194  }
195  }
196  for(j = 0; j < X->Def.NTEInterAllMax; j++){
197  X->Def.ParaTEInterAll[i][j]=0;
198  X->Def.ParaTEInterAllDiagonal[i][j]=0;
199  X->Def.ParaTEInterAllOffDiagonal[i][j]=0;
200  X->Def.TEChemi[i][j]=0;
201  X->Def.SpinTEChemi[i][j]=0;
202  X->Def.ParaTEChemi[i][j]=0;
203  for(k = 0; k< 4; k++){
204  X->Def.TEInterAllDiagonal[i][j][k]=0;
205  }
206  for(k = 0; k< 8; k++){
207  X->Def.TEInterAll[i][j][k]=0;
208  X->Def.TEInterAllOffDiagonal[i][j][k]=0;
209  }
210  }
211  }
212 
213  }
214 }
215 
216 
223 int setmem_large
224 (
225  struct BindStruct *X
226  ) {
227 
228  unsigned long int j = 0;
229  unsigned long int idim_maxMPI;
230 
231  idim_maxMPI = MaxMPI_li(X->Check.idim_max);
232 
233  if (GetlistSize(X) == TRUE) {
234  lui_malloc1(list_1, X->Check.idim_max + 1);
235 #ifdef MPI
236  lui_malloc1(list_1buf, idim_maxMPI + 1);
237  for (j = 0; j < X->Check.idim_max + 1; j++) {
238  list_1buf[j] = 0;
239  }
240 #endif // MPI
241  lui_malloc1(list_2_1, X->Large.SizeOflist_2_1);
242  lui_malloc1(list_2_2, X->Large.SizeOflist_2_2);
243  if (list_1 == NULL
244  || list_2_1 == NULL
245  || list_2_2 == NULL
246  ) {
247  return -1;
248  }
249  for (j = 0; j < X->Check.idim_max + 1; j++) {
250  list_1[j] = 0;
251  }
252  for (j = 0; j < X->Large.SizeOflist_2_1; j++) {
253  list_2_1[j] = 0;
254  }
255  for (j = 0; j < X->Large.SizeOflist_2_2; j++) {
256  list_2_2[j] = 0;
257  }
258  }
259 
260  d_malloc1(list_Diagonal, X->Check.idim_max + 1);
261  c_malloc1(v0, X->Check.idim_max + 1);
262  c_malloc1(v1, X->Check.idim_max + 1);
263  for (j = 0; j < X->Check.idim_max + 1; j++) {
264  list_Diagonal[j] = 0;
265  v0[j] = 0;
266  v1[j] = 0;
267  }
268  if (X->Def.iCalcType == TimeEvolution) {
269  c_malloc1(v2, X->Check.idim_max + 1);
270  } else {
271  c_malloc1(v2, 1);
272  }
273 #ifdef MPI
274  c_malloc1(v1buf, idim_maxMPI + 1);
275  for (j = 0; j < X->Check.idim_max + 1; j++) {
276  v1buf[j] = 0;
277  }
278 #endif // MPI
279  if (X->Def.iCalcType == TPQCalc) {
280  c_malloc1(vg, 1);
281  vg[0] = 0;
282  } else {
283  c_malloc1(vg, X->Check.idim_max + 1);
284  for (j = 0; j < X->Check.idim_max + 1; j++) {
285  vg[j] = 0;
286  }
287  }
288  d_malloc1(alpha, X->Def.Lanczos_max + 1);
289  d_malloc1(beta, X->Def.Lanczos_max + 1);
290 
291  if (
292  list_Diagonal == NULL
293  || v0 == NULL
294  || v1 == NULL
295  || vg == NULL
296  ) {
297  return -1;
298  }
299 
300  if (X->Def.iCalcType == TPQCalc || X->Def.iFlgCalcSpec != CALCSPEC_NOT) {
301  c_malloc2(vec, X->Def.Lanczos_max + 1, X->Def.Lanczos_max + 1);
302  } else if (X->Def.iCalcType == Lanczos || X->Def.iCalcType == CG) {
303  if (X->Def.LanczosTarget > X->Def.nvec) {
304  c_malloc2(vec, X->Def.LanczosTarget + 1, X->Def.Lanczos_max + 1);
305  } else {
306  c_malloc2(vec, X->Def.nvec + 1, X->Def.Lanczos_max + 1);
307  }
308  }
309 
310  if (X->Def.iCalcType == FullDiag) {
311  d_malloc1(X->Phys.all_num_down, X->Check.idim_max + 1);
312  d_malloc1(X->Phys.all_num_up, X->Check.idim_max + 1);
313  d_malloc1(X->Phys.all_energy, X->Check.idim_max + 1);
314  d_malloc1(X->Phys.all_doublon, X->Check.idim_max + 1);
315  d_malloc1(X->Phys.all_sz, X->Check.idim_max + 1);
316  d_malloc1(X->Phys.all_s2, X->Check.idim_max + 1);
317  c_malloc2(Ham, X->Check.idim_max + 1, X->Check.idim_max + 1);
318  c_malloc2(L_vec, X->Check.idim_max + 1, X->Check.idim_max + 1);
319 
320  if (X->Phys.all_num_down == NULL
321  || X->Phys.all_num_up == NULL
322  || X->Phys.all_energy == NULL
323  || X->Phys.all_doublon == NULL
324  || X->Phys.all_s2 == NULL
325  ) {
326  return -1;
327  }
328  for (j = 0; j < X->Check.idim_max + 1; j++) {
329  if (Ham[j] == NULL || L_vec[j] == NULL) {
330  return -1;
331  }
332  }
333  } else if (X->Def.iCalcType == CG) {
334  d_malloc1(X->Phys.all_num_down, X->Def.k_exct);
335  d_malloc1(X->Phys.all_num_up, X->Def.k_exct);
336  d_malloc1(X->Phys.all_energy, X->Def.k_exct);
337  d_malloc1(X->Phys.all_doublon, X->Def.k_exct);
338  d_malloc1(X->Phys.all_sz, X->Def.k_exct);
339  d_malloc1(X->Phys.all_s2, X->Def.k_exct);
340  }
341  fprintf(stdoutMPI, "%s", cProFinishAlloc);
342  return 0;
343 }
344 
355  (
356  int **InterAllOffDiagonal,
357  double complex *ParaInterAllOffDiagonal,
358  int **InterAllDiagonal,
359  double *ParaInterAllDiagonal,
360  const int NInterAll
361  ) {
362  i_malloc2(InterAllOffDiagonal, NInterAll, 8);
363  c_malloc1(ParaInterAllOffDiagonal, NInterAll);
364  i_malloc2(InterAllDiagonal, NInterAll, 4);
365  d_malloc1(ParaInterAllDiagonal, NInterAll);
366  }
367 
368 
378  int GetlistSize
379  (
380  struct BindStruct *X
381  ) {
382  // unsigned int idim_maxMPI;
383 
384 // idim_maxMPI = MaxMPI_li(X->Check.idim_max);
385 
386  switch (X->Def.iCalcModel) {
387  case Spin:
388  case Hubbard:
389  case HubbardNConserved:
390  case Kondo:
391  case KondoGC:
392  if (X->Def.iFlgGeneralSpin == FALSE) {
393  if (X->Def.iCalcModel == Spin && X->Def.Nsite % 2 == 1) {
394  X->Large.SizeOflist_2_1 = X->Check.sdim * 2 + 2;
395  } else {
396  X->Large.SizeOflist_2_1 = X->Check.sdim + 2;
397  }
398  X->Large.SizeOflist_2_2 = X->Check.sdim + 2;
399  X->Large.SizeOflistjb = X->Check.sdim + 2;
400  } else {//for spin-canonical general spin
401  X->Large.SizeOflist_2_1 = X->Check.sdim + 2;
402  X->Large.SizeOflist_2_2 =
403  X->Def.Tpow[X->Def.Nsite - 1] * X->Def.SiteToBit[X->Def.Nsite - 1] / X->Check.sdim + 2;
404  X->Large.SizeOflistjb =
405  X->Def.Tpow[X->Def.Nsite - 1] * X->Def.SiteToBit[X->Def.Nsite - 1] / X->Check.sdim + 2;
406  }
407  break;
408  default:
409  return FALSE;
410  }
411  return TRUE;
412 }
int ** list_6spin_star
Definition: struct.h:402
double complex * v1
Definition: global.h:35
#define D_FileNameMax
Definition: global.h:23
double complex * v2
Definition: global.h:36
double complex ** vec
Definition: global.h:45
long unsigned int * list_1buf
Definition: global.h:48
long unsigned int num_pivot
Definition: struct.h:397
int *** list_6spin_pair
Definition: struct.h:403
void setmem_def(struct BindStruct *X, struct BoostList *xBoost)
Set size of memories for Def and Phys in BindStruct.
Definition: xsetmem.c:56
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
Bind.
Definition: struct.h:408
double complex ** Ham
Definition: global.h:73
double complex *** arrayJ
Definition: struct.h:400
const char * cProFinishAlloc
unsigned long int MaxMPI_li(unsigned long int idim)
MPI wrapper function to obtain maximum unsigned long integer across processes.
Definition: wrapperMPI.c:171
long unsigned int * list_2_1
Definition: global.h:49
static unsigned long int mfint[7]
Definition: xsetmem.c:34
void setmem_HEAD(struct BindStruct *X)
Set size of memories headers of output files.
Definition: xsetmem.c:42
int GetlistSize(struct BindStruct *X)
Set size of lists for the canonical ensemble.
Definition: xsetmem.c:379
long unsigned int R0
Definition: struct.h:395
double complex ** L_vec
Definition: global.h:74
double * alpha
Definition: global.h:44
#define FALSE
Definition: global.h:25
long unsigned int * list_1
Definition: global.h:47
double complex * vg
Definition: global.h:41
long unsigned int * list_2_2
Definition: global.h:50
double * beta
Definition: global.h:44
For Boost.
Definition: struct.h:393
struct EDMainCalStruct X
Definition: struct.h:431
double * list_Diagonal
Definition: global.h:46
void setmem_IntAll_Diagonal(int **InterAllOffDiagonal, double complex *ParaInterAllOffDiagonal, int **InterAllDiagonal, double *ParaInterAllDiagonal, const int NInterAll)
Set the size of memories for InterAllDiagonal and InterAllOffDiagonal arrays.
Definition: xsetmem.c:355
unsigned int NumarrayJ
Definition: struct.h:399
double complex * v0
Definition: global.h:34
double complex * v1buf
Definition: global.h:37
FILE * stdoutMPI
File pointer to the standard output defined in InitializeMPI()
Definition: global.h:164