CHROMA
inline_prop_distillation_stoch_w.cc
Go to the documentation of this file.
1 /*! \file
2  * \brief Compute propagators from distillation
3  *
4  * Propagator calculation in distillation
5  */
6 
7 #include "fermact.h"
10 #include "meas/glue/mesplq.h"
11 #include "qdp_map_obj.h"
12 #include "qdp_map_obj_disk.h"
13 #include "qdp_disk_map_slice.h"
15 #include "util/ferm/transf.h"
16 #include "util/ferm/spin_rep.h"
17 #include "util/ferm/diractodr.h"
19 #include "util/ft/time_slice_set.h"
20 #include "util/ft/sftmom.h"
21 #include "util/info/proginfo.h"
25 
27 
28 namespace Chroma
29 {
30  //----------------------------------------------------------------------------------
31  // Utility functions
32  namespace
33  {
34  //! Error output
35  StandardOutputStream& operator<<(StandardOutputStream& os, const multi1d<int>& d)
36  {
37  if (d.size() > 0)
38  {
39  os << d[0];
40 
41  for(int i=1; i < d.size(); ++i)
42  os << " " << d[i];
43  }
44 
45  return os;
46  }
47  }
48 
49  //----------------------------------------------------------------------------
50  namespace InlinePropDistillationStochEnv
51  {
52  //! Propagator input
53  void read(XMLReader& xml, const std::string& path, Params::NamedObject_t& input)
54  {
55  XMLReader inputtop(xml, path);
56 
57  read(inputtop, "gauge_id", input.gauge_id);
58  read(inputtop, "src_file", input.src_file);
59  read(inputtop, "soln_file", input.soln_file);
60  read(inputtop, "prop_file", input.prop_file);
61  }
62 
63  //! Propagator output
64  void write(XMLWriter& xml, const std::string& path, const Params::NamedObject_t& input)
65  {
66  push(xml, path);
67 
68  write(xml, "gauge_id", input.gauge_id);
69  write(xml, "src_file", input.src_file);
70  write(xml, "soln_file", input.soln_file);
71  write(xml, "prop_file", input.prop_file);
72 
73  pop(xml);
74  }
75 
76 
77  //! Propagator input
78  void read(XMLReader& xml, const std::string& path, Params::Param_t::Contract_t& input)
79  {
80  XMLReader inputtop(xml, path);
81 
82  read(inputtop, "num_vecs", input.num_vecs);
83  read(inputtop, "t_sources", input.t_sources);
84  read(inputtop, "decay_dir", input.decay_dir);
85  read(inputtop, "Nt_forward", input.Nt_forward);
86  read(inputtop, "Nt_backward", input.Nt_backward);
87  read(inputtop, "mass_label", input.mass_label);
88  read(inputtop, "spatial_mask_size", input.spatial_mask_size);
89  }
90 
91  //! Propagator output
92  void write(XMLWriter& xml, const std::string& path, const Params::Param_t::Contract_t& input)
93  {
94  push(xml, path);
95 
96  write(xml, "num_vecs", input.num_vecs);
97  write(xml, "t_sources", input.t_sources);
98  write(xml, "decay_dir", input.decay_dir);
99  write(xml, "Nt_forward", input.Nt_forward);
100  write(xml, "Nt_backward", input.Nt_backward);
101  write(xml, "mass_label", input.mass_label);
102  write(xml, "spatial_mask_size", input.spatial_mask_size);
103 
104  pop(xml);
105  }
106 
107 
108  //! Propagator input
109  void read(XMLReader& xml, const std::string& path, Params::Param_t& input)
110  {
111  XMLReader inputtop(xml, path);
112 
113  read(inputtop, "Propagator", input.prop);
114  read(inputtop, "Contractions", input.contract);
115  }
116 
117  //! Propagator output
118  void write(XMLWriter& xml, const std::string& path, const Params::Param_t& input)
119  {
120  push(xml, path);
121 
122  write(xml, "Propagator", input.prop);
123  write(xml, "Contractions", input.contract);
124 
125  pop(xml);
126  }
127 
128 
129  //! Propagator input
130  void read(XMLReader& xml, const std::string& path, Params& input)
131  {
132  Params tmp(xml, path);
133  input = tmp;
134  }
135 
136  //! Propagator output
137  void write(XMLWriter& xml, const std::string& path, const Params& input)
138  {
139  push(xml, path);
140 
141  write(xml, "Param", input.param);
142  write(xml, "NamedObject", input.named_obj);
143 
144  pop(xml);
145  }
146  } // namespace InlinePropDistillationStochEnv
147 
148 
149  //----------------------------------------------------------------------------
150  //----------------------------------------------------------------------------
151  namespace InlinePropDistillationStochEnv
152  {
153  // Anonymous namespace
154  namespace
155  {
156  // Convenience type
157  typedef QDP::MapObjectDisk<KeyPropDistillation_t, TimeSliceIO<LatticeColorVectorF> > MOD_t;
158 
159  //----------------------------------------------------------------------------
160  //! Get source key
161  KeyPropDistillation_t getSrcKey(int t_source, int colorvec_src)
162  {
164 
165  key.t_source = t_source;
166  key.t_slice = t_source;
167  key.colorvec_src = colorvec_src;
168  key.spin_src = -1;
169  key.spin_snk = -1;
170  // key.mass NOTE: no mass key in source
171 
172  return key;
173  }
174 
175 
176  //----------------------------------------------------------------------------
177  //! Read a source std::vector
178  LatticeColorVector getSrc(MOD_t& source_obj, int t_source, int colorvec_src)
179  {
180  QDPIO::cout << __func__ << ": on t_source= " << t_source << " colorvec_src= " << colorvec_src << std::endl;
181 
182  // Get the source std::vector
183  KeyPropDistillation_t src_key = getSrcKey(t_source, colorvec_src);
184  LatticeColorVectorF vec_srce = zero;
185 
186  TimeSliceIO<LatticeColorVectorF> time_slice_io(vec_srce, t_source);
187  source_obj.get(src_key, time_slice_io);
188 
189  return vec_srce;
190  }
191 
192 
193  //----------------------------------------------------------------------------
194  //! Get active time-slices
195  std::vector<bool> getActiveTSlices(int t_source, int Nt_forward, int Nt_backward)
196  {
197  // Initialize the active time slices
198  const int decay_dir = Nd-1;
199  const int Lt = Layout::lattSize()[decay_dir];
200 
201  std::vector<bool> active_t_slices(Lt);
202  for(int t=0; t < Lt; ++t)
203  {
204  active_t_slices[t] = false;
205  }
206 
207  // Forward
208  for(int dt=0; dt < Nt_forward; ++dt)
209  {
210  int t = t_source + dt;
211  active_t_slices[t % Lt] = true;
212  }
213 
214  // Backward
215  for(int dt=0; dt < Nt_backward; ++dt)
216  {
217  int t = t_source - dt;
218  while (t < 0) {t += Lt;}
219 
220  active_t_slices[t % Lt] = true;
221  }
222 
223  return active_t_slices;
224  }
225 
226 
227  //----------------------------------------------------------------------------
228  //! Get source keys
229  std::list<KeyPropDistillation_t> getSrcKeys(int t_source, int colorvec_src)
230  {
231  std::list<KeyPropDistillation_t> keys;
232 
233  keys.push_back(getSrcKey(t_source,colorvec_src));
234 
235  return keys;
236  }
237 
238 
239  //----------------------------------------------------------------------------
240  //! Get sink keys
241  std::list<KeyPropDistillation_t> getSnkKeys(int t_source, int colorvec_src, int Nt_forward, int Nt_backward, const std::string mass)
242  {
243  std::list<KeyPropDistillation_t> keys;
244 
245  std::vector<bool> active_t_slices = getActiveTSlices(t_source, Nt_forward, Nt_backward);
246 
247  const int decay_dir = Nd-1;
248  const int Lt = Layout::lattSize()[decay_dir];
249 
250  for(int spin_source=0; spin_source < Ns; ++spin_source)
251  {
252  for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
253  {
254  for(int t=0; t < Lt; ++t)
255  {
256  if (! active_t_slices[t]) {continue;}
257 
259 
260  key.t_source = t_source;
261  key.t_slice = t;
262  key.colorvec_src = colorvec_src;
263  key.spin_src = spin_source;
264  key.spin_snk = spin_sink;
265  key.mass = mass;
266 
267  // QDPIO::cout << key << std::flush;
268 
269  keys.push_back(key);
270  } // for t
271  } // for spin_sink
272  } // for spin_source
273 
274  return keys;
275  }
276 
277  //----------------------------------------------------------------------------
278  //! Get soln
279  multi2d<LatticeColorVector> getSoln(MOD_t& source_obj, int t_source, int colorvec_src, const std::string mass)
280  {
281  multi2d<LatticeColorVector> ferm_out(Ns,Ns);
282 
283  const int decay_dir = Nd-1;
284  const int Lt = Layout::lattSize()[decay_dir];
285 
286  for(int spin_source=0; spin_source < Ns; ++spin_source)
287  {
288  for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
289  {
291 
292  key.t_source = t_source;
293  key.t_slice = t_source;
294  key.colorvec_src = colorvec_src;
295  key.spin_src = spin_source;
296  key.spin_snk = spin_sink;
297  key.mass = mass;
298 
299  // Get the source std::vector
300  LatticeColorVectorF vec_srce = zero;
301 
302  TimeSliceIO<LatticeColorVectorF> time_slice_io(vec_srce, t_source);
303  source_obj.get(key, time_slice_io);
304 
305  ferm_out(spin_sink, spin_source) = vec_srce;
306  } // for spin_sink
307  } // for spin_source
308 
309  return ferm_out;
310  }
311 
312 
313  //----------------------------------------------------------------------------
314  //! Get soln
315  void doTrace(multi2d<ComplexD>& trace_mom,
316  const LatticeColorVector& vec_srce,
317  const multi2d<LatticeColorVector>& fred_out, int t_source, int colorvec_src,
318  const std::vector<MatrixSpinRep_t>& diracToDrMatPlus,
319  const std::vector<MatrixSpinRep_t>& diracToDrMatMinus,
320  const SftMom& phases,
321  int num_vecs)
322  {
323  // Check
324  if (1)
325  {
326  Double ddd = norm2(vec_srce); // check for debugging
327  QDPIO::cout << __func__ << ": colorvec_src= " << colorvec_src << " norm(left_vec)= " << ddd << "\n";
328  }
329 
330  // Rotate from DeGrand-Rossi (DR) to Dirac-Pauli (DP)
331  multi2d<LatticeColorVector> ferm_out;
332  {
333  multi2d<LatticeColorVector> ferm_tmp;
334 
335  multiplyRep(ferm_tmp, diracToDrMatMinus, fred_out);
336  multiplyRep(ferm_out, ferm_tmp, diracToDrMatPlus);
337  }
338 
339  for(int spin_source = 0; spin_source < Ns; ++spin_source)
340  {
341  LatticeComplex lop = localInnerProduct(vec_srce, ferm_out(spin_source,spin_source));
342 
343  // Slow fourier-transform
344  for(int mom_num = 0; mom_num < phases.numMom(); ++mom_num)
345  {
346  multi1d<ComplexD> op_sum = sumMulti(phases[mom_num] * lop, phases.getSet());
347 
348 #if 0
349  for(int ttt = 0; ttt < op_sum.size(); ++ttt)
350  {
351  QDPIO::cout << "OP_SUM: "
352  << " t_source= " << ttt
353  << " colorvec_src= " << colorvec_src
354  << " spin_source= " << spin_source
355  << " mom_num= " << mom_num
356  << " mom= " << phases.numToMom(mom_num)
357  << " op_sum[ " << ttt << " ]= " << op_sum[ttt]
358  << std::endl;
359  }
360 #endif
361 
362  trace_mom(t_source, mom_num) += op_sum[t_source];
363  }
364  }
365 
366  for(int mom_num = 0 ; mom_num < phases.numMom() ; ++mom_num)
367  {
368  QDPIO::cout << "TRACE: "
369  << " t_source= " << t_source
370  << " colorvec_src= " << colorvec_src
371  << " num_vecs= " << num_vecs
372  << " mom= " << phases.numToMom(mom_num)
373  << " trace( " << mom_num << " )= " << trace_mom(t_source,mom_num)
374  << std::endl;
375  }
376  }
377 
378  } // end anonymous
379  } // end namespace
380 
381 
382  //----------------------------------------------------------------------------
383  namespace InlinePropDistillationStochEnv
384  {
385  namespace
386  {
387  AbsInlineMeasurement* createMeasurement(XMLReader& xml_in,
388  const std::string& path)
389  {
390  return new InlineMeas(Params(xml_in, path));
391  }
392 
393  //! Local registration flag
394  bool registered = false;
395  }
396 
397  const std::string name = "PROP_DISTILLATION_STOCH";
398 
399  //! Register all the factories
400  bool registerAll()
401  {
402  bool success = true;
403  if (! registered)
404  {
406  success &= TheInlineMeasurementFactory::Instance().registerObject(name, createMeasurement);
407  registered = true;
408  }
409  return success;
410  }
411 
412 
413  //----------------------------------------------------------------------------
414  // Param stuff
416 
417  Params::Params(XMLReader& xml_in, const std::string& path)
418  {
419  try
420  {
421  XMLReader paramtop(xml_in, path);
422 
423  if (paramtop.count("Frequency") == 1)
424  read(paramtop, "Frequency", frequency);
425  else
426  frequency = 1;
427 
428  // Parameters for source construction
429  read(paramtop, "Param", param);
430 
431  // Read in the output propagator/source configuration info
432  read(paramtop, "NamedObject", named_obj);
433 
434  // Possible alternate XML file pattern
435  if (paramtop.count("xml_file") != 0)
436  {
437  read(paramtop, "xml_file", xml_file);
438  }
439  }
440  catch(const std::string& e)
441  {
442  QDPIO::cerr << __func__ << ": Caught Exception reading XML: " << e << std::endl;
443  QDP_abort(1);
444  }
445  }
446 
447 
448  //----------------------------------------------------------------------------
449  //----------------------------------------------------------------------------
450  // Function call
451  void
452  InlineMeas::operator()(unsigned long update_no,
453  XMLWriter& xml_out)
454  {
455  // If xml file not empty, then use alternate
456  if (params.xml_file != "")
457  {
458  std::string xml_file = makeXMLFileName(params.xml_file, update_no);
459 
460  push(xml_out, "PropDistillation");
461  write(xml_out, "update_no", update_no);
462  write(xml_out, "xml_file", xml_file);
463  pop(xml_out);
464 
465  XMLFileWriter xml(xml_file);
466  func(update_no, xml);
467  }
468  else
469  {
470  func(update_no, xml_out);
471  }
472  }
473 
474 
475  // Real work done here
476  void
477  InlineMeas::func(unsigned long update_no,
478  XMLWriter& xml_out)
479  {
480  START_CODE();
481 
482  StopWatch snoop;
483  snoop.reset();
484  snoop.start();
485 
486  // Test and grab a reference to the gauge field
487  multi1d<LatticeColorMatrix> u;
488  XMLBufferWriter gauge_xml;
489  try
490  {
491  u = TheNamedObjMap::Instance().getData< multi1d<LatticeColorMatrix> >(params.named_obj.gauge_id);
492  TheNamedObjMap::Instance().get(params.named_obj.gauge_id).getRecordXML(gauge_xml);
493  }
494  catch( std::bad_cast )
495  {
496  QDPIO::cerr << name << ": caught dynamic cast error" << std::endl;
497  QDP_abort(1);
498  }
499  catch (const std::string& e)
500  {
501  QDPIO::cerr << name << ": std::map call failed: " << e << std::endl;
502  QDP_abort(1);
503  }
504 
505  push(xml_out, "PropDistillation");
506  write(xml_out, "update_no", update_no);
507 
508  QDPIO::cout << name << ": propagator calculation" << std::endl;
509 
510  proginfo(xml_out); // Print out basic program info
511 
512  // Write out the input
513  write(xml_out, "Input", params);
514 
515  // Write out the config header
516  write(xml_out, "Config_info", gauge_xml);
517 
518  push(xml_out, "Output_version");
519  write(xml_out, "out_version", 1);
520  pop(xml_out);
521 
522  // Calculate some gauge invariant observables just for info.
523  MesPlq(xml_out, "Observables", u);
524 
525  // Will use TimeSliceSet-s a lot
526  const int decay_dir = params.param.contract.decay_dir;
527  const int Lt = Layout::lattSize()[decay_dir];
528 
529  // A sanity check
530  if (decay_dir != Nd-1)
531  {
532  QDPIO::cerr << __func__ << ": TimeSliceIO only supports decay_dir= " << Nd-1 << "\n";
533  QDP_abort(1);
534  }
535 
536  // The time-slice set
537  TimeSliceSet time_slice_set(decay_dir);
538 
539  //
540  // Sanity checks
541  //
542  if (params.param.contract.spatial_mask_size.size() != Nd-1)
543  {
544  QDPIO::cerr << name << ": param.contract.spatial mask size incorrect 1" << std::endl;
545  QDP_abort(1);
546  }
547 
548  for(int mu=0; mu < Nd-1; ++mu)
549  {
551  {
552  QDPIO::cerr << name << ": not valid param.contract.spatial mask" << std::endl;
553  QDP_abort(1);
554  }
555 
556  if ((QDP::Layout::lattSize()[mu] % params.param.contract.spatial_mask_size[mu]) != 0)
557  {
558  QDPIO::cerr << name << ": lattice size not divisible by param.contract.spatial mask size" << std::endl;
559  QDP_abort(1);
560  }
561  }
562 
563 
564  //
565  // Build a list of all the colorvec sources
566  //
567  std::list<int> source_list;
568 
569  if (1)
570  {
571 #if 0
572  // The distillation sources
573  for(int colorvec_src = 0; colorvec_src < num_vecs; ++colorvec_src)
574  source_list.push_back(colorvec_src);
575 #endif
576 
577  // The number of spatial noises
578  int num_sites = 1;
579  for(int mu = 0; mu < params.param.contract.spatial_mask_size.size(); ++mu)
580  {
581  num_sites *= params.param.contract.spatial_mask_size[mu];
582  }
583 
584  // Loop over all sites with a spatial mask
585  // These will form the sources
586 #if 0
587  for(int ipos=1; ipos <= num_sites; ++ipos)
588 #endif
589  for(int ipos = num_sites; ipos > 0; --ipos)
590  source_list.push_back(-ipos);
591  }
592 
593 
594  //
595  // Map-object-disk storage of the source file
596  //
597  QDP::MapObjectDisk<KeyPropDistillation_t, TimeSliceIO<LatticeColorVectorF> > source_obj;
598  source_obj.setDebug(0);
599 
600  QDPIO::cout << "Open source file" << std::endl;
601 
602  if (! source_obj.fileExists(params.named_obj.src_file))
603  {
604  QDPIO::cerr << name << ": source file does not exist: src_file= " << params.named_obj.src_file << std::endl;
605  QDP_abort(1);
606  }
607  else
608  {
609  source_obj.open(params.named_obj.src_file, std::ios_base::in);
610  }
611 
612  QDPIO::cout << "Finished opening solution file" << std::endl;
613 
614 
615  //
616  // Map-object-disk storage
617  //
618  QDP::MapObjectDisk<KeyPropDistillation_t, TimeSliceIO<LatticeColorVectorF> > soln_obj;
619  soln_obj.setDebug(0);
620 
621  QDPIO::cout << "Open solution file" << std::endl;
622 
623  if (! soln_obj.fileExists(params.named_obj.soln_file))
624  {
625  QDPIO::cerr << name << ": soln file does not exist: soln_file= " << params.named_obj.soln_file << std::endl;
626  QDP_abort(1);
627  }
628  else
629  {
631  }
632 
633  QDPIO::cout << "Finished opening solution file" << std::endl;
634 
635 
636  //
637  // Map-object-disk storage
638  //
639  QDP::MapObjectDisk<KeyPropDistillation_t, TimeSliceIO<LatticeColorVectorF> > prop_obj;
640  prop_obj.setDebug(0);
641 
642  QDPIO::cout << "Open solution file" << std::endl;
643 
644  if (! prop_obj.fileExists(params.named_obj.prop_file))
645  {
646  XMLBufferWriter file_xml;
647 
648  push(file_xml, "MODMetaData");
649  write(file_xml, "id", std::string("propDistillation"));
650  write(file_xml, "lattSize", QDP::Layout::lattSize());
651  write(file_xml, "decay_dir", decay_dir);
652  write(file_xml, "num_vecs", params.param.contract.num_vecs);
653  write(file_xml, "Nt_backward", params.param.contract.Nt_backward);
654  write(file_xml, "Nt_backward", params.param.contract.Nt_backward);
655  write(file_xml, "mass_label", params.param.contract.mass_label);
656  proginfo(file_xml); // Print out basic program info
657  write(file_xml, "Params", params.param);
658  write(file_xml, "Config_info", gauge_xml);
659  pop(file_xml);
660 
661  std::string file_str(file_xml.str());
662 
663  prop_obj.insertUserdata(file_xml.str());
664  prop_obj.open(params.named_obj.prop_file, std::ios_base::in | std::ios_base::out | std::ios_base::trunc);
665  }
666  else
667  {
668  prop_obj.open(params.named_obj.prop_file);
669  }
670 
671  QDPIO::cout << "Finished opening prop file" << std::endl;
672 
673 
674  // Total number of iterations
675  int ncg_had = 0;
676 
677 
678 #if 1
679  // NOTE: not doing this spin rotation, but leaving the code here for reference
680  // Rotation from DR to DP
681  SpinMatrix diracToDRMat(DiracToDRMat());
682  std::vector<MatrixSpinRep_t> diracToDrMatPlus = convertTwoQuarkSpin(diracToDRMat);
683  std::vector<MatrixSpinRep_t> diracToDrMatMinus = convertTwoQuarkSpin(adj(diracToDRMat));
684 #endif
685 
686 
687 #if 1
688  SftMom phases(2, false, decay_dir);
689  multi2d<ComplexD> trace_mom(phases.numSubsets(), phases.numMom());
690  trace_mom = zero;
691 #endif
692 
693  //
694  // Try the factories
695  //
696  try
697  {
698  StopWatch swatch;
699  swatch.reset();
700  QDPIO::cout << "Try the various factories" << std::endl;
701 
702  // Typedefs to save typing
703  typedef LatticeFermion T;
704  typedef multi1d<LatticeColorMatrix> P;
705  typedef multi1d<LatticeColorMatrix> Q;
706 
707  //
708  // Initialize fermion action
709  //
710  std::istringstream xml_s(params.param.prop.fermact.xml);
711  XMLReader fermacttop(xml_s);
712  QDPIO::cout << "FermAct = " << params.param.prop.fermact.id << std::endl;
713 
714  // Generic Wilson-Type stuff
717  fermacttop,
719 
720  Handle< FermState<T,P,Q> > state(S_f->createState(u));
721 
724 
725  QDPIO::cout << "Suitable factory found: compute all the quark props" << std::endl;
726  swatch.start();
727 
728 
729  //
730  // Loop over the source color and spin, creating the source
731  // and calling the relevant propagator routines.
732  //
733  const int num_vecs = params.param.contract.num_vecs;
734  const multi1d<int>& t_sources = params.param.contract.t_sources;
735 
736  // Loop over each time source
737  for(int tt=0; tt < t_sources.size(); ++tt)
738  {
739  int t_source = t_sources[tt]; // This is the actual time-slice.
740  QDPIO::cout << "t_source = " << t_source << std::endl;
741 
742 #if 1
743  if (1)
744  {
745  for(int colorvec_src = 0; colorvec_src < num_vecs; ++colorvec_src)
746  {
747  doTrace(trace_mom,
748  getSrc(source_obj, t_source, colorvec_src),
749  getSoln(soln_obj, t_source, colorvec_src, params.param.contract.mass_label),
750  t_source, colorvec_src,
751  diracToDrMatPlus,
752  diracToDrMatMinus,
753  phases,
754  num_vecs);
755 
756  }
757  }
758 #endif
759 
760  //
761  // The space distillation loop
762  //
763  for(auto col_src = source_list.begin(); col_src != source_list.end(); ++col_src)
764  {
765  const int colorvec_src = *col_src;
766 
767  StopWatch sniss1;
768  sniss1.reset();
769  sniss1.start();
770  QDPIO::cout << "colorvec_src = " << colorvec_src << std::endl;
771 
772  // Get the source std::vector
773  LatticeColorVector vec_srce = getSrc(source_obj, t_source, colorvec_src);
774 
775  // Check
776  if (1)
777  {
778  Double ddd = norm2(vec_srce); // check for debugging
779  QDPIO::cout << __func__ << ": colorvec_src= " << colorvec_src << " norm(source_vec)= " << ddd << "\n";
780  }
781 
782  //
783  // Loop over each spin source and invert.
784  // Use the same colorstd::vector source. No spin dilution will be used.
785  //
786  multi2d<LatticeColorVector> ferm_out(Ns,Ns);
787 
788  for(int spin_source=0; spin_source < Ns; ++spin_source)
789  {
790  QDPIO::cout << "spin_source = " << spin_source << std::endl;
791 
792  // Insert a ColorVector into spin index spin_source
793  // This only overwrites sections, so need to initialize first
794  LatticeFermion chi = zero;
795  CvToFerm(vec_srce, chi, spin_source);
796 
797  LatticeFermion quark_soln = zero;
798 
799  // Do the propagator inversion
800  SystemSolverResults_t res = (*PP)(quark_soln, chi);
801  ncg_had = res.n_count;
802 
803  // Extract into the temporary output array
804  for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
805  {
806  ferm_out(spin_sink,spin_source) = peekSpin(quark_soln, spin_sink);
807  }
808  } // for spin_source
809 
810 
811 #if 0
812  // NOTE: not doing this spin rotation, but leaving the code here for reference
813  // Instead, leave the code in the DR basis - the convention for all the perambulators
814  // written
815 
816  // Rotate from DeGrand-Rossi (DR) to Dirac-Pauli (DP)
817  {
818  multi2d<LatticeColorVector> ferm_tmp;
819 
820  multiplyRep(ferm_tmp, diracToDrMatMinus, ferm_out);
821  multiplyRep(ferm_out, ferm_tmp, diracToDrMatPlus);
822  }
823 #endif
824 
825  sniss1.stop();
826  QDPIO::cout << "Time to assemble and transmogrify propagators for colorvec_src= " << colorvec_src << " time = "
827  << sniss1.getTimeInSeconds()
828  << " secs" << std::endl;
829 
830 
831  // Write out each time-slice chunk of a lattice colorvec soln to disk
832  QDPIO::cout << "Write propagator solutions to disk" << std::endl;
833  StopWatch sniss2;
834  sniss2.reset();
835  sniss2.start();
836 
837  // Write the solutions
838  QDPIO::cout << "Write propagator solution to disk" << std::endl;
839  std::list<KeyPropDistillation_t> snk_keys(getSnkKeys(t_source, colorvec_src,
843 
844  for(std::list<KeyPropDistillation_t>::const_iterator key= snk_keys.begin();
845  key != snk_keys.end();
846  ++key)
847  {
848  LatticeColorVectorF tmptmp = ferm_out(key->spin_snk,key->spin_src);
849 
850  prop_obj.insert(*key, TimeSliceIO<LatticeColorVectorF>(tmptmp, key->t_slice));
851  } // for key
852 
853 
854 #if 1
855  //
856  // Compute some time-sliced 1pt traces
857  //
858  if (1)
859  {
860  // Get the source std::vector
861  if (colorvec_src < 0)
862  {
863  int orig_vec = (-colorvec_src) | (1 << 15);
864  vec_srce = getSrc(source_obj, t_source, -orig_vec);
865  }
866 
867  doTrace(trace_mom,
868  vec_srce,
869  ferm_out,
870  t_source, colorvec_src,
871  diracToDrMatPlus,
872  diracToDrMatMinus,
873  phases,
874  num_vecs);
875  }
876 #endif
877 
878  sniss2.stop();
879  QDPIO::cout << "Time to write propagators for colorvec_src= " << colorvec_src << " time = "
880  << sniss2.getTimeInSeconds()
881  << " secs" << std::endl;
882 
883  } // for colorvec_src
884  } // for tt
885 
886  swatch.stop();
887  QDPIO::cout << "Propagators computed: time= "
888  << swatch.getTimeInSeconds()
889  << " secs" << std::endl;
890  }
891  catch (const std::string& e)
892  {
893  QDPIO::cout << name << ": caught exception around qprop: " << e << std::endl;
894  QDP_abort(1);
895  }
896 
897  push(xml_out,"Relaxation_Iterations");
898  write(xml_out, "ncg_had", ncg_had);
899  pop(xml_out);
900 
901  pop(xml_out); // prop_dist
902 
903  snoop.stop();
904  QDPIO::cout << name << ": total time = "
905  << snoop.getTimeInSeconds()
906  << " secs" << std::endl;
907 
908  QDPIO::cout << name << ": ran successfully" << std::endl;
909 
910  END_CODE();
911  }
912 
913  }
914 
915 } // namespace Chroma
Inline measurement factory.
Class for counted reference semantics.
Definition: handle.h:33
void func(const unsigned long update_no, XMLWriter &xml_out)
Do the measurement.
void operator()(const unsigned long update_no, XMLWriter &xml_out)
Do the measurement.
Fourier transform phase factor support.
Definition: sftmom.h:35
int numSubsets() const
Number of subsets - length in decay direction.
Definition: sftmom.h:63
multi1d< int > numToMom(int mom_num) const
Convert momenta id to actual array of momenta.
Definition: sftmom.h:78
int numMom() const
Number of momenta.
Definition: sftmom.h:60
const Set & getSet() const
The set to be used in sumMulti.
Definition: sftmom.h:57
static T & Instance()
Definition: singleton.h:432
Builds time slice subsets.
EXTERN Real dt
int mu
Definition: cool.cc:24
Basis rotation matrix from Dirac to Degrand-Rossi (and reverse)
Class structure for fermion actions.
Fermion action factories.
All Wilson-type fermion actions.
SpinMatrixD DiracToDRMat()
The Dirac to Degrand-Rossi spin transformation matrix (and reverse)
Definition: diractodr.cc:22
void CvToFerm(const LatticeColorVectorF &a, LatticeFermionF &b, int spin_index)
Convert (insert) a LatticeColorVector into a LatticeFermion.
Definition: transf.cc:18
void proginfo(XMLWriter &xml)
Print out basic information about this program.
Definition: proginfo.cc:24
std::string makeXMLFileName(std::string xml_file, unsigned long update_no)
Return a xml file name for inline measurements.
TimeSliceSet time_slice_set
Compute the propagator from distillation.
Key for vanilla distillation propagator sources and solutions.
Make xml file writer.
int t
Definition: meslate.cc:37
Nd
Definition: meslate.cc:74
Named object function std::map.
static bool registered
Local registration flag.
multi1d< LatticeColorMatrix > P
QDPSubTypeTrait< typename BinaryReturn< C1, C2, FnLocalInnerProduct >::Type_t >::Type_t localInnerProduct(const QDPSubType< T1, C1 > &l, const QDPType< T2, C2 > &r)
void read(XMLReader &xml, const std::string &path, Params::NamedObject_t &input)
Propagator input.
void write(XMLWriter &xml, const std::string &path, const Params::NamedObject_t &input)
Propagator output.
bool registerAll()
Register all the factories.
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
QDP::StandardOutputStream & operator<<(QDP::StandardOutputStream &s, const multi1d< int > &d)
Definition: npr_vertex_w.cc:12
static multi1d< LatticeColorMatrix > u
LinOpSysSolverMGProtoClover::Q Q
LatticeFermion tmp
Definition: mespbg5p_w.cc:36
push(xml_out,"Condensates")
LinOpSysSolverMGProtoClover::T T
int i
Definition: pbg5p_w.cc:55
Double mass
Definition: pbg5p_w.cc:54
std::vector< MatrixSpinRep_t > convertTwoQuarkSpin(const SpinMatrix &in)
Convert a generic spin matrix into a sparse spin representation.
Definition: spin_rep.cc:67
multi1d< LatticeFermion > chi(Ncb)
QDP::MapObjectDisk< KeyPropDistillution_t, TimeSliceIO< LatticeColorVectorF > > MOD_t
pop(xml_out)
void MesPlq(const multi1d< LatticeColorMatrixF3 > &u, multi2d< Double > &plane_plaq, Double &link)
Definition: mesplq.cc:83
DComplex d
Definition: invbicg.cc:99
static QDP_ColorVector * out
Constructor.
START_CODE()
const WilsonTypeFermAct< multi1d< LatticeFermion > > Handle< const ConnectState > state
Definition: pbg5p_w.cc:28
const WilsonTypeFermAct< multi1d< LatticeFermion > > & S_f
Definition: pbg5p_w.cc:27
Double zero
Definition: invbicg.cc:106
void multiplyRep(multi2d< LatticeColorVector > &dist_rep, const std::vector< MatrixSpinRep_t > &spin, const multi2d< LatticeColorVector > &prop_rep)
Dist(t2) = SpinMatrix*Prop(t2)
static QDP_ColorVector * in
FloatingPoint< double > Double
Definition: gtest.h:7351
::std::string string
Definition: gtest.h:1979
Print out basic info about this program.
Fourier transform phase factor support.
Sparse matrix representation of spin matrices.
GroupXML_t invParam
Definition: qprop_io.h:84
GroupXML_t fermact
Definition: qprop_io.h:80
Distillation propagators.
Holds return info from SystemSolver call.
Definition: syssolver.h:17
Convenience for building time-slice subsets.
Contraction operators for two quarks.