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 46KB

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