You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

ABACUS_Scan.h 37KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844
  1. /**********************************************************
  2. This software is part of J.-S. Caux's ABACUS library.
  3. Copyright (c) J.-S. Caux.
  4. -----------------------------------------------------------
  5. File: ABACUS_Scan.h
  6. Purpose: Declares all classes and functions used in the
  7. ABACUS logic of scanning with threads.
  8. ***********************************************************/
  9. #ifndef ABACUS_SCAN_H
  10. #define ABACUS_SCAN_H
  11. #include "ABACUS.h"
  12. namespace ABACUS {
  13. const int MAX_STATE_LIST_SIZE = 10000;
  14. // Functions in src/UTILS/Data_File_Name.cc:
  15. void Data_File_Name (std::stringstream& name, char whichDSF, DP c_int, DP L, int N,
  16. int iKmin, int iKmax, DP kBT, DP L2, std::string defaultname);
  17. void Data_File_Name (std::stringstream& name, char whichDSF, int iKmin, int iKmax, DP kBT,
  18. LiebLin_Bethe_State& State, LiebLin_Bethe_State& RefScanState, std::string defaultname);
  19. void Data_File_Name (std::stringstream& name, char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  20. DP kBT, int N2, std::string defaultname);
  21. void Data_File_Name (std::stringstream& name, char whichDSF, int iKmin, int iKmax, DP kBT,
  22. Heis_Bethe_State& State, Heis_Bethe_State& RefScanState, std::string defaultname);
  23. void ODSLF_Data_File_Name (std::stringstream& name, char whichDSF, DP Delta, int N, int M,
  24. int iKmin, int iKmax, DP kBT, int N2, std::string defaultname);
  25. void Data_File_Name (std::stringstream& name, char whichDSF, int iKmin, int iKmax, DP kBT,
  26. ODSLF_Bethe_State& State, ODSLF_Bethe_State& RefScanState, std::string defaultname);
  27. // Coding to convert ints to strings: for application in reduced labels
  28. //const int ABACUScodingsize = 64; // use a multiple of 2 to accelerate divisions in labeling.
  29. //const char ABACUScoding[] = {'a', 'b', 'c', 'd', 'e', 'f', 'g', 'h', 'i', 'j', 'k', 'l', 'm', 'n', 'o', 'p', 'q', 'r', 's', 't', 'u', 'v', 'w', 'x', 'y', 'z', 'A', 'B', 'C', 'D', 'E', 'F', 'G', 'H', 'I', 'J', 'K', 'L', 'M', 'N', 'O', 'P', 'Q', 'R', 'S', 'T', 'U', 'V', 'W', 'X', 'Y', 'Z', '0', '1', '2', '3', '4', '5', '6', '7', '8', '9', '!', '?'};
  30. // From ABACUS++T_8 onwards: forbid special characters as |, :, !, ? and all capital letters in labels.
  31. // This is due to the dumb capitalization-preserving but capitalization-insensitive HFS+ filesystem on Mac OS X.
  32. const int ABACUScodingsize = 32;
  33. const char ABACUScoding[] = {'0', '1', '2', '3', '4', '5', '6', '7', '8', '9', 'a', 'b', 'c', 'd', 'e', 'f', 'g', 'h', 'i', 'j', 'k', 'l', 'm', 'n', 'o', 'p', 'q', 'r', 's', 't', 'u', 'v'};
  34. const char LABELSEP = '_'; // was _
  35. const char TYPESEP = 'x'; // was |
  36. const char EXCSEP = 'y'; // was :
  37. const char INEXCSEP = 'z'; // was @
  38. struct State_Label_Data {
  39. Vect<int> type; // integer type labels of the types present
  40. Vect<int> M; // how many particles of each type
  41. Vect<int> nexc; // how many excitations as compared to the reference state used
  42. Vect<Vect<int> > Ix2old; // which Ix2 will be excited
  43. Vect<Vect<int> > Ix2exc; // which Ix2 the excitation has shifted to
  44. State_Label_Data (const Vect<int>& type_ref, const Vect<int>& M_ref,
  45. const Vect<int>& nexc_ref, const Vect<Vect<int> >& Ix2old_ref, const Vect<Vect<int> >& Ix2exc_ref)
  46. {
  47. type = type_ref; M = M_ref; nexc = nexc_ref; Ix2old = Ix2old_ref; Ix2exc = Ix2exc_ref;
  48. }
  49. };
  50. std::string Extract_Base_Label (std::string label); // works for labels and complabels
  51. std::string Extract_nexc_Label (std::string label);
  52. // For compressed labels: conversions between integers and char/strings.
  53. std::string Convert_POSINT_to_STR (int int_to_convert);
  54. int Convert_CHAR_to_POSINT (char char_to_convert);
  55. int Convert_STR_to_POSINT (std::string str_to_convert);
  56. State_Label_Data Read_Base_Label (std::string label);
  57. State_Label_Data Read_State_Label (std::string label, const Vect<Vect<int> >& OriginIx2);
  58. State_Label_Data Read_State_Label (std::string label, const Vect<int>& OriginIx2); // if there is only one type
  59. std::string Return_State_Label (State_Label_Data data, const Vect<Vect<int> >& OriginIx2);
  60. std::string Return_State_Label (State_Label_Data data, const Vect<int>& OriginIx2); // if there is only one type
  61. std::string Return_State_Label (const Vect<Vect<int> >& ScanIx2, const Vect<Vect<int> >& OriginIx2);
  62. std::string Return_State_Label (const Vect<int>& ScanIx2, const Vect<int>& OriginIx2); // if there is only one type
  63. Vect<Vect<int> > Return_Ix2_from_Label (std::string label_ref, const Vect<Vect<int> >& OriginIx2);
  64. Vect<int> Return_Ix2_from_Label (std::string label_ref, const Vect<int>& OriginIx2); // specialization to Lieb-Liniger
  65. // Functions for descending states: in SCAN/Descendents.cc
  66. Vect<std::string> Descendent_States_with_iK_Stepped_Up (std::string ScanIx2_label, const LiebLin_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  67. Vect<std::string> Descendent_States_with_iK_Stepped_Down (std::string ScanIx2_label, const LiebLin_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  68. Vect<std::string> Descendent_States_with_iK_Preserved (std::string ScanIx2_label, const LiebLin_Bethe_State& OriginState, bool disperse_only_current_exc_up, bool preserve_nexc_up, bool disperse_only_current_exc_down, bool preserve_nexc_down);
  69. Vect<std::string> Descendent_States_with_iK_Stepped_Up (std::string ScanIx2_label, const Heis_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  70. Vect<std::string> Descendent_States_with_iK_Stepped_Down (std::string ScanIx2_label, const Heis_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  71. Vect<std::string> Descendent_States_with_iK_Preserved (std::string ScanIx2_label, const Heis_Bethe_State& OriginState, bool disperse_only_current_exc_up, bool preserve_nexc_up, bool disperse_only_current_exc_down, bool preserve_nexc_down);
  72. // For symmetric state scanning:
  73. Vect<std::string> Descendent_States_with_iK_Stepped_Up_rightIx2only
  74. (std::string ScanIx2_label, const LiebLin_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  75. Vect<std::string> Descendent_States_with_iK_Stepped_Down_rightIx2only
  76. (std::string ScanIx2_label, const LiebLin_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  77. Vect<std::string> Descendent_States_with_iK_Stepped_Up_rightIx2only
  78. (std::string ScanIx2_label, const Heis_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  79. Vect<std::string> Descendent_States_with_iK_Stepped_Down_rightIx2only
  80. (std::string ScanIx2_label, const Heis_Bethe_State& OriginState, bool disperse_only_current_exc, bool preserve_nexc);
  81. // Functions to prepare and wrapup parallel scans:
  82. void Prepare_Parallel_Scan_LiebLin (char whichDSF, DP c_int, DP L, int N, int iKmin, int iKmax, DP kBT,
  83. std::string defaultname, int paralevel, Vect<int> rank_lower_paralevels,
  84. Vect<int> nr_processors_lower_paralevels, int nr_processors_at_newlevel);
  85. void Wrapup_Parallel_Scan_LiebLin (char whichDSF, DP c_int, DP L, int N, int iKmin, int iKmax, DP kBT,
  86. std::string defaultname, int paralevel, Vect<int> rank_lower_paralevels,
  87. Vect<int> nr_processors_lower_paralevels, int nr_processors_at_newlevel);
  88. void Prepare_Parallel_Scan_Heis (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  89. int paralevel, Vect<int> rank_lower_paralevels,
  90. Vect<int> nr_processors_lower_paralevels, int nr_processors_at_newlevel);
  91. void Wrapup_Parallel_Scan_Heis (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  92. int paralevel, Vect<int> rank_lower_paralevels,
  93. Vect<int> nr_processors_lower_paralevels, int nr_processors_at_newlevel);
  94. void Filter_RAW_File_for_iK (const char ff_file[], int iKneeded);
  95. void Sort_RAW_File (const char ffsq_file[], char optionchar);
  96. void Sort_RAW_File (const char ffsq_file[], char optionchar, char whichDSF);
  97. // Functions for data interpretation:
  98. DP Smoothen_RAW_into_SF (std::string prefix, int iKmin, int iKmax, int DiK, DP devmax,
  99. DP ommin, DP ommax, int Nom, DP gwidth, DP normalization, DP denom_sum_K);
  100. DP Smoothen_RAW_into_SF (std::string prefix, Vect<std::string> rawfilename, Vect<DP> weight, int iKmin, int iKmax, int DiK, DP devmax,
  101. DP ommin, DP ommax, int Nom, DP gwidth, DP normalization, DP denom_sum_K);
  102. void Write_K_File (DP Length, int iKmin, int iKmax);
  103. void Write_Omega_File (int Nout_omega, DP omegamin, DP omegamax);
  104. void Write_Time_File (int Nt, DP tmin, DP tmax);
  105. // Smoothen with gaussian width scaled with two-particle bandwidth
  106. DP Smoothen_RAW_into_SF_LiebLin_Scaled (std::string prefix, DP L, int N, int iKmin, int iKmax, int DiK, DP ommin, DP ommax, int Nom, DP width, DP normalization);
  107. //****************************************************************************
  108. struct Scan_Info {
  109. DP sumrule_obtained;
  110. DP Nfull; // dimensionality of (sub)Hilbert space considered
  111. long long int Ninadm;
  112. long long int Ndata;
  113. long long int Ndata_conv;
  114. long long int Ndata_conv0;
  115. double TT; // total computation time in seconds
  116. public:
  117. Scan_Info(); // constructor, puts everything to zero
  118. Scan_Info (DP sr, DP Nf, long long int Ni, long long int Nd, long long int Ndc, long long int Ndc0, double t);
  119. void Save (const char* outfile_Cstr);
  120. void Load (const char* infile_Cstr);
  121. inline Scan_Info& operator = (const Scan_Info& ref_info)
  122. {
  123. sumrule_obtained = ref_info.sumrule_obtained;
  124. Nfull = ref_info.Nfull;
  125. Ninadm = ref_info.Ninadm;
  126. Ndata = ref_info.Ndata;
  127. Ndata_conv = ref_info.Ndata_conv;
  128. Ndata_conv0 = ref_info.Ndata_conv0;
  129. TT = ref_info.TT;
  130. return(*this);
  131. }
  132. inline Scan_Info& operator+= (const Scan_Info& ref_info)
  133. {
  134. if (this != &ref_info) {
  135. sumrule_obtained += ref_info.sumrule_obtained;
  136. Nfull += ref_info.Nfull;
  137. Ninadm += ref_info.Ninadm;
  138. Ndata += ref_info.Ndata;
  139. Ndata_conv += ref_info.Ndata_conv;
  140. Ndata_conv0 += ref_info.Ndata_conv0;
  141. TT += ref_info.TT;
  142. }
  143. return(*this);
  144. }
  145. inline Scan_Info& operator-= (const Scan_Info& ref_info)
  146. {
  147. if (this != &ref_info) {
  148. sumrule_obtained -= ref_info.sumrule_obtained;
  149. Nfull -= ref_info.Nfull;
  150. Ninadm -= ref_info.Ninadm;
  151. Ndata -= ref_info.Ndata;
  152. Ndata_conv -= ref_info.Ndata_conv;
  153. Ndata_conv0 -= ref_info.Ndata_conv0;
  154. TT -= ref_info.TT;
  155. }
  156. return(*this);
  157. }
  158. };
  159. std::ostream& operator<< (std::ostream& s, const Scan_Info& info);
  160. // Functions in src/SCAN/General_Scan.cc:
  161. template<class Tstate>
  162. Scan_Info General_Scan (char whichDSF, int iKmin, int iKmax, int iKmod, DP kBT, Tstate& AveragingState, Tstate& SeedScanState,
  163. std::string defaultScanStatename, int Max_Secs, DP target_sumrule, bool refine, int paralevel, Vect<int> rank, Vect<int> nr_processors);
  164. Scan_Info Scan_LiebLin (char whichDSF, DP c_int, DP L, int N, int iKmin, int iKmax, DP kBT,
  165. int Max_Secs, DP target_sumrule, bool refine,
  166. int paralevel, Vect<int> rank, Vect<int> nr_processors);
  167. Scan_Info Scan_LiebLin (char whichDSF, DP c_int, DP L, int N, int iKmin, int iKmax, DP kBT,
  168. int Max_Secs, DP target_sumrule, bool refine);
  169. Scan_Info Scan_LiebLin (char whichDSF, LiebLin_Bethe_State AveragingState, std::string defaultScanStatename,
  170. int iKmin, int iKmax,
  171. int Max_Secs, DP target_sumrule, bool refine,
  172. int paralevel, Vect<int> rank, Vect<int> nr_processors);
  173. Scan_Info Scan_LiebLin (char whichDSF, LiebLin_Bethe_State AveragingState, std::string defaultname,
  174. int iKmin, int iKmax,
  175. int Max_Secs, DP target_sumrule, bool refine);
  176. Scan_Info Scan_LiebLin_Geometric_Quench (DP c_int, DP L_1, int type_id_1, long long int id_1, DP L_2, int N,
  177. int iK_UL, int Max_Secs, DP target_sumrule, bool refine);
  178. Scan_Info Scan_Heis (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  179. int Max_Secs, DP target_sumrule, bool refine,
  180. int paralevel, Vect<int> rank, Vect<int> nr_processors);
  181. Scan_Info Scan_Heis (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  182. int Max_Secs, DP target_sumrule, bool refine);
  183. Scan_Info Scan_Heis (char whichDSF, XXZ_Bethe_State& AveragingState, std::string defaultScanStatename,
  184. int iKmin, int iKmax,
  185. int Max_Secs, DP target_sumrule, bool refine,
  186. int paralevel, Vect<int> rank, Vect<int> nr_processors);
  187. Scan_Info Scan_Heis (char whichDSF, XXX_Bethe_State& AveragingState, std::string defaultScanStatename,
  188. int iKmin, int iKmax,
  189. int Max_Secs, DP target_sumrule, bool refine,
  190. int paralevel, Vect<int> rank, Vect<int> nr_processors);
  191. Scan_Info Scan_Heis (char whichDSF, XXZ_gpd_Bethe_State& AveragingState, std::string defaultScanStatename,
  192. int iKmin, int iKmax,
  193. int Max_Secs, DP target_sumrule, bool refine,
  194. int paralevel, Vect<int> rank, Vect<int> nr_processors);
  195. Scan_Info Scan_ODSLF (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  196. int Max_Secs, DP target_sumrule, bool refine, int rank, int nr_processors);
  197. Scan_Info Scan_ODSLF (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax, int Max_Secs, bool refine);
  198. Scan_Info Scan_ODSLF (char whichDSF, DP Delta, int N, int M, int iKneeded, int Max_Secs, bool refine);
  199. Scan_Info Scan_ODSLF (char whichDSF, DP Delta, int N, int M, int Max_Secs, bool refine);
  200. //****************************************************************************
  201. // Functions in src/SCAN/Descendents.cc:
  202. Vect<std::string> Descendents (const LiebLin_Bethe_State& ScanState, const LiebLin_Bethe_State& OriginState, int type_required);
  203. Vect<std::string> Descendents (const Heis_Bethe_State& ScanState, const Heis_Bethe_State& OriginState, int type_required);
  204. struct Scan_Thread {
  205. std::string label;
  206. int type;
  207. Scan_Thread ();
  208. Scan_Thread (std::string label_ref, int type_ref) {
  209. label = label_ref;
  210. type = type_ref;
  211. }
  212. Scan_Thread& operator= (const Scan_Thread& RefThread);
  213. };
  214. struct Scan_Thread_Data {
  215. // By convention, a Scan_Thread_Data object handles a list of threads which are yet to be descended.
  216. // Improvement on Scan_Thread_Set used up to ABACUS++G_7, saving data to disk instead of holding it in memory.
  217. int nlists = 6400; // number of threads lists, fixed to this number by convention.
  218. DP logscale = (1.0/64) * log(2.0); // each separate list contains threads differing by a scale factor of 2^{1/64} \approx 1.01.
  219. std::string thrdir_name; // directory in which threads files are saved.
  220. Vect<int> nthreads_total;
  221. Vect<int> nthreads_on_disk;
  222. int lowest_il_with_nthreads_neq_0;
  223. // In-memory storage, for adding threads efficiently without constantly writing to disk
  224. // These objects are saved to disk when Next_Scan_Threads are called.
  225. Vect<int> dim;
  226. Vect<int> nthreads_in_memory;
  227. Vect<Vect<std::string> > label;
  228. Vect<Vect<int> > type; // which type of descendent is needed
  229. Vect<std::string> filename;
  230. Scan_Thread_Data ();
  231. Scan_Thread_Data (std::string thrdir_name_ref, bool refine);
  232. ~Scan_Thread_Data ();
  233. bool Increase_Memory_Size (int il, int nr_to_add);
  234. void Include_Thread (DP abs_data_value_ref, std::string label_ref, int type_ref);
  235. void Include_Thread (int il, std::string label_ref, int type_ref);
  236. Vect<Scan_Thread> Extract_Next_Scan_Threads (); // returns a vector of the threads that are next in line. By defn, all threads with index il == lowest_il_with_nthreads_neq_0. These are removed from the object.
  237. Vect<Scan_Thread> Extract_Next_Scan_Threads (int min_nr); // as above, but returns a minimum of min_nr threads.
  238. void Flush_to_Disk (int il);
  239. void Save ();
  240. void Load ();
  241. };
  242. //****************************************************************************
  243. // To populate a list of states for scanning:
  244. inline void Scan_for_Possible_Bases (const Vect<int> SeedNrap, const Vect<int> Str_L,
  245. int Mdown_remaining, Vect<std::string>& possible_base_label, int& nfound, int nexc_max_used,
  246. int base_level_to_scan, Vect<int>& Nrapidities)
  247. {
  248. if (Mdown_remaining < 0) { ABACUSerror("Scan_for_Possible_Bases: shouldn't be here..."); } // reached inconsistent point
  249. if (base_level_to_scan == 0) {
  250. if (Str_L[0] != 1) ABACUSerror("Str_L[0] != 1 in ABACUS_Scan.h Scan_for_Possible_Bases.");
  251. Nrapidities[0] = Mdown_remaining;
  252. // Set label:
  253. std::stringstream M0out;
  254. M0out << Nrapidities[0];
  255. possible_base_label[nfound] = M0out.str();
  256. for (int itype = 1; itype < Nrapidities.size(); ++itype)
  257. if (Nrapidities[itype] > 0) {
  258. possible_base_label[nfound] += TYPESEP;
  259. std::stringstream typeout;
  260. typeout << itype;
  261. possible_base_label[nfound] += typeout.str();
  262. possible_base_label[nfound] += EXCSEP;
  263. std::stringstream Mout;
  264. Mout << Nrapidities[itype];
  265. possible_base_label[nfound] += Mout.str();
  266. }
  267. nfound++;
  268. }
  269. else {
  270. // Preserve the number of strings at this level as compared to SeedState:
  271. Nrapidities[base_level_to_scan] = SeedNrap[base_level_to_scan];
  272. if (Mdown_remaining - Str_L[base_level_to_scan] * Nrapidities[base_level_to_scan] >= 0)
  273. Scan_for_Possible_Bases (SeedNrap, Str_L, Mdown_remaining - Str_L[base_level_to_scan] * Nrapidities[base_level_to_scan],
  274. possible_base_label, nfound, nexc_max_used, base_level_to_scan - 1, Nrapidities);
  275. // Reduce number of strings at this level as compared to SeedState:
  276. for (int i = 1; i <= ABACUS::min(SeedNrap[base_level_to_scan], nexc_max_used/Str_L[base_level_to_scan]); ++i) {
  277. Nrapidities[base_level_to_scan] = SeedNrap[base_level_to_scan] - i;
  278. if (Mdown_remaining - Str_L[base_level_to_scan] * Nrapidities[base_level_to_scan] >= 0)
  279. Scan_for_Possible_Bases (SeedNrap, Str_L, Mdown_remaining - Str_L[base_level_to_scan] * Nrapidities[base_level_to_scan],
  280. possible_base_label, nfound, nexc_max_used - i*Str_L[base_level_to_scan], base_level_to_scan - 1, Nrapidities);
  281. }
  282. // Increase the number of strings at this level as compared to SeedState:
  283. for (int i = 1; i <= ABACUS::min(Mdown_remaining/Str_L[base_level_to_scan], nexc_max_used/Str_L[base_level_to_scan]); ++i) {
  284. Nrapidities[base_level_to_scan] = SeedNrap[base_level_to_scan] + i;
  285. if (Mdown_remaining - Str_L[base_level_to_scan] * Nrapidities[base_level_to_scan] >= 0)
  286. Scan_for_Possible_Bases (SeedNrap, Str_L, Mdown_remaining - Str_L[base_level_to_scan] * Nrapidities[base_level_to_scan],
  287. possible_base_label, nfound, nexc_max_used - i*Str_L[base_level_to_scan], base_level_to_scan - 1, Nrapidities);
  288. }
  289. }
  290. }
  291. inline Vect<std::string> Possible_Bases (const Vect<int> SeedNrap, const Vect<int> Str_L, int Mdown)//const Heis_Bethe_State& SeedState)
  292. {
  293. int nexc_max_used = NEXC_MAX_HEIS;
  294. Vect<std::string> possible_base_label (1000);
  295. int nfound = 0;
  296. Vect<int> Nrapidities = SeedNrap;
  297. int Mdown_remaining = Mdown;
  298. Scan_for_Possible_Bases (SeedNrap, Str_L, Mdown_remaining, possible_base_label, nfound, nexc_max_used, SeedNrap.size() - 1, Nrapidities);
  299. // Copy results into a clean vector:
  300. Vect<std::string> possible_base_label_found (nfound);
  301. for (int i = 0; i < nfound; ++i) possible_base_label_found[i] = possible_base_label[i];
  302. return(possible_base_label_found);
  303. }
  304. //****************************************************************************
  305. template<class Tstate>
  306. class Scan_State_List {
  307. public:
  308. int ndef;
  309. Vect<Tstate> State;
  310. Vect<std::string> base_label;
  311. Vect<Scan_Info> info; // info for base and type of State[n]
  312. Vect<bool> flag_for_scan; // set to true, next round of scanning will use this base/type
  313. Vect<bool> scan_attempted; // whether this has already been attempted
  314. public:
  315. inline Scan_State_List (char whichDSF, const Tstate& SeedScanState);
  316. public:
  317. inline Tstate& Return_State (std::string base_label_ref); // returns a state corresponding to same base and type
  318. inline void Populate_List (char whichDSF, const Tstate& SeedScanState); // creates all types of states containing up to nexc_max excitations
  319. inline void Include_Info (Scan_Info& info_to_add, std::string base_label_ref);
  320. inline void Raise_Scanning_Flags (DP threshold); // checks whether base/type should be scanned based on simpler base/type combinations
  321. inline void Order_in_SRC ();
  322. inline void Save_Info (const char* sumfile_Cstr);
  323. inline void Load_Info (const char* sumfile_Cstr);
  324. };
  325. // Do the explicit class specializations:
  326. template<>
  327. inline Scan_State_List<LiebLin_Bethe_State>::Scan_State_List (char whichDSF, const LiebLin_Bethe_State& SeedScanState)
  328. : ndef(0), State(Vect<LiebLin_Bethe_State>(MAX_STATE_LIST_SIZE)), base_label(Vect<std::string>(MAX_STATE_LIST_SIZE)),
  329. info(Vect<Scan_Info>(MAX_STATE_LIST_SIZE)), flag_for_scan(Vect<bool>(false, MAX_STATE_LIST_SIZE)),
  330. scan_attempted(Vect<bool>(false, MAX_STATE_LIST_SIZE))
  331. {
  332. State[0] = SeedScanState;
  333. }
  334. template<>
  335. inline Scan_State_List<XXZ_Bethe_State>::Scan_State_List (char whichDSF, const XXZ_Bethe_State& SeedScanState)
  336. : ndef(0), State(Vect<XXZ_Bethe_State>(MAX_STATE_LIST_SIZE)), base_label(Vect<std::string>(MAX_STATE_LIST_SIZE)),
  337. info(Vect<Scan_Info>(MAX_STATE_LIST_SIZE)), flag_for_scan(Vect<bool>(false, MAX_STATE_LIST_SIZE)),
  338. scan_attempted(Vect<bool>(false, MAX_STATE_LIST_SIZE))
  339. {
  340. State[0] = SeedScanState;
  341. }
  342. template<>
  343. inline Scan_State_List<XXX_Bethe_State>::Scan_State_List (char whichDSF, const XXX_Bethe_State& SeedScanState)
  344. : ndef(0), State(Vect<XXX_Bethe_State>(MAX_STATE_LIST_SIZE)), base_label(Vect<std::string>(MAX_STATE_LIST_SIZE)),
  345. info(Vect<Scan_Info>(MAX_STATE_LIST_SIZE)), flag_for_scan(Vect<bool>(false, MAX_STATE_LIST_SIZE)),
  346. scan_attempted(Vect<bool>(false, MAX_STATE_LIST_SIZE))
  347. {
  348. State[0] = SeedScanState;
  349. }
  350. template<>
  351. inline Scan_State_List<XXZ_gpd_Bethe_State>::Scan_State_List (char whichDSF, const XXZ_gpd_Bethe_State& SeedScanState)
  352. : ndef(0), State(Vect<XXZ_gpd_Bethe_State>(MAX_STATE_LIST_SIZE)), base_label(Vect<std::string>(MAX_STATE_LIST_SIZE)),
  353. info(Vect<Scan_Info>(MAX_STATE_LIST_SIZE)), flag_for_scan(Vect<bool>(false, MAX_STATE_LIST_SIZE)),
  354. scan_attempted(Vect<bool>(false, MAX_STATE_LIST_SIZE))
  355. {
  356. State[0] = SeedScanState;
  357. }
  358. /* IN_DEVELOPMENT
  359. template<>
  360. inline Scan_State_List<ODSLF_XXZ_Bethe_State>::Scan_State_List (char whichDSF, const ODSLF_XXZ_Bethe_State& RefState)
  361. : ndef(0), State(Vect<ODSLF_XXZ_Bethe_State>(MAX_STATE_LIST_SIZE)), base_label(Vect<std::string>(MAX_STATE_LIST_SIZE)),
  362. info(Vect<Scan_Info>(MAX_STATE_LIST_SIZE)), flag_for_scan(Vect<bool>(false, MAX_STATE_LIST_SIZE)),
  363. scan_attempted(Vect<bool>(false, MAX_STATE_LIST_SIZE))
  364. {
  365. if (whichDSF == 'Z' || whichDSF == 'z') State[0] = ODSLF_XXZ_Bethe_State(RefState.chain, RefState.base.Mdown);
  366. else if (whichDSF == 'm') State[0] = ODSLF_XXZ_Bethe_State(RefState.chain, RefState.base.Mdown - 1);
  367. else if (whichDSF == 'p') State[0] = ODSLF_XXZ_Bethe_State(RefState.chain, RefState.base.Mdown + 1);
  368. else ABACUSerror("Unknown whichDSF in Scan_State_List<ODSLF_XXZ... constructor.");
  369. }
  370. */
  371. template<>
  372. inline LiebLin_Bethe_State& Scan_State_List<LiebLin_Bethe_State>::Return_State (std::string base_label_ref)
  373. {
  374. int n = 0;
  375. while (n < ndef && base_label_ref.compare(base_label[n]) != 0) n++;
  376. if (n == ndef) {
  377. State[n] = State[0];
  378. base_label[n] = base_label_ref;
  379. info[n].Nfull = 1LL; // Nfull not definable for LiebLin
  380. ndef++;
  381. }
  382. return(State[n]);
  383. }
  384. template<>
  385. inline XXZ_Bethe_State& Scan_State_List<XXZ_Bethe_State>::Return_State (std::string base_label_ref)
  386. {
  387. int n = 0;
  388. while (n < ndef && base_label_ref.compare(base_label[n]) != 0) n++;
  389. if (n == ndef) {
  390. Heis_Base checkbase (State[0].chain, base_label_ref);
  391. State[n] = XXZ_Bethe_State (State[0].chain, checkbase);
  392. info[n].Nfull = checkbase.dimH;
  393. ndef++;
  394. }
  395. return(State[n]);
  396. }
  397. template<>
  398. inline XXX_Bethe_State& Scan_State_List<XXX_Bethe_State>::Return_State (std::string base_label_ref)
  399. {
  400. int n = 0;
  401. while (n < ndef && base_label_ref.compare(base_label[n]) != 0) n++;
  402. if (n == ndef) {
  403. Heis_Base checkbase (State[0].chain, base_label_ref);
  404. State[n] = XXX_Bethe_State (State[0].chain, checkbase);
  405. info[n].Nfull = checkbase.dimH;
  406. ndef++;
  407. }
  408. return(State[n]);
  409. }
  410. template<>
  411. inline XXZ_gpd_Bethe_State& Scan_State_List<XXZ_gpd_Bethe_State>::Return_State (std::string base_label_ref)
  412. {
  413. int n = 0;
  414. while (n < ndef && base_label_ref.compare(base_label[n]) != 0) n++;
  415. if (n == ndef) {
  416. Heis_Base checkbase (State[0].chain, base_label_ref);
  417. State[n] = XXZ_gpd_Bethe_State (State[0].chain, checkbase);
  418. info[n].Nfull = checkbase.dimH;
  419. ndef++;
  420. }
  421. return(State[n]);
  422. }
  423. /* IN DEVELOPMENT
  424. template<>
  425. inline ODSLF_XXZ_Bethe_State& Scan_State_List<ODSLF_XXZ_Bethe_State>::Return_State (long long int base_id_ref, long long int type_id_ref)
  426. {
  427. int n = 0;
  428. while (n < ndef && !(base_id_ref == State[n].base_id && type_id_ref == State[n].type_id)) n++;
  429. if (n == ndef) {
  430. State[n] = ODSLF_XXZ_Bethe_State (State[0].chain, base_id_ref, type_id_ref);
  431. info[n].Nfull = State[n].maxid + 1LL;
  432. ndef++;
  433. }
  434. return(State[n]);
  435. }
  436. */
  437. template<>
  438. inline void Scan_State_List<LiebLin_Bethe_State>::Populate_List (char whichDSF, const LiebLin_Bethe_State& SeedScanState)
  439. {
  440. // For LiebLin_Bethe_State: only one base is used, so there is only one state here.
  441. if (ndef != 0) ABACUSerror("Please only populate a virgin Scan_State_List.");
  442. std::stringstream baselabel;
  443. baselabel << State[0].N;
  444. base_label[0] = baselabel.str();
  445. std::stringstream label0;
  446. label0 << State[0].N << LABELSEP << ABACUScoding[0] << LABELSEP;//"_0_";
  447. State[0].Set_to_Label (label0.str(), SeedScanState.Ix2);
  448. info[0].Nfull = 1LL; // Nfull not definable for LiebLin
  449. ndef = 1;
  450. }
  451. template<>
  452. inline void Scan_State_List<XXZ_Bethe_State>::Populate_List (char whichDSF, const XXZ_Bethe_State& SeedScanState)
  453. {
  454. // creates all types of states containing up to nexc_max excitations
  455. if (ndef != 0) ABACUSerror("Please only populate a virgin Scan_State_List.");
  456. // We assume that SeedScanState has quantum numbers which are set according to the relevant AveragingState.
  457. // This function creates a list of states with other bases in the vicinity of that of SeedScanState,
  458. // matching the quantum numbers as closely as possible.
  459. Vect<int> Str_L(SeedScanState.chain.Nstrings);
  460. for (int i = 0; i < SeedScanState.chain.Nstrings; ++i) Str_L[i] = SeedScanState.chain.Str_L[i];
  461. // First of all, we create a list of the possible bases themselves.
  462. Vect<std::string> bases_label = Possible_Bases (SeedScanState.base.Nrap, Str_L, SeedScanState.base.Mdown); // returns a vector of possible bases
  463. for (int ib = 0; ib < bases_label.size(); ++ib) {
  464. Heis_Base checkbase (State[0].chain, bases_label[ib]);
  465. State[ndef] = XXZ_Bethe_State (State[0].chain, checkbase);
  466. State[ndef].Set_to_Closest_Matching_Ix2_fixed_Base (SeedScanState);
  467. State[ndef].Set_Label_from_Ix2 (State[ndef].Ix2); // sets to trivial label for this base
  468. base_label[ndef] = bases_label[ib];
  469. info[ndef].Nfull = State[ndef].base.dimH;
  470. ndef++;
  471. if (ndef >= MAX_STATE_LIST_SIZE) ABACUSerror("Increase number of elements in ScanStateList.");
  472. }
  473. }
  474. template<>
  475. inline void Scan_State_List<XXX_Bethe_State>::Populate_List (char whichDSF, const XXX_Bethe_State& SeedScanState)
  476. {
  477. // creates all types of states containing up to nexc_max excitations
  478. if (ndef != 0) ABACUSerror("Please only populate a virgin Scan_State_List.");
  479. // We assume that SeedScanState has quantum numbers which are set according to the relevant AveragingState.
  480. // This function creates a list of states with other bases in the vicinity of that of SeedScanState,
  481. // matching the quantum numbers as closely as possible.
  482. Vect<int> Str_L(SeedScanState.chain.Nstrings);
  483. for (int i = 0; i < SeedScanState.chain.Nstrings; ++i) Str_L[i] = SeedScanState.chain.Str_L[i];
  484. // To take infinite rapidities into account, we use intermediate states with up to 2 less finite rapidities (1 for Szz, 2 for Spm)
  485. int nrinfrapmax = 0;
  486. if (whichDSF == 'z') nrinfrapmax = 1;
  487. else if (whichDSF == 'p') nrinfrapmax = ABACUS::min(2, SeedScanState.base.Mdown);
  488. Vect<int> Nrapmod = SeedScanState.base.Nrap;
  489. for (int nrinfrap = 0; nrinfrap <= nrinfrapmax; ++nrinfrap) {
  490. Nrapmod[0] = SeedScanState.base.Nrap[0] - nrinfrap;
  491. if (Nrapmod[0] < 0) ABACUSerror("Putting too many rapidities at infinity in ABACUS_Scan.h: Possible_Bases.");
  492. Vect<std::string> bases_label = Possible_Bases (Nrapmod, Str_L, SeedScanState.base.Mdown-nrinfrap); // returns a vector of possible bases
  493. for (int ib = 0; ib < bases_label.size(); ++ib) {
  494. Heis_Base checkbase (State[0].chain, bases_label[ib]);
  495. State[ndef] = XXX_Bethe_State (State[0].chain, checkbase);
  496. State[ndef].Set_to_Closest_Matching_Ix2_fixed_Base (SeedScanState);
  497. base_label[ndef] = bases_label[ib];
  498. info[ndef].Nfull = State[ndef].base.dimH;
  499. ndef++;
  500. if (ndef >= MAX_STATE_LIST_SIZE) ABACUSerror("Increase number of elements in ScanStateList.");
  501. }
  502. } // for nrinfrap
  503. }
  504. template<>
  505. inline void Scan_State_List<XXZ_gpd_Bethe_State>::Populate_List (char whichDSF, const XXZ_gpd_Bethe_State& SeedScanState)
  506. {
  507. // creates all types of states containing up to nexc_max excitations
  508. if (ndef != 0) ABACUSerror("Please only populate a virgin Scan_State_List.");
  509. // We assume that SeedScanState has quantum numbers which are set according to the relevant AveragingState.
  510. // This function creates a list of states with other bases in the vicinity of that of SeedScanState,
  511. // matching the quantum numbers as closely as possible.
  512. Vect<int> Str_L(SeedScanState.chain.Nstrings);
  513. for (int i = 0; i < SeedScanState.chain.Nstrings; ++i) Str_L[i] = SeedScanState.chain.Str_L[i];
  514. // First of all, we create a list of the possible bases themselves.
  515. Vect<std::string> bases_label = Possible_Bases (SeedScanState.base.Nrap, Str_L, SeedScanState.base.Mdown); // returns a vector of possible bases
  516. for (int ib = 0; ib < bases_label.size(); ++ib) {
  517. Heis_Base checkbase (State[0].chain, bases_label[ib]);
  518. State[ndef] = XXZ_gpd_Bethe_State (State[0].chain, checkbase);
  519. State[ndef].Set_to_Closest_Matching_Ix2_fixed_Base (SeedScanState);
  520. base_label[ndef] = bases_label[ib];
  521. info[ndef].Nfull = State[ndef].base.dimH;
  522. ndef++;
  523. if (ndef >= MAX_STATE_LIST_SIZE) ABACUSerror("Increase number of elements in ScanStateList.");
  524. }
  525. }
  526. /* IN DEVELOPMENT
  527. template<>
  528. inline void Scan_State_List<ODSLF_XXZ_Bethe_State>::Populate_List ()
  529. {
  530. // creates all types of states containing up to nexc_max excitations
  531. if (ndef != 0) ABACUSerror("Please only populate a virgin Scan_State_List.");
  532. //std::cout << "In Populate_List: " << State[0] << std::endl;
  533. Vect<long long int> bases_id = State[0].chain.Possible_Bases (State[0].base.Mdown); // returns a vector of possible bases
  534. //std::cout << "Mdown = " << State[0].base.Mdown << "\tPossible bases size: " << bases_id.size() << "\tPossible bases: " << bases_id << std::endl;
  535. for (int ib = 0; ib < bases_id.size(); ++ib) {
  536. ODSLF_Base checkbase (State[0].chain, bases_id[ib]);
  537. Vect<long long int> types_id = checkbase.Possible_Types (); // returns a vector of possible types
  538. //std::cout << "For base_id " << bases_id[ib] << "\t found types " << types_id << std::endl;
  539. for (int it = 0; it < types_id.size(); ++it) {
  540. if (bases_id[ib] < 1000000) { // FUDGE: consider only one-strings
  541. //std::cout << "Populate list: constructing state: " << bases_id[ib] << "\t" << types_id[it] << std::endl;
  542. State[ndef] = ODSLF_XXZ_Bethe_State (State[0].chain, bases_id[ib], types_id[it]);
  543. //std::cout << "Populate list: before setting id: " << std::endl << State[ndef] << std::endl;
  544. State[ndef].Set_to_id(0LL);
  545. //std::cout << "Populate list: after setting id: " << std::endl << State[ndef] << std::endl;
  546. info[ndef].Nfull = State[ndef].maxid + 1LL;
  547. ndef++;
  548. if (ndef >= MAX_STATE_LIST_SIZE) ABACUSerror("Increase number of elements in ScanStateList.");
  549. }
  550. }
  551. }
  552. }
  553. */
  554. template<class Tstate>
  555. inline void Scan_State_List<Tstate>::Include_Info (Scan_Info& info_to_add, std::string base_label_ref)
  556. {
  557. int n = 0;
  558. while (n < ndef && base_label_ref.compare(base_label[n]) != 0) n++;
  559. if (n == ndef) {
  560. std::cout << "ndef = " << ndef << std::endl;
  561. for (int i = 0; i < ndef; ++i) std::cout << base_label[i] << "\t";
  562. std::cout << std::endl;
  563. std::cout << "base_label_ref " << base_label_ref << std::endl;
  564. ABACUSerror("Did not find base_label_ref in Scan_State_List::Include_Info.");
  565. }
  566. info[n] += info_to_add;
  567. return;
  568. }
  569. template<class Tstate>
  570. inline void Scan_State_List<Tstate>::Raise_Scanning_Flags (DP threshold)
  571. {
  572. flag_for_scan = true;
  573. }
  574. template<class Tstate>
  575. inline void Scan_State_List<Tstate>::Order_in_SRC ()
  576. {
  577. if (ndef > 0) {
  578. Vect_INT index(ndef);
  579. for (int i = 0; i < ndef; ++i) index[i] = i;
  580. Vect<DP> sr (ndef);
  581. for (int i = 0; i < ndef; ++i) sr[i] = info[i].sumrule_obtained;
  582. sr.QuickSort(index, 0, ndef - 1);
  583. Vect<Tstate> State_ordered(ndef);
  584. Vect<std::string> base_label_ordered(ndef);
  585. Vect<Scan_Info> info_ordered(ndef);
  586. Vect<bool> flag_for_scan_ordered(ndef);
  587. // Put data in proper order
  588. for (int i = 0; i < ndef; ++i) {
  589. State_ordered[i] = State[index[ndef - 1 - i] ];
  590. base_label_ordered[i] = base_label[index[ndef - 1 - i] ];
  591. info_ordered[i] = info[index[ndef - 1 - i] ];
  592. flag_for_scan_ordered[i] = flag_for_scan[index[ndef - 1 - i] ];
  593. }
  594. // Put back in *this object:
  595. for (int i = 0; i < ndef; ++i) {
  596. State[i] = State_ordered[i];
  597. base_label[i] = base_label_ordered[i];
  598. info[i] = info_ordered[i];
  599. flag_for_scan[i] = flag_for_scan_ordered[i];
  600. } // The rest are all simply 0.
  601. }
  602. }
  603. template<class Tstate>
  604. inline void Scan_State_List<Tstate>::Save_Info (const char* sumfile_Cstr)
  605. {
  606. std::ofstream outfile;
  607. outfile.open(sumfile_Cstr);
  608. if (outfile.fail()) ABACUSerror("Could not open outfile... ");
  609. outfile.setf(std::ios::fixed);
  610. outfile.precision(16);
  611. outfile << std::setw(20) << "base" << std::setw(25) << "sumrule_obtained" << std::setw(25) << "Nfull" << std::setw(10) << "Ninadm" << std::setw(10) << "Ndata" << std::setw(10) << "conv" << std::setw(10) << "conv0" << std::setw(10) << "TT.";
  612. for (int i = 0; i < ndef; ++i)
  613. if (info[i].Nfull > 0.0) {
  614. int TT_hr = int(info[i].TT/3600);
  615. int TT_min = int((info[i].TT - 3600.0*TT_hr)/60);
  616. outfile << std::endl << std::setw(20) << base_label[i] << std::setw(25);
  617. if (info[i].sumrule_obtained < 1.0) outfile << std::fixed;
  618. else outfile << std::scientific;
  619. outfile << std::setprecision(16) << info[i].sumrule_obtained;
  620. if (info[i].Nfull < 1.0e+10) outfile << std::setw(25) << std::fixed << std::setprecision(0) << info[i].Nfull;
  621. else outfile << std::setw(25) << std::scientific << std::setprecision(16) << info[i].Nfull;
  622. outfile << std::setw(10) << info[i].Ninadm << std::setw(10) << info[i].Ndata << std::setw(10) << info[i].Ndata_conv << std::setw(10) << info[i].Ndata_conv0 << std::setw(10) << TT_hr << " h " << TT_min << " m " << std::fixed << std::showpoint << std::setprecision(3) << info[i].TT - 3600.0*TT_hr - 60.0*TT_min << " s";
  623. }
  624. outfile.close();
  625. }
  626. template<class Tstate>
  627. inline void Scan_State_List<Tstate>::Load_Info (const char* sumfile_Cstr)
  628. {
  629. std::ifstream infile;
  630. infile.open(sumfile_Cstr);
  631. if(infile.fail()) {
  632. std::cout << std::endl << sumfile_Cstr << std::endl;
  633. ABACUSerror("Could not open input file in Scan_State_List::Load_Info.");
  634. }
  635. // Load first line, containing informative text:
  636. char junk[256];
  637. infile.getline(junk, 256);
  638. // Now load the previous info's:
  639. std::string base_label_ref;
  640. DP sr_ref;
  641. DP Nfull_ref;
  642. long long int Ninadm_ref, Ndata_ref, conv_ref, conv0_ref;
  643. DP TT_ref;
  644. int TT_hr, TT_min;
  645. DP TT_sec;
  646. char a;
  647. while (infile.peek() != EOF) {
  648. infile >> base_label_ref >> sr_ref >> Nfull_ref >> Ninadm_ref >> Ndata_ref >> conv_ref >> conv0_ref >> TT_hr >> a >> TT_min >> a >> TT_sec >> a;
  649. TT_ref = 3600.0 * TT_hr + 60.0* TT_min + TT_sec;
  650. Scan_Info info_ref (sr_ref, Nfull_ref, Ninadm_ref, Ndata_ref, conv_ref, conv0_ref, TT_ref);
  651. (*this).Include_Info (info_ref, base_label_ref);
  652. }
  653. infile.close();
  654. return;
  655. }
  656. } // namespace ABACUS
  657. #endif