CHROMA
inline_prop_and_matelem_colorvec_w.cc
Go to the documentation of this file.
1 /*! \file
2  * \brief Compute the perambulators
3  *
4  * Compute propagator elements M^-1 * multi1d<LatticeColorVector>
5  * and the matrix element of LatticeColorVector*M^-1*LatticeColorVector
6  */
7 
8 #include "fermact.h"
11 #include "meas/glue/mesplq.h"
17 #include "util/ferm/key_val_db.h"
18 #include "util/ferm/transf.h"
19 #include "util/ft/sftmom.h"
20 #include "util/info/proginfo.h"
24 
26 
27 namespace Chroma
28 {
29  namespace InlinePropAndMatElemColorVecEnv
30  {
31  //! Propagator input
32  void read(XMLReader& xml, const std::string& path, Params::NamedObject_t& input)
33  {
34  XMLReader inputtop(xml, path);
35 
36  read(inputtop, "gauge_id", input.gauge_id);
37  read(inputtop, "colorvec_id", input.colorvec_id);
38  read(inputtop, "prop_id", input.prop_id);
39  read(inputtop, "prop_op_file", input.prop_op_file);
40 
41  // User Specified MapObject tags
42  input.prop_obj = readXMLGroup(inputtop, "PropMapObject", "MapObjType");
43  }
44 
45  //! Propagator output
46  void write(XMLWriter& xml, const std::string& path, const Params::NamedObject_t& input)
47  {
48  push(xml, path);
49 
50  write(xml, "gauge_id", input.gauge_id);
51  write(xml, "colorvec_id", input.colorvec_id);
52  write(xml, "prop_id", input.prop_id);
53  write(xml, "prop_op_file", input.prop_op_file);
54  xml << input.prop_obj.xml;
55 
56  pop(xml);
57  }
58 
59 
60  //! Propagator input
61  void read(XMLReader& xml, const std::string& path, Params::Param_t::Contract_t& input)
62  {
63  XMLReader inputtop(xml, path);
64 
65  read(inputtop, "num_vecs", input.num_vecs);
66  read(inputtop, "t_sources", input.t_sources);
67  read(inputtop, "decay_dir", input.decay_dir);
68  read(inputtop, "mass_label", input.mass_label);
69  }
70 
71  //! Propagator output
72  void write(XMLWriter& xml, const std::string& path, const Params::Param_t::Contract_t& input)
73  {
74  push(xml, path);
75 
76  write(xml, "num_vecs", input.num_vecs);
77  write(xml, "t_sources", input.t_sources);
78  write(xml, "decay_dir", input.decay_dir);
79  write(xml, "mass_label", input.mass_label);
80 
81  pop(xml);
82  }
83 
84 
85  //! Propagator input
86  void read(XMLReader& xml, const std::string& path, Params::Param_t& input)
87  {
88  XMLReader inputtop(xml, path);
89 
90  read(inputtop, "Propagator", input.prop);
91  read(inputtop, "Contractions", input.contract);
92  }
93 
94  //! Propagator output
95  void write(XMLWriter& xml, const std::string& path, const Params::Param_t& input)
96  {
97  push(xml, path);
98 
99  write(xml, "Propagator", input.prop);
100  write(xml, "Contractions", input.contract);
101 
102  pop(xml);
103  }
104 
105 
106  //! Propagator input
107  void read(XMLReader& xml, const std::string& path, Params& input)
108  {
109  Params tmp(xml, path);
110  input = tmp;
111  }
112 
113  //! Propagator output
114  void write(XMLWriter& xml, const std::string& path, const Params& input)
115  {
116  push(xml, path);
117 
118  write(xml, "Param", input.param);
119  write(xml, "NamedObject", input.named_obj);
120 
121  pop(xml);
122  }
123 
124 
125  //---------------------------------------------------------------------------------------------
126  namespace
127  {
128  AbsInlineMeasurement* createMeasurement(XMLReader& xml_in,
129  const std::string& path)
130  {
131  return new InlineMeas(Params(xml_in, path));
132  }
133 
134  //! Local registration flag
135  bool registered = false;
136 
137  const std::string name = "PROP_AND_MATELEM_COLORVEC";
138  }
139 
140  //! Register all the factories
141  bool registerAll()
142  {
143  bool success = true;
144  if (! registered)
145  {
147  success &= TheInlineMeasurementFactory::Instance().registerObject(name, createMeasurement);
149  registered = true;
150  }
151  return success;
152  }
153 
154 
155  //----------------------------------------------------------------------------
156  // Param stuff
158 
159  Params::Params(XMLReader& xml_in, const std::string& path)
160  {
161  try
162  {
163  XMLReader paramtop(xml_in, path);
164 
165  if (paramtop.count("Frequency") == 1)
166  read(paramtop, "Frequency", frequency);
167  else
168  frequency = 1;
169 
170  // Parameters for source construction
171  read(paramtop, "Param", param);
172 
173  // Read in the output propagator/source configuration info
174  read(paramtop, "NamedObject", named_obj);
175 
176  // Possible alternate XML file pattern
177  if (paramtop.count("xml_file") != 0)
178  {
179  read(paramtop, "xml_file", xml_file);
180  }
181  }
182  catch(const std::string& e)
183  {
184  QDPIO::cerr << name << ": Caught Exception reading XML: " << e << std::endl;
185  QDP_abort(1);
186  }
187  }
188 
189 
190 
191  // Function call
192  void
193  InlineMeas::operator()(unsigned long update_no,
194  XMLWriter& xml_out)
195  {
196  // If xml file not empty, then use alternate
197  if (params.xml_file != "")
198  {
199  std::string xml_file = makeXMLFileName(params.xml_file, update_no);
200 
201  push(xml_out, "PropColorVec");
202  write(xml_out, "update_no", update_no);
203  write(xml_out, "xml_file", xml_file);
204  pop(xml_out);
205 
206  XMLFileWriter xml(xml_file);
207  func(update_no, xml);
208  }
209  else
210  {
211  func(update_no, xml_out);
212  }
213  }
214 
215 
216  // Real work done here
217  void
218  InlineMeas::func(unsigned long update_no,
219  XMLWriter& xml_out)
220  {
221  START_CODE();
222 
223  StopWatch snoop;
224  snoop.reset();
225  snoop.start();
226 
227  // Test and grab a reference to the gauge field
228  multi1d<LatticeColorMatrix> u;
229  XMLBufferWriter gauge_xml;
230  try
231  {
232  u = TheNamedObjMap::Instance().getData< multi1d<LatticeColorMatrix> >(params.named_obj.gauge_id);
233  TheNamedObjMap::Instance().get(params.named_obj.gauge_id).getRecordXML(gauge_xml);
234  }
235  catch( std::bad_cast )
236  {
237  QDPIO::cerr << name << ": caught dynamic cast error" << std::endl;
238  QDP_abort(1);
239  }
240  catch (const std::string& e)
241  {
242  QDPIO::cerr << name << ": std::map call failed: " << e << std::endl;
243  QDP_abort(1);
244  }
245 
246  push(xml_out, "PropColorVec");
247  write(xml_out, "update_no", update_no);
248 
249  QDPIO::cout << name << ": propagator calculation" << std::endl;
250 
251  proginfo(xml_out); // Print out basic program info
252 
253  // Write out the input
254  write(xml_out, "Input", params);
255 
256  // Write out the config header
257  write(xml_out, "Config_info", gauge_xml);
258 
259  push(xml_out, "Output_version");
260  write(xml_out, "out_version", 1);
261  pop(xml_out);
262 
263  // Calculate some gauge invariant observables just for info.
264  MesPlq(xml_out, "Observables", u);
265 
266  //
267  // Read in the source along with relevant information.
268  //
269  XMLReader source_file_xml, source_record_xml;
270 
271  QDPIO::cout << "Snarf the source from a named buffer" << std::endl;
272  try
273  {
275 
276  // Snarf the source info. This is will throw if the colorvec_id is not there
277 
278  TheNamedObjMap::Instance().get(params.named_obj.colorvec_id).getFileXML(source_file_xml);
279  TheNamedObjMap::Instance().get(params.named_obj.colorvec_id).getRecordXML(source_record_xml);
280 
281  // Write out the source header
282  write(xml_out, "Source_file_info", source_file_xml);
283  write(xml_out, "Source_record_info", source_record_xml);
284  }
285  catch (std::bad_cast) {
286  QDPIO::cerr << name << ": caught dynamic cast error" << std::endl;
287  QDP_abort(1);
288  }
289  catch (const std::string& e) {
290  QDPIO::cerr << name << ": error extracting source_header: " << e << std::endl;
291  QDP_abort(1);
292  }
293  catch( const char* e) {
294  QDPIO::cerr << name <<": Caught some char* exception:" << std::endl;
295  QDPIO::cerr << e << std::endl;
296  QDPIO::cerr << "Rethrowing" << std::endl;
297  throw;
298  }
299 
300  // Cast should be valid now
301  const QDP::MapObject<int,EVPair<LatticeColorVector> >& eigen_source =
303 
304  QDPIO::cout << "Source successfully read and parsed" << std::endl;
305 
306  //
307  // Create the output files
308  //
309  try
310  {
311  // Generate a metadata
312  std::string file_str;
313  if (1)
314  {
315  XMLBufferWriter file_xml;
316 
317  push(file_xml, "MODMetaData");
318  write(file_xml, "id", std::string("propColorVec"));
319  write(file_xml, "lattSize", QDP::Layout::lattSize());
320  write(file_xml, "decay_dir", params.param.contract.decay_dir);
321  write(file_xml, "num_vecs", params.param.contract.num_vecs);
322  // write(file_xml, "mass_label", params.param.contract.mass_label); // no simple kind of mass label available. That is unfortunate...
323  proginfo(file_xml); // Print out basic program info
324  write(file_xml, "Propagator", params.param.prop);
325  write(file_xml, "Contractions", params.param.contract);
326  write(file_xml, "Config_info", gauge_xml);
327  pop(file_xml);
328 
329  file_str = file_xml.str();
330  }
331 
332  // Create the map object
333  std::istringstream xml_s(params.named_obj.prop_obj.xml);
334  XMLReader MapObjReader(xml_s);
335 
336  // Create the entry
340  MapObjReader,
342  file_str);
343  }
344  catch (std::bad_cast)
345  {
346  QDPIO::cerr << name << ": caught dynamic cast error" << std::endl;
347  QDP_abort(1);
348  }
349  catch (const std::string& e)
350  {
351  QDPIO::cerr << name << ": error creating prop: " << e << std::endl;
352  QDP_abort(1);
353  }
354 
355  // Cast should be valid now
356  QDP::MapObject<KeyPropColorVec_t,LatticeFermion>& prop_obj =
358 
359  // Sanity check - write out the norm2 of the source in the Nd-1 direction
360  // Use this for any possible verification
361  {
362  // Initialize the slow Fourier transform phases
363  SftMom phases(0, true, Nd-1);
364 
365  EVPair<LatticeColorVector> tmpvec; eigen_source.get(0, tmpvec);
366  multi1d<Double> source_corrs = sumMulti(localNorm2(tmpvec.eigenVector), phases.getSet());
367 
368  push(xml_out, "Source_correlators");
369  write(xml_out, "source_corrs", source_corrs);
370  pop(xml_out);
371  }
372 
373  // Another sanity check
375  {
376  QDPIO::cerr << __func__ << ": num_vecs= " << params.param.contract.num_vecs
377  << " is greater than the number of available colorvectors= "
378  << eigen_source.size() << std::endl;
379  QDP_abort(1);
380  }
381 
382 
383  //
384  // DB storage
385  //
386  BinaryStoreDB< SerialDBKey<KeyPropElementalOperator_t>, SerialDBData<ValPropElementalOperator_t> > qdp_db;
387 
388  // Open the file, and write the meta-data and the binary for this operator
389  if (! qdp_db.fileExists(params.named_obj.prop_op_file))
390  {
391  XMLBufferWriter file_xml;
392 
393  push(file_xml, "DBMetaData");
394  write(file_xml, "id", std::string("propElemOp"));
395  write(file_xml, "lattSize", QDP::Layout::lattSize());
396  write(file_xml, "decay_dir", params.param.contract.decay_dir);
397  proginfo(file_xml); // Print out basic program info
398  write(file_xml, "Params", params.param);
399  write(file_xml, "Config_info", gauge_xml);
401  pop(file_xml);
402 
403  std::string file_str(file_xml.str());
404  qdp_db.setMaxUserInfoLen(file_str.size());
405 
406  qdp_db.open(params.named_obj.prop_op_file, O_RDWR | O_CREAT, 0664);
407 
408  qdp_db.insertUserdata(file_str);
409  }
410  else
411  {
412  qdp_db.open(params.named_obj.prop_op_file, O_RDWR, 0664);
413  }
414 
415 
416  // Total number of iterations
417  int ncg_had = 0;
418 
419  //
420  // Try the factories
421  //
422  try
423  {
424  StopWatch swatch;
425  swatch.reset();
426  QDPIO::cout << "Try the various factories" << std::endl;
427 
428  // Typedefs to save typing
429  typedef LatticeFermion T;
430  typedef multi1d<LatticeColorMatrix> P;
431  typedef multi1d<LatticeColorMatrix> Q;
432 
433  //
434  // Initialize fermion action
435  //
436  std::istringstream xml_s(params.param.prop.fermact.xml);
437  XMLReader fermacttop(xml_s);
438  QDPIO::cout << "FermAct = " << params.param.prop.fermact.id << std::endl;
439 
440  // Generic Wilson-Type stuff
443  fermacttop,
445 
446  Handle< FermState<T,P,Q> > state(S_f->createState(u));
447 
450 
451  QDPIO::cout << "Suitable factory found: compute all the quark props" << std::endl;
452  swatch.start();
453 
454  //
455  // Loop over the source color and spin, creating the source
456  // and calling the relevant propagator routines.
457  //
458  const int num_vecs = params.param.contract.num_vecs;
459  const int decay_dir = params.param.contract.decay_dir;
460  const multi1d<int>& t_sources = params.param.contract.t_sources;
461 
462  // Initialize the slow Fourier transform phases
463  SftMom phases(0, true, decay_dir);
464 
465 
466  // Loop over each operator
467  for(int tt=0; tt < t_sources.size(); ++tt)
468  {
469  int t_source = t_sources[tt];
470  QDPIO::cout << "t_source = " << t_source << std::endl;
471 
472  // All the loops
473  for(int spin_source=0; spin_source < Ns; ++spin_source)
474  {
475  QDPIO::cout << "spin_source = " << spin_source << std::endl;
476 
477  //
478  // Initialize all the keys for this spin spin and all spin sinks on all t-slices
479  //
480  multi2d<KeyValPropElementalOperator_t> buf(phases.numSubsets(),Ns);
481  for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
482  {
483  for(int t=0; t < phases.numSubsets(); ++t)
484  {
485  buf(t,spin_sink).key.key().t_slice = t;
486  buf(t,spin_sink).key.key().t_source = t_source;
487  buf(t,spin_sink).key.key().spin_src = spin_source;
488  buf(t,spin_sink).key.key().spin_snk = spin_sink;
489  buf(t,spin_sink).key.key().mass_label = params.param.contract.mass_label;
490  buf(t,spin_sink).val.data().mat.resize(num_vecs,num_vecs);
491  }
492  }
493 
494  // Propagator solution part
495  for(int colorvec_source=0; colorvec_source < num_vecs; ++colorvec_source)
496  {
497  QDPIO::cout << "colorvec_source = " << colorvec_source << std::endl;
498 
499  // Pull out a time-slice of the color std::vector source
500  LatticeColorVector vec_srce = zero;
501  {
503  eigen_source.get(colorvec_source, tmpvec);
504  vec_srce[phases.getSet()[t_source]] = tmpvec.eigenVector;
505  }
506 
507  // Insert a ColorVector into spin index spin_source
508  // This only overwrites sections, so need to initialize first
509  LatticeFermion chi = zero;
510  CvToFerm(vec_srce, chi, spin_source);
511 
512  LatticeFermion quark_soln = zero;
513 
514  // Do the propagator inversion
515  SystemSolverResults_t res = (*PP)(quark_soln, chi);
516  ncg_had = res.n_count;
517 
518  // Insert into the fermion solution std::map object
519  KeyPropColorVec_t key;
520  key.t_source = t_source;
521  key.colorvec_src = colorvec_source;
522  key.spin_src = spin_source;
523 
524  prop_obj.insert(key, quark_soln);
525 
526  //
527  // The perambulator part
528  //
529  for(int colorvec_sink=0; colorvec_sink < num_vecs; ++colorvec_sink)
530  {
531  EVPair<LatticeColorVector> vec_sink; eigen_source.get(colorvec_sink,vec_sink);
532 
533  for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
534  {
535  multi1d<ComplexD> hsum(sumMulti(localInnerProduct(vec_sink.eigenVector, peekSpin(quark_soln, spin_sink)),
536  phases.getSet()));
537 
538  for(int t=0; t < hsum.size(); ++t)
539  {
540  buf(t,spin_sink).val.data().mat(colorvec_sink,colorvec_source) = hsum[t];
541  }
542  } // for spin_sink
543  } // for colorvec_sink
544  } // for colorvec_source
545 
546  for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
547  {
548  QDPIO::cout << "insert: spin_source= " << spin_source << " spin_sink= " << spin_sink << std::endl;
549  for(int t=0; t < phases.numSubsets(); ++t)
550  {
551  qdp_db.insert(buf(t,spin_sink).key, buf(t,spin_sink).val);
552  }
553  } // for spin_sink
554  } // for spin_source
555  } // for t_source
556 
557  // Finished writing
558  swatch.stop();
559  QDPIO::cout << "Propagators computed: time= "
560  << swatch.getTimeInSeconds()
561  << " secs" << std::endl;
562  }
563  catch (const std::string& e)
564  {
565  QDPIO::cout << name << ": caught exception around qprop: " << e << std::endl;
566  QDP_abort(1);
567  }
568 
569  push(xml_out,"Relaxation_Iterations");
570  write(xml_out, "ncg_had", ncg_had);
571  pop(xml_out);
572 
573  pop(xml_out); // prop_colorvec
574 
575  // Flush the object-std::map
576  prop_obj.flush();
577 
578  // Write the meta-data for this operator
579  {
580  XMLBufferWriter file_xml;
581 
582  push(file_xml, "PropColorVectors");
583  write(file_xml, "num_records", prop_obj.size());
584  write(file_xml, "Params", params.param);
585  write(file_xml, "Config_info", gauge_xml);
586  pop(file_xml);
587 
588  XMLBufferWriter record_xml;
589  push(record_xml, "PropColorVector");
590  write(record_xml, "num_records", prop_obj.size());
591  pop(record_xml);
592 
593  // Write the propagator xml info
594  TheNamedObjMap::Instance().get(params.named_obj.prop_id).setFileXML(file_xml);
595  TheNamedObjMap::Instance().get(params.named_obj.prop_id).setRecordXML(record_xml);
596  }
597 
598  snoop.stop();
599  QDPIO::cout << name << ": total time = "
600  << snoop.getTimeInSeconds()
601  << " secs" << std::endl;
602 
603  QDPIO::cout << name << ": ran successfully" << std::endl;
604 
605  END_CODE();
606  }
607 
608  }
609 
610 } // namespace Chroma
Inline measurement factory.
Class for counted reference semantics.
Definition: handle.h:33
void operator()(const unsigned long update_no, XMLWriter &xml_out)
Do the measurement.
void func(const unsigned long update_no, XMLWriter &xml_out)
Do the measurement.
Serializable value harness.
Definition: key_val_db.h:69
Fourier transform phase factor support.
Definition: sftmom.h:35
int numSubsets() const
Number of subsets - length in decay direction.
Definition: sftmom.h:63
const Set & getSet() const
The set to be used in sumMulti.
Definition: sftmom.h:57
static T & Instance()
Definition: singleton.h:432
Class structure for fermion actions.
Fermion action factories.
All Wilson-type fermion actions.
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.
GroupXML_t readXMLGroup(XMLReader &xml_in, const std::string &path, const std::string &type_name)
Read group and return as a std::string.
MODS_t & eigen_source
Eigenvectors.
Compute the perambulators.
Key for propagator colorstd::vector sources.
Key for propagator colorstd::vector matrix elements.
Key and values for DB.
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.
const std::string name
Name to be used.
multi1d< LatticeColorMatrix > P
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.
QDPSubTypeTrait< typename BinaryReturn< C1, C2, FnLocalInnerProduct >::Type_t >::Type_t localInnerProduct(const QDPSubType< T1, C1 > &l, const QDPType< T2, C2 > &r)
bool registerAll()
aggregate everything
bool registerAll()
Register all the factories.
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
static multi1d< LatticeColorMatrix > u
LinOpSysSolverMGProtoClover::Q Q
LatticeFermion tmp
Definition: mespbg5p_w.cc:36
push(xml_out,"Condensates")
LinOpSysSolverMGProtoClover::T T
multi1d< LatticeFermion > chi(Ncb)
multi1d< SubsetVectorWeight_t > getEigenValues(const MapObject< int, EVPair< LatticeColorVector > > &eigen_source, int num_vecs)
pop(xml_out)
void MesPlq(const multi1d< LatticeColorMatrixF3 > &u, multi2d< Double > &plane_plaq, Double &link)
Definition: mesplq.cc:83
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
::std::string string
Definition: gtest.h:1979
Print out basic info about this program.
Fourier transform phase factor support.
GroupXML_t invParam
Definition: qprop_io.h:84
GroupXML_t fermact
Definition: qprop_io.h:80
A Pair type.
Holds return info from SystemSolver call.
Definition: syssolver.h:17
Holds of vectors and weights.