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.

General_Scan.cc 47KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229
  1. /**********************************************************
  2. This software is part of J.-S. Caux's ABACUS library.
  3. Copyright (c) J.-S. Caux.
  4. -----------------------------------------------------------
  5. File: src/SCAN/General_Scan.cc
  6. Purpose: universal implementation of state scanning:
  7. functions to descend down hierarchy of intermediate states.
  8. NOTE: since templated functions have to be in the same file,
  9. we put all scanning functions here. The externally-accessible
  10. functions are defined at the end of this file.
  11. ***********************************************************/
  12. #include <omp.h>
  13. #include "ABACUS.h"
  14. using namespace std;
  15. using namespace ABACUS;
  16. string LABEL_TO_CHECK = "bla";
  17. //string LABEL_TO_CHECK = "6_0_";
  18. //string LABEL_TO_CHECK = "6_2_22y32";
  19. namespace ABACUS {
  20. // Types of descendents:
  21. // 14 == iK stepped up, leading exc one step further (can lead to ph recombination)
  22. // 13 == iK up, next exc (nr of ph preserved, not taking possible ph recombination into account)
  23. // 12 == iK up, next ext (nr ph increased, not taking possible ph recombination into account)
  24. // 11 == iK down, leading exc one step further (can lead to ph recombination)
  25. // 10 == iK down, next exc (nr of ph preserved, not taking possible ph recombination into account)
  26. // 9 == iK down, next exc (nr ph increased, not taking possible ph recombination into account)
  27. // 8 == iK preserved, 14 and 11 (up to 2 ph recombinations)
  28. // 7 == iK preserved, 14 and 10
  29. // 6 == iK preserved, 14 and 9
  30. // 5 == iK preserved, 13 and 11
  31. // 4 == iK preserved, 13 and 10
  32. // 3 == iK preserved, 13 and 9
  33. // 2 == iK preserved, 12 and 11
  34. // 1 == iK preserved, 12 and 10
  35. // 0 == iK preserved, 12 and 9
  36. // For scanning over symmetric states, the interpretation is slightly different.
  37. // Types 14, 13 and 12 step iK up using the Ix2 on the right only,
  38. // and mirrors the change on the left Ix2.
  39. // Types 11, 10 and 9 step iK down using the Ix2 on the right only,
  40. // and mirrors the change on the left Ix2.
  41. // There is then no need of scanning over types 0 - 8.
  42. // By convention, types 9, 10 and 11 can call types 9 - 14; types 12-14 can only call types 12-14.
  43. bool Expect_ph_Recombination_iK_Up (string ScanIx2_label, const Vect<Vect<int> >& OriginIx2,
  44. const Vect<Vect<int> >& BaseScanIx2)
  45. {
  46. // This function returns true if descending further can lead to a particle-hole recombination.
  47. // The criteria which are used are:
  48. // - the active excitation has moved at least one step (so it has already created its p-h pair)
  49. // - there exists an OriginIx2 between the active Ix2 and the next Ix2
  50. // (to right or left depending on type of descendent)
  51. Vect<Vect<int> > ScanIx2 = Return_Ix2_from_Label (ScanIx2_label, OriginIx2);
  52. // Determine the level and index of the bottom-most left-most right-moving quantum number sits:
  53. int exclevel = -1;
  54. int excindex = 0;
  55. bool excfound = false;
  56. do {
  57. exclevel++;
  58. if (exclevel == ScanIx2.size()) { // there is no right-moving quantum number in ScanIx2
  59. break;
  60. }
  61. for (int alpha = 0; alpha < ScanIx2[exclevel].size(); ++alpha)
  62. if (ScanIx2[exclevel][alpha] > BaseScanIx2[exclevel][alpha]) {
  63. excindex = alpha;
  64. excfound = true;
  65. break;
  66. }
  67. } while (!excfound);
  68. // If we haven't found an excitation, then exclevel == ScanIx2.size() and excindex = 0;
  69. if (excfound && !BaseScanIx2[exclevel].includes(ScanIx2[exclevel][excindex])) {
  70. // there exists an already dispersing excitation which isn't in Origin
  71. // Is there a possible recombination?
  72. if (excindex < ScanIx2[exclevel].size() - 1) {
  73. // a particle to the right of excitation has already move right, so there is a hole
  74. // check that there exists an occupied Ix2 in Origin sitting between the excitation
  75. // and the next Ix2 to its right in ScanIx2
  76. for (int alpha = BaseScanIx2[exclevel].size() - 1; alpha >= 0; --alpha)
  77. if (BaseScanIx2[exclevel][alpha] > ScanIx2[exclevel][excindex]
  78. && BaseScanIx2[exclevel][alpha] < ScanIx2[exclevel][excindex + 1]) {
  79. return(true);
  80. }
  81. }
  82. } // if (excfound)
  83. return(false);
  84. }
  85. // Specialization for Lieb-Liniger:
  86. bool Expect_ph_Recombination_iK_Up (string ScanIx2_label, const LiebLin_Bethe_State& OriginState)
  87. {
  88. Vect<Vect<int> > OriginIx2here(1);
  89. OriginIx2here[0] = OriginState.Ix2;
  90. Vect<Vect<int> > BaseScanIx2here(1);
  91. BaseScanIx2here[0] = OriginState.Ix2;
  92. return(Expect_ph_Recombination_iK_Up (ScanIx2_label, OriginIx2here, BaseScanIx2here));
  93. }
  94. // Specialization for Heis
  95. bool Expect_ph_Recombination_iK_Up (string ScanIx2_label, const Heis_Bethe_State& OriginState)
  96. {
  97. return(Expect_ph_Recombination_iK_Up (ScanIx2_label, OriginState.Ix2, OriginState.Ix2));
  98. }
  99. bool Expect_ph_Recombination_iK_Down (string ScanIx2_label, const Vect<Vect<int> >& OriginIx2,
  100. const Vect<Vect<int> >& BaseScanIx2)
  101. {
  102. // This function returns true if descending further can lead to a particle-hole recombination.
  103. // The criteria which are used are:
  104. // - the active excitation has moved at least one step (so it has already created its p-h pair)
  105. // - there exists an OriginIx2 between the active Ix2 and the next Ix2
  106. // (to right or left depending on type of descendent)
  107. Vect<Vect<int> > ScanIx2 = Return_Ix2_from_Label (ScanIx2_label, OriginIx2);
  108. // Determine the level and index of the bottom-most right-most left-moving quantum number sits:
  109. int exclevel = -1;
  110. int excindex = 0;
  111. bool excfound = false;
  112. do {
  113. exclevel++;
  114. if (exclevel == ScanIx2.size()) { // there isn't a single left-moving quantum number in ScanIx2
  115. break;
  116. }
  117. for (int alpha = ScanIx2[exclevel].size() - 1; alpha >= 0; --alpha) {
  118. if (ScanIx2[exclevel][alpha] < BaseScanIx2[exclevel][alpha]) {
  119. excindex = alpha;
  120. excfound = true;
  121. break;
  122. }
  123. }
  124. } while (!excfound);
  125. // If we haven't found an excitation, then exclevel == ScanIx2.size() and excindex = 0;
  126. if (!excfound) excindex = ScanIx2[exclevel].size() - 1;
  127. if (excfound && !BaseScanIx2[exclevel].includes(ScanIx2[exclevel][excindex])) {
  128. // there exists an already dispersing excitation which isn't in Origin
  129. // Is there a possible recombination?
  130. if (excindex > 0) {
  131. // a particle to the left of excitation has already moved left, so there is a hole
  132. // check that there exists an occupied Ix2 in Origin sitting between the excitation
  133. // and the next Ix2 to its left in ScanIx2
  134. for (int alpha = 0; alpha < BaseScanIx2[exclevel].size(); ++alpha)
  135. if (BaseScanIx2[exclevel][alpha] > ScanIx2[exclevel][excindex - 1]
  136. && BaseScanIx2[exclevel][alpha] < ScanIx2[exclevel][excindex]) {
  137. return(true);
  138. }
  139. }
  140. } // if (excfound)
  141. return(false);
  142. }
  143. // Specialization for Lieb-Liniger:
  144. bool Expect_ph_Recombination_iK_Down (string ScanIx2_label, const LiebLin_Bethe_State& OriginState)
  145. {
  146. Vect<Vect<int> > OriginIx2here(1);
  147. OriginIx2here[0] = OriginState.Ix2;
  148. Vect<Vect<int> > BaseScanIx2here(1);
  149. BaseScanIx2here[0] = OriginState.Ix2;
  150. return(Expect_ph_Recombination_iK_Down (ScanIx2_label, OriginIx2here, BaseScanIx2here));
  151. }
  152. // Specialization for Heis
  153. bool Expect_ph_Recombination_iK_Down (string ScanIx2_label, const Heis_Bethe_State& OriginState)
  154. {
  155. return(Expect_ph_Recombination_iK_Down (ScanIx2_label, OriginState.Ix2, OriginState.Ix2));
  156. }
  157. template<class Tstate>
  158. Scan_Info General_Scan (char whichDSF, int iKmin, int iKmax, int iKmod, DP kBT,
  159. Tstate& AveragingState, Tstate& SeedScanState, string defaultScanStatename,
  160. int Max_Secs, DP target_sumrule, bool refine,
  161. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  162. {
  163. // Performs the scan over excited states, writing data to file.
  164. // AveragingState is the state on which the correlations are calculated.
  165. // SeedScanState is the originator of all scan states.
  166. // This distinction is kept to allow for quenches and finite temperatures.
  167. // This function is also called by the parallel implementation of ABACUS.
  168. // In this case, file names carry a rank and nr_processors suffix.
  169. // In fact, the parallelization can be done in incremental levels.
  170. // If paralevel == 0, the run is serial.
  171. // If paralevel == n, the run is parallelized in a tree with n levels of branching.
  172. // A paralevel == 1 branching's files have a suffix of the form "_3_8", meaning that this
  173. // is the rank 3 out of 8 processors.
  174. // A paralevel == 2 branching's files have a suffix of the form "_3_8_2_8", meaning that this
  175. // is the rank 2 out of 8 subscan of the _3_8 scan.
  176. bool in_parallel = (paralevel > 0);
  177. if (in_parallel && (rank.size() != paralevel || nr_processors.size() != paralevel)) {
  178. cout << "paralevel = " << paralevel << "\trank.size() = " << rank.size()
  179. << "\tnr_processors.size() = " << nr_processors.size() << endl;
  180. cout << "rank = " << rank << endl;
  181. cout << "nr_processors = " << nr_processors << endl;
  182. ABACUSerror("Inconsistent paralevel, rank or nr_processors in General_Scan.");
  183. }
  184. if (in_parallel && !refine) ABACUSerror("Must refine when using parallel ABACUS");
  185. // expected cost on data_value of adding a particle-hole excitation.
  186. DP ph_cost = Particle_Hole_Excitation_Cost (whichDSF, AveragingState);
  187. int Max_Secs_used = int(0.9 * Max_Secs); // we don't start any new ithread loop beyond this point
  188. int Max_Secs_alert = int(0.95 * Max_Secs); // we break any ongoing ithread loop beyond this point
  189. stringstream filenameprefix;
  190. Data_File_Name (filenameprefix, whichDSF, iKmin, iKmax, kBT,
  191. AveragingState, SeedScanState, defaultScanStatename);
  192. if (in_parallel)
  193. for (int r = 0; r < paralevel; ++r)
  194. filenameprefix << "_" << rank[r] << "_" << nr_processors[r];
  195. string prefix = filenameprefix.str();
  196. stringstream filenameprefix_prevparalevel;
  197. // without the rank and nr_processors of the highest paralevel
  198. Data_File_Name (filenameprefix_prevparalevel, whichDSF, iKmin, iKmax, kBT,
  199. AveragingState, SeedScanState, defaultScanStatename);
  200. if (in_parallel) for (int r = 0; r < paralevel - 1; ++r)
  201. filenameprefix << "_" << rank[r] << "_" << nr_processors[r];
  202. string prefix_prevparalevel = filenameprefix_prevparalevel.str();
  203. stringstream RAW_stringstream; string RAW_string;
  204. stringstream INADM_stringstream; string INADM_string;
  205. stringstream CONV0_stringstream; string CONV0_string;
  206. stringstream STAT_stringstream; string STAT_string;
  207. stringstream LOG_stringstream; string LOG_string;
  208. stringstream THR_stringstream; string THR_string;
  209. stringstream THRDIR_stringstream; string THRDIR_string;
  210. stringstream SRC_stringstream; string SRC_string;
  211. stringstream SUM_stringstream; string SUM_string;
  212. RAW_stringstream << prefix << ".raw";
  213. INADM_stringstream << prefix << ".inadm";
  214. CONV0_stringstream << prefix << ".conv0";
  215. STAT_stringstream << prefix << ".stat";
  216. LOG_stringstream << prefix << ".log";
  217. THR_stringstream << prefix << ".thr";
  218. THRDIR_stringstream << prefix << "_thrdir";
  219. SRC_stringstream << prefix << ".src";
  220. SUM_stringstream << prefix << ".sum";
  221. RAW_string = RAW_stringstream.str(); const char* RAW_Cstr = RAW_string.c_str();
  222. INADM_string = INADM_stringstream.str(); const char* INADM_Cstr = INADM_string.c_str();
  223. CONV0_string = CONV0_stringstream.str(); const char* CONV0_Cstr = CONV0_string.c_str();
  224. STAT_string = STAT_stringstream.str(); const char* STAT_Cstr = STAT_string.c_str();
  225. LOG_string = LOG_stringstream.str(); const char* LOG_Cstr = LOG_string.c_str();
  226. THR_string = THR_stringstream.str(); const char* THR_Cstr = THR_string.c_str();
  227. SRC_string = SRC_stringstream.str(); const char* SRC_Cstr = SRC_string.c_str();
  228. SUM_string = SUM_stringstream.str(); const char* SUM_Cstr = SUM_string.c_str();
  229. THRDIR_string = THRDIR_stringstream.str();
  230. fstream RAW_outfile;
  231. if (!refine || in_parallel) RAW_outfile.open(RAW_Cstr, fstream::out | fstream::trunc);
  232. else RAW_outfile.open(RAW_Cstr, fstream::out | fstream::app);
  233. if (RAW_outfile.fail()) {
  234. cout << RAW_Cstr << endl;
  235. ABACUSerror("Could not open RAW_outfile... ");
  236. }
  237. RAW_outfile.precision(16);
  238. fstream INADM_outfile;
  239. if (!refine || in_parallel) INADM_outfile.open(INADM_Cstr, fstream::out | fstream::trunc);
  240. else INADM_outfile.open(INADM_Cstr, fstream::out | fstream::app);
  241. if (INADM_outfile.fail()) ABACUSerror("Could not open INADM_outfile... ");
  242. INADM_outfile.precision(16);
  243. fstream CONV0_outfile;
  244. if (!refine || in_parallel) CONV0_outfile.open(CONV0_Cstr, fstream::out | fstream::trunc);
  245. else CONV0_outfile.open(CONV0_Cstr, fstream::out | fstream::app);
  246. if (CONV0_outfile.fail()) ABACUSerror("Could not open CONV0_outfile... ");
  247. CONV0_outfile.precision(16);
  248. fstream STAT_outfile;
  249. if (!refine || in_parallel) STAT_outfile.open(STAT_Cstr, fstream::out | fstream::trunc);
  250. else STAT_outfile.open(STAT_Cstr, fstream::out | fstream::app);
  251. if (STAT_outfile.fail()) ABACUSerror("Could not open STAT_outfile... ");
  252. STAT_outfile.precision(8);
  253. ofstream LOG_outfile;
  254. if (!in_parallel) {
  255. if (!refine) LOG_outfile.open(LOG_Cstr, fstream::out | fstream::trunc);
  256. else LOG_outfile.open(LOG_Cstr, fstream::out | fstream::app);
  257. if (LOG_outfile.fail()) ABACUSerror("Could not open LOG_outfile... ");
  258. LOG_outfile.precision(16);
  259. }
  260. else { // in_parallel
  261. LOG_outfile.open(LOG_Cstr, fstream::out | fstream::trunc);
  262. if (LOG_outfile.fail()) ABACUSerror("Could not open LOG_outfile... ");
  263. LOG_outfile.precision(16);
  264. }
  265. Scan_Info scan_info;
  266. if (!refine) mkdir(THRDIR_string.c_str(), S_IRWXU | S_IRWXG | S_IRWXO);
  267. Scan_Thread_Data paused_thread_data (THRDIR_string, refine);
  268. if (refine) {
  269. paused_thread_data.Load();
  270. if (!in_parallel) scan_info.Load(SRC_Cstr);
  271. }
  272. Scan_Info scan_info_before = scan_info; // for LOG file
  273. Scan_Info scan_info_before_descent;
  274. Scan_Info scan_info_obtained_in_descent;
  275. Scan_State_List<Tstate> ScanStateList (whichDSF, SeedScanState);
  276. ScanStateList.Populate_List(whichDSF, SeedScanState);
  277. if (refine && !in_parallel) ScanStateList.Load_Info (SUM_Cstr);
  278. else if (in_parallel && rank.sum() == 0) {}; // do nothing, keep info in the higher .sum file!
  279. DP Chem_Pot = Chemical_Potential (AveragingState);
  280. DP sumrule_factor = Sumrule_Factor (whichDSF, AveragingState, Chem_Pot, iKmin, iKmax);
  281. // Now go for it !
  282. bool at_least_one_new_flag_raised = false;
  283. int Ndata_previous_cycle = 0;
  284. int ninadm = 0; // nr of inadm states for which we save info in .inadm file. Save first 1000.
  285. int nconv0 = 0; // nr of unconv states for which we save info in .conv0 file. Save first 1000.
  286. double start_time_omp = omp_get_wtime();
  287. double current_time_omp = omp_get_wtime();
  288. #pragma omp parallel
  289. do {
  290. int omp_thread_nr = omp_get_thread_num();
  291. if ((paused_thread_data.lowest_il_with_nthreads_neq_0 == paused_thread_data.nlists - 1)
  292. && omp_thread_nr > 0) {
  293. double start_time_wait = omp_get_wtime();
  294. double stop_time_wait;
  295. do {
  296. for (int i = 0; i < 100000; ++i) { }
  297. stop_time_wait = omp_get_wtime();
  298. } while (stop_time_wait - start_time_wait < 5.0);
  299. }
  300. double start_time_cycle_omp = omp_get_wtime();
  301. at_least_one_new_flag_raised = false;
  302. #pragma omp master
  303. {
  304. double start_time_flags = omp_get_wtime();
  305. // First flag the new base/type 's that we need to include:
  306. ScanStateList.Raise_Scanning_Flags
  307. (exp(-paused_thread_data.logscale * paused_thread_data.lowest_il_with_nthreads_neq_0));
  308. // Get these base/type started:
  309. for (int i = 0; i < ScanStateList.ndef; ++i) {
  310. if (ScanStateList.flag_for_scan[i] && ScanStateList.info[i].Ndata == 0
  311. && !ScanStateList.scan_attempted[i]) {
  312. Scan_Info scan_info_flags;
  313. at_least_one_new_flag_raised = true;
  314. ScanStateList.scan_attempted[i] = true;
  315. Tstate ScanState;
  316. ScanState = ScanStateList.State[i];
  317. DP data_value = -1.0;
  318. bool admissible = ScanState.Check_Admissibility(whichDSF);
  319. if (admissible) {
  320. ScanState.Compute_All(true);
  321. if (ScanState.conv) {
  322. // Put momentum in fundamental window, if possible:
  323. int iKexc = ScanState.iK - AveragingState.iK;
  324. while (iKexc > iKmax && iKexc - iKmod >= iKmin) iKexc -= iKmod;
  325. while (iKexc < iKmin && iKexc + iKmod <= iKmax) iKexc += iKmod;
  326. stringstream rawfile_entry;
  327. data_value = Compute_Matrix_Element_Contrib (whichDSF, iKmin, iKmax, ScanState, AveragingState, Chem_Pot, rawfile_entry);
  328. {
  329. #pragma omp critical
  330. RAW_outfile << rawfile_entry.str();
  331. }
  332. {
  333. #pragma omp critical
  334. if (iKexc >= iKmin && iKexc <= iKmax) {
  335. scan_info_flags.Ndata++;
  336. scan_info_flags.Ndata_conv++;
  337. scan_info_flags.sumrule_obtained += data_value*sumrule_factor;
  338. }
  339. }
  340. // If we force descent: modify data_value by hand so that descent is forced on next scanning pass
  341. for (int itype = 0; itype < 15; ++itype) {
  342. DP data_value_used = 0.1* exp(-paused_thread_data.logscale * ABACUS::min(0, paused_thread_data.lowest_il_with_nthreads_neq_0));
  343. if (Force_Descent(whichDSF, ScanState, AveragingState, itype, iKmod, Chem_Pot))
  344. data_value = data_value_used;
  345. }
  346. Vect<bool> allowed(false, 15);
  347. if (whichDSF == 'B') { // symmetric state scanning
  348. allowed[9] = true; allowed[10] = true; allowed[11] = true;
  349. allowed[12] = true; allowed[13] = true; allowed[14] = true;
  350. }
  351. else {
  352. allowed[0] = (iKexc >= iKmin && iKexc <= iKmax);
  353. allowed[1] = allowed[0]; allowed[2] = allowed[0];
  354. allowed[3] = allowed[0]; allowed[4] = allowed[0];
  355. allowed[5] = allowed[0]; allowed[6] = allowed[0];
  356. allowed[7] = allowed[0]; allowed[8] = allowed[0];
  357. allowed[9] = (iKexc > iKmin);
  358. allowed[10] = allowed[9]; allowed[11] = allowed[9];
  359. allowed[12] = (iKexc < iKmax);
  360. allowed[13] = allowed[12]; allowed[14] = allowed[12];
  361. }
  362. for (int type_required_here = 0; type_required_here < 15; ++type_required_here) {
  363. if (!allowed[type_required_here]) continue;
  364. // All cases here are such that the ScanState hasn't been descended yet,
  365. // so we simply use data_value as expected data value:
  366. {
  367. #pragma omp critical
  368. paused_thread_data.Include_Thread (abs(data_value), ScanState.label, type_required_here);
  369. }
  370. }
  371. } // if (ScanState.conv)
  372. else {
  373. if (nconv0++ < 1000)
  374. CONV0_outfile << setw(25) << ScanState.label << setw(25) << ScanState.diffsq
  375. << setw(5) << ScanState.Check_Rapidities()
  376. << setw(25) << ScanState.String_delta() << endl;
  377. scan_info_flags.Ndata++;
  378. scan_info_flags.Ndata_conv0++;
  379. }
  380. } // if admissible
  381. else { // if inadmissible, modify data_value by hand so that descent is forced on next scanning pass
  382. if (ninadm++ < 10000000) INADM_outfile << ScanState.label << endl;
  383. scan_info_flags.Ndata++;
  384. scan_info_flags.Ninadm++;
  385. // Put momentum in fundamental window, if possible:
  386. int iKexc = ScanState.iK - AveragingState.iK;
  387. while (iKexc > iKmax && iKexc - iKmod >= iKmin) iKexc -= iKmod;
  388. while (iKexc < iKmin && iKexc + iKmod <= iKmax) iKexc += iKmod;
  389. DP data_value = 1.0e-32;
  390. for (int itype = 0; itype < 15; ++itype)
  391. if (Force_Descent(whichDSF, ScanState, AveragingState, itype, iKmod, Chem_Pot))
  392. data_value = 0.1* exp(-paused_thread_data.logscale * paused_thread_data.lowest_il_with_nthreads_neq_0);
  393. Vect<bool> allowed(false, 15);
  394. if (whichDSF == 'B') {
  395. // We scan over symmetric states. Only types 14 down to 9 are allowed.
  396. allowed[9] = true; allowed[10] = true; allowed[11] = true;
  397. allowed[12] = true; allowed[13] = true; allowed[14] = true;
  398. }
  399. else {
  400. allowed[0] = (iKexc >= iKmin && iKexc <= iKmax);
  401. allowed[1] = allowed[0]; allowed[2] = allowed[0];
  402. allowed[3] = allowed[0]; allowed[4] = allowed[0];
  403. allowed[5] = allowed[0]; allowed[6] = allowed[0];
  404. allowed[7] = allowed[0]; allowed[8] = allowed[0];
  405. allowed[9] = (iKexc > iKmin);
  406. allowed[10] = allowed[9]; allowed[11] = allowed[9];
  407. allowed[12] = (iKexc < iKmax);
  408. allowed[13] = allowed[12]; allowed[14] = allowed[12];
  409. }
  410. for (int type_required_here = 0; type_required_here < 15; ++type_required_here) {
  411. if (!allowed[type_required_here]) continue;
  412. {
  413. #pragma omp critical
  414. paused_thread_data.Include_Thread (abs(data_value), ScanState.label, type_required_here);
  415. }
  416. }
  417. } // inadmissible
  418. scan_info_flags.TT += omp_get_wtime() - start_time_flags;
  419. // Put this info into the appropriate ScanStateList.info
  420. {
  421. #pragma omp critical
  422. ScanStateList.Include_Info(scan_info_flags, ScanStateList.base_label[i]);
  423. scan_info += scan_info_flags;
  424. }
  425. } // if flag_for_scan
  426. } // for i
  427. } // #pragma omp master
  428. // Now we deal with the previously existing paused threads:
  429. Vect<Scan_Thread> threads_to_do;
  430. int il_to_do = paused_thread_data.lowest_il_with_nthreads_neq_0; // for resaving threads in case we're out of time
  431. DP current_threshold = exp(-paused_thread_data.logscale * il_to_do);
  432. {
  433. #pragma omp critical
  434. threads_to_do = paused_thread_data.Extract_Next_Scan_Threads();
  435. }
  436. int ithread;
  437. {
  438. for (ithread = 0; ithread < threads_to_do.size(); ++ithread) {
  439. Scan_Info scan_info_this_ithread;
  440. double start_time_this_ithread = omp_get_wtime();
  441. // If we don't have time anymore, resave the threads instead of computing them:
  442. if (start_time_this_ithread - start_time_omp > Max_Secs_alert) {
  443. for (int ith = ithread; ith < threads_to_do.size(); ++ith) {
  444. #pragma omp critical
  445. paused_thread_data.Include_Thread (il_to_do, threads_to_do[ith].label, threads_to_do[ith].type);
  446. }
  447. break; // jump out of ithread loop
  448. }
  449. Tstate ScanState;
  450. {
  451. #pragma omp critical
  452. ScanState = ScanStateList.Return_State(Extract_Base_Label(threads_to_do[ithread].label));
  453. }
  454. Tstate BaseScanState; BaseScanState = ScanState;
  455. ScanState.Set_to_Label(threads_to_do[ithread].label, BaseScanState.Ix2);
  456. Tstate ScanStateBeingDescended = ScanState;
  457. // STARTING Descend_and_Compute block:
  458. int type_required = threads_to_do[ithread].type;
  459. ScanState.Compute_Momentum();
  460. Vect<string> desc_label;
  461. bool disperse_only_current_exc_up = false;
  462. if (type_required == 14 || type_required == 8 || type_required == 7 || type_required == 6)
  463. disperse_only_current_exc_up = true;
  464. bool preserve_nexc_up = false;
  465. if (type_required == 13 || type_required == 5 || type_required == 4 || type_required == 3)
  466. preserve_nexc_up = true;
  467. bool disperse_only_current_exc_down = false;
  468. if (type_required == 11 || type_required == 8 || type_required == 5 || type_required == 2)
  469. disperse_only_current_exc_down = true;
  470. bool preserve_nexc_down = false;
  471. if (type_required == 10 || type_required == 7 || type_required == 4 || type_required == 1)
  472. preserve_nexc_down = true;
  473. if (whichDSF == 'B') { // symmetric state scanning
  474. if (type_required >= 9 && type_required <= 11)
  475. desc_label = Descendent_States_with_iK_Stepped_Down_rightIx2only
  476. (ScanState.label, BaseScanState, disperse_only_current_exc_down, preserve_nexc_down);
  477. else if (type_required >= 12 && type_required <= 14)
  478. desc_label = Descendent_States_with_iK_Stepped_Up_rightIx2only
  479. (ScanState.label, BaseScanState, disperse_only_current_exc_up, preserve_nexc_up);
  480. }
  481. else {
  482. if (type_required >= 0 && type_required <= 8) {
  483. desc_label = Descendent_States_with_iK_Preserved
  484. (ScanState.label, BaseScanState, disperse_only_current_exc_up, preserve_nexc_up,
  485. disperse_only_current_exc_down, preserve_nexc_down);
  486. }
  487. else if (type_required >= 9 && type_required <= 11)
  488. desc_label = Descendent_States_with_iK_Stepped_Down
  489. (ScanState.label, BaseScanState, disperse_only_current_exc_down, preserve_nexc_down);
  490. else if (type_required >= 12 && type_required <= 14)
  491. desc_label = Descendent_States_with_iK_Stepped_Up
  492. (ScanState.label, BaseScanState, disperse_only_current_exc_up, preserve_nexc_up);
  493. }
  494. string label_here = ScanState.label;
  495. for (int idesc = 0; idesc < desc_label.size(); ++idesc) {
  496. ScanState.Set_to_Label (desc_label[idesc], BaseScanState.Ix2);
  497. bool admissible = ScanState.Check_Admissibility(whichDSF);
  498. DP data_value = 0.0;
  499. ScanState.conv = false;
  500. ScanState.Compute_Momentum(); // since momentum is used as forced descent criterion
  501. if (admissible) {
  502. ScanState.Compute_All (idesc == 0);
  503. if (ScanState.conv) {
  504. // Put momentum in fundamental window, if possible:
  505. int iKexc = ScanState.iK - AveragingState.iK;
  506. while (iKexc > iKmax && iKexc - iKmod >= iKmin) iKexc -= iKmod;
  507. while (iKexc < iKmin && iKexc + iKmod <= iKmax) iKexc += iKmod;
  508. stringstream rawfile_entry;
  509. data_value = Compute_Matrix_Element_Contrib (whichDSF, iKmin, iKmax, ScanState, AveragingState,
  510. Chem_Pot, rawfile_entry);
  511. {
  512. #pragma omp critical
  513. RAW_outfile << rawfile_entry.str();
  514. if (iKexc >= iKmin && iKexc <= iKmax) {
  515. scan_info_this_ithread.Ndata++;
  516. scan_info_this_ithread.Ndata_conv++;
  517. scan_info_this_ithread.sumrule_obtained += data_value*sumrule_factor;
  518. }
  519. }
  520. // Uncomment line below if .stat file is desired:
  521. // STAT_outfile << setw(20) << label_here << "\t" << setw(5) << type_required
  522. // << "\t" << setw(16) << std::scientific
  523. // << exp(-paused_thread_data.logscale * il_to_do)
  524. // << "\t" << setw(20) << ScanState.label << "\t" << setw(16) << data_value
  525. // << "\t" << setw(16) << std::fixed << setprecision(8)
  526. // << data_value/exp(-paused_thread_data.logscale * il_to_do) << endl;
  527. // Uncomment below if alerts for unexpectedly high data_value (as compared to threshold) are desired
  528. // if (fabs(data_value) > 10.0* current_threshold) {
  529. // cout << "\nAlert: data_value > 10* threshold, " << data_value << "\t" << current_threshold << endl;
  530. // cout << " for state " << ScanState.label << " descendent of type " << type_required
  531. // << " of state " << ScanStateBeingDescended.label << endl;
  532. // cout << AveragingState.Ix2 << endl;
  533. // cout << ScanStateBeingDescended.Ix2 << endl;
  534. // cout << ScanState.Ix2 << endl;
  535. // }
  536. } // if (ScanState.conv)
  537. else {
  538. if (nconv0++ < 1000)
  539. CONV0_outfile << setw(25) << ScanState.label << setw(25)
  540. << ScanState.diffsq << setw(5) << ScanState.Check_Rapidities()
  541. << setw(25) << ScanState.String_delta() << endl;
  542. scan_info_this_ithread.Ndata++;
  543. scan_info_this_ithread.Ndata_conv0++;
  544. }
  545. } // if (admissible)
  546. else {
  547. if (ninadm++ < 1000000) INADM_outfile << ScanState.label << endl;
  548. scan_info_this_ithread.Ndata++;
  549. scan_info_this_ithread.Ninadm++;
  550. }
  551. Tstate state_to_descend; state_to_descend = ScanState; // for checking
  552. ScanState.Compute_Momentum();
  553. // Put momentum in fundamental window, if possible:
  554. int iKexc = ScanState.iK - AveragingState.iK;
  555. while (iKexc > iKmax && iKexc - iKmod >= iKmin) iKexc -= iKmod;
  556. while (iKexc < iKmin && iKexc + iKmod <= iKmax) iKexc += iKmod;
  557. // Momentum-preserving are only descended to momentum-preserving.
  558. // Momentum-increasing are only descended to momentum-preserving and momentum-increasing.
  559. // Momentum-decreasing are only descended to momentum-preserving and momentum-decreasing.
  560. Vect<bool> allowed(false, 15);
  561. if (whichDSF == 'B') {
  562. // We scan over symmetric states. Only types 14 down to 9 are allowed.
  563. if (type_required >= 9 && type_required <= 11) { // iK stepped down on rightIx2; step further up or down
  564. allowed[9] = true; allowed[10] = true; allowed[11] = true;
  565. allowed[12] = true; allowed[13] = true; allowed[14] = true;
  566. }
  567. else if (type_required >= 12 && type_required <= 14) { // iK stepped up on rightIx2; only step further up
  568. allowed[12] = true; allowed[13] = true; allowed[14] = true;
  569. }
  570. }
  571. else {
  572. if (type_required >= 0 && type_required <= 8) { // momentum-preserving
  573. allowed[0] = (iKexc >= iKmin && iKexc <= iKmax);
  574. allowed[9] = false;
  575. allowed[12] = false;
  576. }
  577. if (type_required >= 9 && type_required <= 11) { // momentum-decreasing
  578. allowed[0] = (iKexc >= iKmin && iKexc <= iKmax);
  579. allowed[9] = (iKexc > iKmin);
  580. allowed[12] = false;
  581. }
  582. if (type_required >= 12 && type_required <= 14) { // momentum-increasing
  583. allowed[0] = (iKexc >= iKmin && iKexc <= iKmax);
  584. allowed[9] = false;
  585. allowed[12] = (iKexc < iKmax);
  586. }
  587. // The others are just copies of the ones above:
  588. allowed[1] = allowed[0]; allowed[2] = allowed[0]; allowed[3] = allowed[0]; allowed[4] = allowed[0]; allowed[5] = allowed[0]; allowed[6] = allowed[0]; allowed[7] = allowed[0]; allowed[8] = allowed[0];
  589. allowed[10] = allowed[9]; allowed[11] = allowed[9];
  590. allowed[13] = allowed[12]; allowed[14] = allowed[12];
  591. }
  592. for (int type_required_here = 0; type_required_here < 15; ++type_required_here) {
  593. if (!allowed[type_required_here]) continue;
  594. // Reset ScanState to what it was, if change on first pass
  595. if (type_required_here > 0) ScanState = state_to_descend;
  596. // We determine if we carry on scanning based on the data_value obtained, or forcing conditions:
  597. DP expected_abs_data_value = abs(data_value);
  598. //++G_7 logic
  599. if ((type_required_here == 14 || type_required_here == 8
  600. || type_required_here == 7 || type_required_here == 6)
  601. && Expect_ph_Recombination_iK_Up (ScanState.label, BaseScanState))
  602. expected_abs_data_value /= ph_cost;
  603. if (type_required_here == 12 || type_required_here == 2
  604. || type_required_here == 1 || type_required_here == 0)
  605. expected_abs_data_value *= ph_cost;
  606. if ((type_required_here == 11 || type_required_here == 8
  607. || type_required_here == 5 || type_required_here == 2)
  608. && Expect_ph_Recombination_iK_Down (ScanState.label, BaseScanState))
  609. expected_abs_data_value /= ph_cost;
  610. if (type_required_here == 9 || type_required_here == 6
  611. || type_required_here == 3 || type_required_here == 0)
  612. expected_abs_data_value *= ph_cost;
  613. {
  614. #pragma omp critical
  615. paused_thread_data.Include_Thread (expected_abs_data_value, ScanState.label, type_required_here);
  616. }
  617. } // for type_required_here
  618. } // for idesc
  619. // FINISHED Descend_and_Compute block
  620. scan_info_this_ithread.TT += omp_get_wtime() - start_time_this_ithread;
  621. #pragma omp critical
  622. {
  623. scan_info += scan_info_this_ithread;
  624. ScanStateList.Include_Info(scan_info_this_ithread, Extract_Base_Label(threads_to_do[ithread].label));
  625. }
  626. } // for ithread
  627. } // omp parallel region
  628. #pragma omp master
  629. {
  630. if (!in_parallel)
  631. LOG_outfile << "Master cycling. Ndata_conv " << scan_info.Ndata_conv
  632. << ". Threshold " << paused_thread_data.lowest_il_with_nthreads_neq_0 << " "
  633. << setw(9) << setprecision(3)
  634. << exp(-paused_thread_data.logscale * paused_thread_data.lowest_il_with_nthreads_neq_0)
  635. << ". " << setw(12) << scan_info.Ndata - Ndata_previous_cycle << " new data. Nr of threads: "
  636. << setw(14) << paused_thread_data.nthreads_total.sum()
  637. << ". Saturation: " << setprecision(12) << scan_info.sumrule_obtained << endl;
  638. Ndata_previous_cycle = scan_info.Ndata;
  639. }
  640. current_time_omp = omp_get_wtime();
  641. } while (current_time_omp - start_time_omp < Max_Secs_used
  642. && scan_info.sumrule_obtained < target_sumrule
  643. );
  644. // This closes the #pragram omp parallel block
  645. RAW_outfile.close();
  646. INADM_outfile.close();
  647. CONV0_outfile.close();
  648. STAT_outfile.close();
  649. scan_info.Save(SRC_Cstr);
  650. Scan_Info scan_info_refine = scan_info;
  651. scan_info_refine -= scan_info_before;
  652. if (!in_parallel) {
  653. if (scan_info.sumrule_obtained >= target_sumrule)
  654. LOG_outfile << endl << "Achieved sumrule saturation of " << scan_info.sumrule_obtained
  655. << "\t(target was " << target_sumrule << ")." << endl << endl;
  656. if (!refine) {
  657. LOG_outfile << "Main run info: " << scan_info << endl;
  658. LOG_outfile << "Latest threshold level " << paused_thread_data.lowest_il_with_nthreads_neq_0
  659. << " " << std::scientific << setprecision(3)
  660. << exp(-paused_thread_data.logscale * paused_thread_data.lowest_il_with_nthreads_neq_0) << endl;
  661. }
  662. else if (refine) {
  663. LOG_outfile << "Refining info: " << scan_info_refine << endl;
  664. LOG_outfile << "Latest threshold level " << paused_thread_data.lowest_il_with_nthreads_neq_0
  665. << " " << std::scientific << setprecision(3)
  666. << exp(-paused_thread_data.logscale * paused_thread_data.lowest_il_with_nthreads_neq_0) << endl;
  667. LOG_outfile << "Resulting info: " << scan_info << endl;
  668. }
  669. time_t current_time = time(nullptr);
  670. char timestr[100];
  671. strftime(timestr, sizeof(timestr), "%Y-%m-%d %H:%M:%S", gmtime(&current_time));
  672. LOG_outfile << "Run completion timestamp: " << timestr << " UTC" << endl;
  673. LOG_outfile << "ABACUS version " << ABACUS_VERSION << ", copyright J.-S. Caux." << endl << endl;
  674. LOG_outfile.close();
  675. }
  676. else { // in_parallel
  677. LOG_outfile << "rank " << rank << " out of " << nr_processors << " processors: "
  678. << "run info: " << scan_info << endl << "Latest threshold = "
  679. << exp(-paused_thread_data.logscale * paused_thread_data.lowest_il_with_nthreads_neq_0) << endl;
  680. }
  681. paused_thread_data.Save();
  682. ScanStateList.Order_in_SRC ();
  683. ScanStateList.Save_Info (SUM_Cstr);
  684. // Evaluate f-sumrule:
  685. if (!in_parallel ) if (whichDSF != 'q')
  686. Evaluate_F_Sumrule (prefix_prevparalevel, whichDSF, AveragingState, Chem_Pot, iKmin, iKmax);
  687. return(scan_info);
  688. }
  689. //******************************************************
  690. // Functions to initiate scans:
  691. // General version for equilibrium correlators at generic (possibly finite) temperature:
  692. Scan_Info Scan_LiebLin (char whichDSF, DP c_int, DP L, int N, int iKmin, int iKmax, DP kBT,
  693. int Max_Secs, DP target_sumrule, bool refine,
  694. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  695. {
  696. // This function scans the Hilbert space of the LiebLin gas,
  697. // for the function identified by whichDSF.
  698. // whichDSF == 'Z': canonical partition function
  699. // whichDSF == 'd': density-density correlation function
  700. // whichDSF == 'g': Green's function < \Psi \Psi^{\dagger}>
  701. // whichDSF == 'o': one-body function < \Psi^{\dagger} \Psi >
  702. // Delta is the number of sites involved in the smoothing of the entropy
  703. //int Delta = int(sqrt(N))/2;//6;//N/20;
  704. //DP epsilon = log(L)/L; // using Gaussian for density in entropy.
  705. //DP epsilon = 1.0/L; // using Lorentzian for density in entropy.
  706. // Construct the finite-size saddle-point state:
  707. // if we refine, read the quantum numbers of the saddle point state (and seed sps) from the sps file:
  708. stringstream SPS_stringstream; string SPS_string;
  709. Data_File_Name (SPS_stringstream, whichDSF, c_int, L, N, iKmin, iKmax, kBT, 0.0, "");
  710. SPS_stringstream << ".sps";
  711. SPS_string = SPS_stringstream.str(); const char* SPS_Cstr = SPS_string.c_str();
  712. fstream spsfile;
  713. if (refine) spsfile.open(SPS_Cstr, fstream::in);
  714. else spsfile.open(SPS_Cstr, fstream::out | fstream::trunc);
  715. if (spsfile.fail()) {
  716. cout << SPS_Cstr << endl; ABACUSerror("Could not open spsfile.");
  717. }
  718. LiebLin_Bethe_State spstate;
  719. if (!refine) { // obtain the sps from discretized TBA
  720. //spstate = Canonical_Saddle_Point_State (c_int, L, N, whichDSF == 'Z' ? 0.0 : kBT);
  721. spstate = Canonical_Saddle_Point_State (c_int, L, N, kBT);
  722. }
  723. else { // read it from the sps file
  724. // Check that the sps has the right number of Ix2:
  725. int Nspsread;
  726. spsfile >> Nspsread;
  727. if (Nspsread != N) {
  728. cout << Nspsread << "\t" << N << endl;
  729. ABACUSerror("Wrong number of Ix2 in saddle-point state.");
  730. }
  731. spstate = LiebLin_Bethe_State (c_int, L, N);
  732. for (int i = 0; i < N; ++i) spsfile >> spstate.Ix2[i];
  733. }
  734. spstate.Compute_All(true);
  735. int Nscan = N;
  736. if (whichDSF == 'o') Nscan = N - 1;
  737. if (whichDSF == 'g') Nscan = N + 1;
  738. // Now construct or read off the seed scan state:
  739. // TO MODIFY: this is not a good idea, since this might construct a state with many p-h w/r to the AveragingState.
  740. LiebLin_Bethe_State SeedScanState;
  741. if (whichDSF != 'o' && whichDSF != 'g') SeedScanState = spstate;
  742. else if (whichDSF == 'o' || whichDSF == 'g') {
  743. if (!refine) {
  744. if (whichDSF == 'o') SeedScanState = Remove_Particle_at_Center (spstate);
  745. else SeedScanState = Add_Particle_at_Center (spstate);
  746. }
  747. else { // read it from the sps file
  748. // Check that the sps has the right number of Ix2:
  749. int Nsspsread;
  750. spsfile >> Nsspsread;
  751. if (Nsspsread != Nscan) {
  752. cout << Nsspsread << "\t" << Nscan << endl;
  753. ABACUSerror("Wrong number of Ix2 in scan saddle-point state.");
  754. }
  755. SeedScanState = LiebLin_Bethe_State (c_int, L, Nscan);
  756. for (int i = 0; i < Nscan; ++i) spsfile >> SeedScanState.Ix2[i];
  757. }
  758. } // if one-body or Green's function
  759. SeedScanState.Compute_All(true);
  760. LiebLin_Bethe_State ScanState = SeedScanState;
  761. DP delta = sqrt(DP(N)) * (spstate.lambdaoc[N-1] - spstate.lambdaoc[0])/N;
  762. if (!refine) { // we write data to the sps file
  763. spsfile << N << endl;
  764. spsfile << spstate.Ix2 << endl;
  765. spsfile << Nscan << endl;
  766. spsfile << SeedScanState.Ix2 << endl;
  767. spsfile << endl << spstate << endl << endl;
  768. for (int i = 1; i < spstate.N - 2; ++i)
  769. spsfile << 0.5 * (spstate.lambdaoc[i] + spstate.lambdaoc[i+1])
  770. << "\t" << 1.0/spstate.L * (0.25/(spstate.lambdaoc[i] - spstate.lambdaoc[i-1])
  771. + 0.5/(spstate.lambdaoc[i+1] - spstate.lambdaoc[i])
  772. + 0.25/(spstate.lambdaoc[i+2] - spstate.lambdaoc[i+1]))
  773. << "\t" << rho_of_lambdaoc_1 (spstate, 0.5 * (spstate.lambdaoc[i] + spstate.lambdaoc[i+1]), delta)
  774. << "\t" << rho_of_lambdaoc_2 (spstate, 0.5 * (spstate.lambdaoc[i] + spstate.lambdaoc[i+1]), delta)
  775. << endl;
  776. }
  777. spsfile.close();
  778. // Perform the scan:
  779. return General_Scan (whichDSF, iKmin, iKmax, 100000000, kBT, spstate, SeedScanState, "",
  780. Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  781. }
  782. Scan_Info Scan_LiebLin (char whichDSF, DP c_int, DP L, int N, int iKmin, int iKmax, DP kBT,
  783. int Max_Secs, DP target_sumrule, bool refine)
  784. {
  785. int paralevel = 0;
  786. Vect<int> rank(0,1);
  787. Vect<int> nr_processors(0,1);
  788. return Scan_LiebLin (whichDSF, c_int, L, N, iKmin, iKmax, kBT, Max_Secs, target_sumrule,
  789. refine, paralevel, rank, nr_processors);
  790. }
  791. // Scanning on an excited state defined by a set of Ix2:
  792. Scan_Info Scan_LiebLin (char whichDSF, LiebLin_Bethe_State AveragingState, string defaultScanStatename,
  793. int iKmin, int iKmax, int Max_Secs, DP target_sumrule, bool refine,
  794. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  795. {
  796. // This function is as Scan_LiebLin for generic T defined above, except that the
  797. // averaging is now done on a state defined by AveragingStateIx2
  798. // PRECONDITIONS:
  799. // - the Ix2 of AveragingState are properly set.
  800. DP c_int = AveragingState.c_int;
  801. DP L = AveragingState.L;
  802. int N = AveragingState.N;
  803. // The label of the Averaging State is by definition the `empty' label
  804. AveragingState.Set_Label_from_Ix2 (AveragingState.Ix2);
  805. AveragingState.Compute_All(true);
  806. int Nscan = N;
  807. if (whichDSF == 'o') Nscan = N - 1;
  808. if (whichDSF == 'g') Nscan = N + 1;
  809. LiebLin_Bethe_State SeedScanState (c_int, L, Nscan);
  810. if (whichDSF == 'd' || whichDSF == 'B') SeedScanState.Ix2 = AveragingState.Ix2;
  811. // If 'o', remove midmost and shift quantum numbers by half-integer towards removed one:
  812. if (whichDSF == 'o') {
  813. for (int i = 0; i < N-1; ++i)
  814. SeedScanState.Ix2[i] = AveragingState.Ix2[i + (i >= N/2)] + 1 - 2*(i >= N/2);
  815. }
  816. // If 'g', add a quantum number in middle (explicitly: to right of index N/2)
  817. // and shift quantum numbers by half-integer away from added one:
  818. if (whichDSF == 'g') {
  819. SeedScanState.Ix2[N/2] = AveragingState.Ix2[N/2] - 1;
  820. for (int i = 0; i < N+1; ++i)
  821. SeedScanState.Ix2[i + (i >= N/2)] = AveragingState.Ix2[i] - 1 + 2*(i >= N/2);
  822. }
  823. SeedScanState.Compute_All(true);
  824. SeedScanState.Set_Label_from_Ix2 (SeedScanState.Ix2);
  825. DP kBT = 0.0;
  826. // Perform the scan:
  827. return General_Scan (whichDSF, iKmin, iKmax, 100000000, kBT,
  828. AveragingState, SeedScanState, defaultScanStatename,
  829. Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  830. }
  831. // Simplified function call of the above:
  832. Scan_Info Scan_LiebLin (char whichDSF, LiebLin_Bethe_State AveragingState, string defaultScanStatename,
  833. int iKmin, int iKmax, int Max_Secs, DP target_sumrule, bool refine)
  834. {
  835. int paralevel = 0;
  836. Vect<int> rank(0,1);
  837. Vect<int> nr_processors(0,1);
  838. return Scan_LiebLin (whichDSF, AveragingState, defaultScanStatename, iKmin, iKmax, Max_Secs,
  839. target_sumrule, refine, paralevel, rank, nr_processors);
  840. }
  841. // Scanning on a previously-defined AveragingState
  842. Scan_Info Scan_Heis (char whichDSF, XXZ_Bethe_State& AveragingState, string defaultScanStatename,
  843. int iKmin, int iKmax,
  844. int Max_Secs, DP target_sumrule, bool refine,
  845. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  846. {
  847. // General state scanning for Heisenberg chains
  848. // PRECONDITIONS:
  849. // - the Ix2 of AveragingState are properly set.
  850. // Prepare the AveragingState:
  851. AveragingState.Compute_All(true);
  852. XXZ_Bethe_State SeedScanState;
  853. if (whichDSF == 'Z' || whichDSF == 'z') SeedScanState = AveragingState;
  854. else if (whichDSF == 'm') SeedScanState = Remove_Particle_at_Center (AveragingState);
  855. else if (whichDSF == 'p') SeedScanState = Add_Particle_at_Center (AveragingState);
  856. else ABACUSerror("Unknown whichDSF in Scan_Heis.");
  857. // Now the scan itself
  858. General_Scan (whichDSF, iKmin, iKmax, AveragingState.chain.Nsites, 0.0, AveragingState, SeedScanState,
  859. defaultScanStatename, Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  860. }
  861. // Scanning on a previously-defined AveragingState
  862. Scan_Info Scan_Heis (char whichDSF, XXX_Bethe_State& AveragingState, string defaultScanStatename,
  863. int iKmin, int iKmax,
  864. int Max_Secs, DP target_sumrule, bool refine,
  865. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  866. {
  867. // General state scanning for Heisenberg chains
  868. // PRECONDITIONS:
  869. // - the Ix2 of AveragingState are properly set.
  870. // Prepare the AveragingState:
  871. AveragingState.Compute_All(true);
  872. XXX_Bethe_State SeedScanState;
  873. if (whichDSF == 'Z' || whichDSF == 'z') SeedScanState = AveragingState;
  874. else if (whichDSF == 'm') SeedScanState = Remove_Particle_at_Center (AveragingState);
  875. else if (whichDSF == 'p') SeedScanState = Add_Particle_at_Center (AveragingState);
  876. else ABACUSerror("Unknown whichDSF in Scan_Heis.");
  877. // Now the scan itself
  878. return General_Scan (whichDSF, iKmin, iKmax, AveragingState.chain.Nsites, 0.0, AveragingState, SeedScanState,
  879. defaultScanStatename, Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  880. }
  881. // Scanning on a previously-defined AveragingState
  882. Scan_Info Scan_Heis (char whichDSF, XXZ_gpd_Bethe_State& AveragingState, string defaultScanStatename,
  883. int iKmin, int iKmax,
  884. int Max_Secs, DP target_sumrule, bool refine,
  885. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  886. {
  887. // General state scanning for Heisenberg chains
  888. // PRECONDITIONS:
  889. // - the Ix2 of AveragingState are properly set.
  890. // Prepare the AveragingState:
  891. AveragingState.Compute_All(true);
  892. XXZ_gpd_Bethe_State SeedScanState;
  893. if (whichDSF == 'Z' || whichDSF == 'z') SeedScanState = AveragingState;
  894. else if (whichDSF == 'm') SeedScanState = Remove_Particle_at_Center (AveragingState);
  895. else if (whichDSF == 'p') SeedScanState = Add_Particle_at_Center (AveragingState);
  896. else ABACUSerror("Unknown whichDSF in Scan_Heis.");
  897. // Now the scan itself
  898. return General_Scan (whichDSF, iKmin, iKmax, AveragingState.chain.Nsites, 0.0, AveragingState, SeedScanState,
  899. defaultScanStatename, Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  900. }
  901. Scan_Info Scan_Heis (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  902. int Max_Secs, DP target_sumrule, bool refine,
  903. int paralevel, Vect<int> rank, Vect<int> nr_processors)
  904. {
  905. // This function scans the Hilbert space of the Heisenberg spin-1/2 chain
  906. // for the function identified by whichDSF.
  907. // whichDSF == 'Z': canonical partition function
  908. // whichDSF == 'm': S^{-+}
  909. // whichDSF == 'z': S^{zz}
  910. // whichDSF == 'p': S^{+-}
  911. // whichDSF == 'a': < S^z_j S^z_{j+1} S^z_l S^z_{l+1} > for RIXS
  912. // whichDSF == 'b': < S^z_j S^-_{j+1} S^-_l S^z_{l+1} > + (m <-> z) for RIXS
  913. // whichDSF == 'c': < S^-_j S^-_{j+1} S^-_l S^-_{l+1} > for RIXS
  914. Heis_Chain BD1(1.0, Delta, 0.0, N);
  915. Vect_INT Nrapidities_groundstate(0, BD1.Nstrings);
  916. Nrapidities_groundstate[0] = M;
  917. Heis_Base baseconfig_groundstate(BD1, Nrapidities_groundstate);
  918. if ((Delta > 0.0) && (Delta < 1.0)) {
  919. XXZ_Bethe_State GroundState(BD1, baseconfig_groundstate);
  920. GroundState.Compute_All(true);
  921. // The ground state is now fully defined.
  922. XXZ_Bethe_State SeedScanState;
  923. if (whichDSF == 'Z' || whichDSF == 'z') SeedScanState = GroundState;
  924. else if (whichDSF == 'm') SeedScanState = XXZ_Bethe_State(GroundState.chain, M - 1);
  925. else if (whichDSF == 'p') SeedScanState = XXZ_Bethe_State(GroundState.chain, M + 1);
  926. else ABACUSerror("Unknown whichDSF in Scan_Heis.");
  927. // Now the scan itself
  928. General_Scan (whichDSF, iKmin, iKmax, N, 0.0, GroundState, SeedScanState, "",
  929. Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  930. }
  931. else if (Delta == 1.0) {
  932. XXX_Bethe_State GroundState(BD1, baseconfig_groundstate);
  933. GroundState.Compute_All(true);
  934. // The ground state is now fully defined.
  935. XXX_Bethe_State SeedScanState;
  936. if (whichDSF == 'Z' || whichDSF == 'z' || whichDSF == 'a' || whichDSF == 'q') SeedScanState = GroundState;
  937. else if (whichDSF == 'm') SeedScanState = XXX_Bethe_State(GroundState.chain, M - 1);
  938. else if (whichDSF == 'p') SeedScanState = XXX_Bethe_State(GroundState.chain, M + 1);
  939. else if (whichDSF == 'c') SeedScanState = XXX_Bethe_State(GroundState.chain, M - 2);
  940. else ABACUSerror("Unknown whichDSF in Scan_Heis.");
  941. // Now the scan itself
  942. General_Scan (whichDSF, iKmin, iKmax, N, 0.0, GroundState, SeedScanState, "",
  943. Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  944. }
  945. else if (Delta > 1.0) {
  946. XXZ_gpd_Bethe_State GroundState(BD1, baseconfig_groundstate);
  947. GroundState.Compute_All(true);
  948. // The ground state is now fully defined.
  949. XXZ_gpd_Bethe_State SeedScanState;
  950. if (whichDSF == 'Z' || whichDSF == 'z') SeedScanState = GroundState;
  951. else if (whichDSF == 'm') SeedScanState = XXZ_gpd_Bethe_State(GroundState.chain, M - 1);
  952. else if (whichDSF == 'p') SeedScanState = XXZ_gpd_Bethe_State(GroundState.chain, M + 1);
  953. else ABACUSerror("Unknown whichDSF in Scan_Heis.");
  954. // Now the scan itself
  955. return General_Scan (whichDSF, iKmin, iKmax, N, 0.0, GroundState, SeedScanState, "",
  956. Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  957. }
  958. else ABACUSerror("Delta out of range in Heis_Structure_Factor");
  959. return Scan_Info();
  960. }
  961. Scan_Info Scan_Heis (char whichDSF, DP Delta, int N, int M, int iKmin, int iKmax,
  962. int Max_Secs, DP target_sumrule, bool refine)
  963. {
  964. int paralevel = 0;
  965. Vect<int> rank(0,1);
  966. Vect<int> nr_processors(0,1);
  967. return Scan_Heis (whichDSF, Delta, N, M, iKmin, iKmax,
  968. Max_Secs, target_sumrule, refine, paralevel, rank, nr_processors);
  969. }
  970. } // namespace ABACUS