inline_annih_prop_matelem_colorvec_w.cc

Go to the documentation of this file.
00001 // $Id: inline_annih_prop_matelem_colorvec_w.cc,v 3.4 2009/09/14 20:50:14 edwards Exp $
00002 /*! \file
00003  * \brief Compute the annihilation diagram propagator elements    M^-1 * multi1d<LatticeColorVector>
00004  *
00005  * Annihilation diagrams version of propagator calculation on a colorvector
00006  */
00007 
00008 #include "fermact.h"
00009 #include "meas/inline/hadron/inline_annih_prop_matelem_colorvec_w.h"
00010 #include "meas/inline/abs_inline_measurement_factory.h"
00011 #include "meas/glue/mesplq.h"
00012 #include "meas/sources/zN_src.h"
00013 #include "util/ferm/subset_vectors.h"
00014 #include "util/ferm/map_obj.h"
00015 #include "util/ferm/key_val_db.h"
00016 #include "util/ferm/key_prop_colorvec.h"
00017 #include "util/ferm/key_prop_matelem.h"
00018 #include "util/ferm/transf.h"
00019 #include "util/ft/sftmom.h"
00020 #include "util/info/proginfo.h"
00021 #include "actions/ferm/fermacts/fermact_factory_w.h"
00022 #include "actions/ferm/fermacts/fermacts_aggregate_w.h"
00023 #include "meas/inline/make_xml_file.h"
00024 
00025 #include "meas/inline/io/named_objmap.h"
00026 
00027 namespace Chroma 
00028 { 
00029   namespace InlineAnnihPropMatElemColorVecEnv 
00030   {
00031     //! Propagator input
00032     void read(XMLReader& xml, const string& path, Params::NamedObject_t& input)
00033     {
00034       XMLReader inputtop(xml, path);
00035 
00036       read(inputtop, "gauge_id", input.gauge_id);
00037       read(inputtop, "colorvec_id", input.colorvec_id);
00038       read(inputtop, "prop_op_file", input.prop_op_file);
00039     }
00040 
00041     //! Propagator output
00042     void write(XMLWriter& xml, const string& path, const Params::NamedObject_t& input)
00043     {
00044       push(xml, path);
00045 
00046       write(xml, "gauge_id", input.gauge_id);
00047       write(xml, "colorvec_id", input.colorvec_id);
00048       write(xml, "prop_op_file", input.prop_op_file);
00049 
00050       pop(xml);
00051     }
00052 
00053 
00054     //! Propagator input
00055     void read(XMLReader& xml, const string& path, Params::Param_t::Contract_t& input)
00056     {
00057       XMLReader inputtop(xml, path);
00058 
00059       read(inputtop, "num_vecs", input.num_vecs);
00060       read(inputtop, "t_sources_start", input.t_sources_start);
00061       read(inputtop, "dt", input.dt);
00062       read(inputtop, "decay_dir", input.decay_dir);
00063       read(inputtop, "N", input.N);
00064       read(inputtop, "ran_seed", input.ran_seed);
00065       read(inputtop, "mass_label", input.mass_label);
00066     }
00067 
00068     //! Propagator output
00069     void write(XMLWriter& xml, const string& path, const Params::Param_t::Contract_t& input)
00070     {
00071       push(xml, path);
00072 
00073       write(xml, "num_vecs", input.num_vecs);
00074       write(xml, "t_sources_start", input.t_sources_start);
00075       write(xml, "dt", input.dt);
00076       write(xml, "decay_dir", input.decay_dir);
00077       write(xml, "N", input.N);
00078       write(xml, "ran_seed", input.ran_seed);
00079       write(xml, "mass_label", input.mass_label);
00080 
00081       pop(xml);
00082     }
00083 
00084 
00085     //! Propagator input
00086     void read(XMLReader& xml, const string& path, Params::Param_t& input)
00087     {
00088       XMLReader inputtop(xml, path);
00089 
00090       read(inputtop, "Propagator", input.prop);
00091       read(inputtop, "Contractions", input.contract);
00092     }
00093 
00094     //! Propagator output
00095     void write(XMLWriter& xml, const string& path, const Params::Param_t& input)
00096     {
00097       push(xml, path);
00098 
00099       write(xml, "Propagator", input.prop);
00100       write(xml, "Contractions", input.contract);
00101 
00102       pop(xml);
00103     }
00104 
00105 
00106     //! Propagator input
00107     void read(XMLReader& xml, const string& path, Params& input)
00108     {
00109       Params tmp(xml, path);
00110       input = tmp;
00111     }
00112 
00113     //! Propagator output
00114     void write(XMLWriter& xml, const string& path, const Params& input)
00115     {
00116       push(xml, path);
00117     
00118       write(xml, "Param", input.param);
00119       write(xml, "NamedObject", input.named_obj);
00120 
00121       pop(xml);
00122     }
00123   } // namespace InlinePropColorVecEnv 
00124 
00125 
00126   namespace InlineAnnihPropMatElemColorVecEnv 
00127   {
00128     namespace
00129     {
00130       AbsInlineMeasurement* createMeasurement(XMLReader& xml_in, 
00131                                               const std::string& path) 
00132       {
00133         return new InlineMeas(Params(xml_in, path));
00134       }
00135 
00136       //! Local registration flag
00137       bool registered = false;
00138     }
00139       
00140     const std::string name = "ANNIH_PROP_MATELEM_COLORVEC";
00141 
00142     //! Register all the factories
00143     bool registerAll() 
00144     {
00145       bool success = true; 
00146       if (! registered)
00147       {
00148         success &= WilsonTypeFermActsEnv::registerAll();
00149         success &= TheInlineMeasurementFactory::Instance().registerObject(name, createMeasurement);
00150         registered = true;
00151       }
00152       return success;
00153     }
00154 
00155 
00156     //----------------------------------------------------------------------------
00157     // Param stuff
00158     Params::Params() { frequency = 0; }
00159 
00160     Params::Params(XMLReader& xml_in, const std::string& path) 
00161     {
00162       try 
00163       {
00164         XMLReader paramtop(xml_in, path);
00165 
00166         if (paramtop.count("Frequency") == 1)
00167           read(paramtop, "Frequency", frequency);
00168         else
00169           frequency = 1;
00170 
00171         // Parameters for source construction
00172         read(paramtop, "Param", param);
00173 
00174         // Read in the output propagator/source configuration info
00175         read(paramtop, "NamedObject", named_obj);
00176 
00177         // Possible alternate XML file pattern
00178         if (paramtop.count("xml_file") != 0) 
00179         {
00180           read(paramtop, "xml_file", xml_file);
00181         }
00182       }
00183       catch(const std::string& e) 
00184       {
00185         QDPIO::cerr << __func__ << ": Caught Exception reading XML: " << e << endl;
00186         QDP_abort(1);
00187       }
00188     }
00189 
00190 
00191 
00192     // Function call
00193     void 
00194     InlineMeas::operator()(unsigned long update_no,
00195                            XMLWriter& xml_out) 
00196     {
00197       // If xml file not empty, then use alternate
00198       if (params.xml_file != "")
00199       {
00200         string xml_file = makeXMLFileName(params.xml_file, update_no);
00201 
00202         push(xml_out, "AnnihPropMatElemColorVec");
00203         write(xml_out, "update_no", update_no);
00204         write(xml_out, "xml_file", xml_file);
00205         pop(xml_out);
00206 
00207         XMLFileWriter xml(xml_file);
00208         func(update_no, xml);
00209       }
00210       else
00211       {
00212         func(update_no, xml_out);
00213       }
00214     }
00215 
00216 
00217     // Real work done here
00218     void 
00219     InlineMeas::func(unsigned long update_no,
00220                      XMLWriter& xml_out) 
00221     {
00222       START_CODE();
00223 
00224       StopWatch snoop;
00225       snoop.reset();
00226       snoop.start();
00227 
00228       // Test and grab a reference to the gauge field
00229       multi1d<LatticeColorMatrix> u;
00230       XMLBufferWriter gauge_xml;
00231       try
00232       {
00233         u = TheNamedObjMap::Instance().getData< multi1d<LatticeColorMatrix> >(params.named_obj.gauge_id);
00234         TheNamedObjMap::Instance().get(params.named_obj.gauge_id).getRecordXML(gauge_xml);
00235       }
00236       catch( std::bad_cast ) 
00237       {
00238         QDPIO::cerr << name << ": caught dynamic cast error" << endl;
00239         QDP_abort(1);
00240       }
00241       catch (const string& e) 
00242       {
00243         QDPIO::cerr << name << ": map call failed: " << e << endl;
00244         QDP_abort(1);
00245       }
00246 
00247       push(xml_out, "AnnihPropMatElemColorVec");
00248       write(xml_out, "update_no", update_no);
00249 
00250       QDPIO::cout << name << ": annihilation propagator matrix element calculation" << endl;
00251 
00252       proginfo(xml_out);    // Print out basic program info
00253 
00254       // Write out the input
00255       write(xml_out, "Input", params);
00256 
00257       // Write out the config header
00258       write(xml_out, "Config_info", gauge_xml);
00259 
00260       push(xml_out, "Output_version");
00261       write(xml_out, "out_version", 1);
00262       pop(xml_out);
00263 
00264       // Calculate some gauge invariant observables just for info.
00265       MesPlq(xml_out, "Observables", u);
00266 
00267       //
00268       // Read in the source along with relevant information.
00269       // 
00270       XMLReader source_file_xml, source_record_xml;
00271 
00272       QDPIO::cout << "Snarf the source from a named buffer" << endl;
00273       try
00274       {
00275         TheNamedObjMap::Instance().getData< SubsetVectors<LatticeColorVector> >(params.named_obj.colorvec_id);
00276 
00277         // Snarf the source info. This is will throw if the colorvec_id is not there
00278         TheNamedObjMap::Instance().get(params.named_obj.colorvec_id).getFileXML(source_file_xml);
00279         TheNamedObjMap::Instance().get(params.named_obj.colorvec_id).getRecordXML(source_record_xml);
00280 
00281         // Write out the source header
00282         write(xml_out, "Source_file_info", source_file_xml);
00283         write(xml_out, "Source_record_info", source_record_xml);
00284       }    
00285       catch (std::bad_cast)
00286       {
00287         QDPIO::cerr << name << ": caught dynamic cast error" << endl;
00288         QDP_abort(1);
00289       }
00290       catch (const string& e) 
00291       {
00292         QDPIO::cerr << name << ": error extracting source_header: " << e << endl;
00293         QDP_abort(1);
00294       }
00295       const SubsetVectors<LatticeColorVector>& eigen_source = 
00296         TheNamedObjMap::Instance().getData< SubsetVectors<LatticeColorVector> >(params.named_obj.colorvec_id);
00297 
00298       QDPIO::cout << "Source successfully read and parsed" << endl;
00299 
00300       // Sanity check - write out the norm2 of the source in the Nd-1 direction
00301       // Use this for any possible verification
00302       {
00303         // Initialize the slow Fourier transform phases
00304         SftMom phases(0, true, Nd-1);
00305 
00306         multi1d< multi1d<Double> > source_corrs(eigen_source.getNumVectors());
00307         for(int m=0; m < source_corrs.size(); ++m)
00308         {
00309           source_corrs[m] = sumMulti(localNorm2(eigen_source.getEvectors()[m]), phases.getSet());
00310         }
00311 
00312         push(xml_out, "Source_correlators");
00313         write(xml_out, "source_corrs", source_corrs);
00314         pop(xml_out);
00315       }
00316 
00317       // Another sanity check
00318       if (params.param.contract.num_vecs > eigen_source.getNumVectors())
00319       {
00320         QDPIO::cerr << __func__ << ": num_vecs= " << params.param.contract.num_vecs
00321                     << " is greater than the number of available colorvectors= "
00322                     << eigen_source.getNumVectors() << endl;
00323         QDP_abort(1);
00324       }
00325 
00326       // Interval of sources should divide into lattice extent
00327       if ((QDP::Layout::lattSize()[params.param.contract.decay_dir] % params.param.contract.dt) != 0)
00328       {
00329         QDPIO::cerr << __func__ << ": dt= " << params.param.contract.dt
00330                     << " does not divide into lattice extent" << endl;
00331         QDP_abort(1);
00332       }
00333 
00334 
00335       //
00336       // DB storage
00337       //
00338       BinaryStoreDB< SerialDBKey<KeyPropElementalOperator_t>, SerialDBData<ValPropElementalOperator_t> > qdp_db;
00339 
00340       // Open the file, and write the meta-data and the binary for this operator
00341       if (! qdp_db.fileExists(params.named_obj.prop_op_file))
00342       {
00343         XMLBufferWriter file_xml;
00344 
00345         push(file_xml, "DBMetaData");
00346         write(file_xml, "id", string("propElemOp"));
00347         write(file_xml, "lattSize", QDP::Layout::lattSize());
00348         write(file_xml, "decay_dir", params.param.contract.decay_dir);
00349         proginfo(file_xml);    // Print out basic program info
00350         write(file_xml, "Params", params.param.contract);
00351         write(file_xml, "Config_info", gauge_xml);
00352         write(file_xml, "Weights", eigen_source.getEvalues());
00353         pop(file_xml);
00354 
00355         std::string file_str(file_xml.str());
00356         qdp_db.setMaxUserInfoLen(file_str.size());
00357 
00358         qdp_db.open(params.named_obj.prop_op_file, O_RDWR | O_CREAT, 0664);
00359 
00360         qdp_db.insertUserdata(file_str);
00361       }
00362       else
00363       {
00364         qdp_db.open(params.named_obj.prop_op_file, O_RDWR, 0664);
00365       }
00366 
00367 
00368 
00369       // Save current seed
00370       Seed ran_seed;
00371       QDP::RNG::savern(ran_seed);
00372         
00373       // Set the seed to desired value
00374       QDP::RNG::setrn(params.param.contract.ran_seed);
00375 
00376       // Total number of iterations
00377       int ncg_had = 0;
00378 
00379       //
00380       // Try the factories
00381       //
00382       try
00383       {
00384         StopWatch swatch;
00385         swatch.reset();
00386         QDPIO::cout << "Try the various factories" << endl;
00387 
00388         // Typedefs to save typing
00389         typedef LatticeFermion               T;
00390         typedef multi1d<LatticeColorMatrix>  P;
00391         typedef multi1d<LatticeColorMatrix>  Q;
00392 
00393         //
00394         // Initialize fermion action
00395         //
00396         std::istringstream  xml_s(params.param.prop.fermact.xml);
00397         XMLReader  fermacttop(xml_s);
00398         QDPIO::cout << "FermAct = " << params.param.prop.fermact.id << endl;
00399 
00400         // Generic Wilson-Type stuff
00401         Handle< FermionAction<T,P,Q> >
00402           S_f(TheFermionActionFactory::Instance().createObject(params.param.prop.fermact.id,
00403                                                                fermacttop,
00404                                                                params.param.prop.fermact.path));
00405 
00406         Handle< FermState<T,P,Q> > state(S_f->createState(u));
00407 
00408         Handle< SystemSolver<LatticeFermion> > PP = S_f->qprop(state,
00409                                                                params.param.prop.invParam);
00410       
00411         QDPIO::cout << "Suitable factory found: compute all the quark props" << endl;
00412         swatch.start();
00413 
00414         //
00415         // Loop over the source color and spin, creating the source
00416         // and calling the relevant propagator routines.
00417         //
00418         const int num_vecs        = params.param.contract.num_vecs;
00419         const int decay_dir       = params.param.contract.decay_dir;
00420         const int dt              = params.param.contract.dt;
00421         const int Lt              = QDP::Layout::lattSize()[decay_dir];
00422         const int num_sources     = Lt / params.param.contract.dt;
00423         const multi1d<int>& t_sources_start = params.param.contract.t_sources_start;
00424 
00425         // Initialize the slow Fourier transform phases
00426         SftMom phases(0, true, decay_dir);
00427 
00428         for(int tt=0; tt < t_sources_start.size(); ++tt)
00429         {
00430           int t_source_start = params.param.contract.t_sources_start[tt];
00431           // Construct all the solution vectors for only the first time source
00432           QDPIO::cout << "t_source_start = " << t_source_start << endl; 
00433 
00434           //
00435           // Have to hold temporarily the solution vectors
00436           //
00437           MapObject<KeyPropColorVec_t,LatticeFermion> map_obj;
00438 
00439           // The random numbers used for the stochastic combination of sources is held here
00440           MapObject<int,Complex> rng_map_obj;
00441 
00442           // Build up list of time sources
00443           multi1d<int> t_sources(num_sources);
00444 
00445           for(int nn=0; nn < t_sources.size(); ++nn)
00446           {
00447             int t = (t_source_start + dt*nn + Lt) % Lt;
00448             t_sources[nn] = t;
00449 
00450             rng_map_obj.insert(t, zN_rng(params.param.contract.N));
00451           }
00452 
00453 
00454           // All the loops
00455           for(int colorvec_source=0; colorvec_source < num_vecs; ++colorvec_source)
00456           {
00457             QDPIO::cout << "colorvec_source = " << colorvec_source << endl; 
00458 
00459             // Accumulate the source from several time-slices
00460             LatticeColorVector vec_srce = zero;
00461             for(int nn=0; nn < t_sources.size(); ++nn)
00462             {
00463               int t = t_sources[nn];
00464 
00465               // Pull out a time-slice of the color vector source, give it a random weight
00466               vec_srce[phases.getSet()[t]] = 
00467                 rng_map_obj[t] * eigen_source.getEvectors()[colorvec_source];
00468             }
00469         
00470             for(int spin_source=0; spin_source < Ns; ++spin_source)
00471             {
00472               QDPIO::cout << "spin_source = " << spin_source << endl; 
00473 
00474               // Insert a ColorVector into spin index spin_source
00475               // This only overwrites sections, so need to initialize first
00476               LatticeFermion chi = zero;
00477               CvToFerm(vec_srce, chi, spin_source);
00478 
00479               LatticeFermion quark_soln = zero;
00480 
00481               // Do the propagator inversion
00482               SystemSolverResults_t res = (*PP)(quark_soln, chi);
00483               ncg_had = res.n_count;
00484 
00485               KeyPropColorVec_t key;
00486               key.t_source     = t_source_start;
00487               key.colorvec_src = colorvec_source;
00488               key.spin_src     = spin_source;
00489             
00490               map_obj.insert(key, quark_soln);
00491             } // for spin_source
00492           } // for colorvec_source
00493 
00494 
00495           swatch.stop();
00496           QDPIO::cout << "Propagators computed: time= " 
00497                       << swatch.getTimeInSeconds() 
00498                       << " secs" << endl;
00499 
00500 
00501           //
00502           // All the loops
00503           //
00504           // NOTE: I pull the spin source and sink loops to the outside intentionally.
00505           // The idea is to create a colorvector index (2d) array. These are not
00506           // too big, but are big enough to make the IO efficient, and the DB efficient
00507           // on reading. For N=32 and Lt=128, the mats are 2MB.
00508           //
00509           swatch.reset();
00510           swatch.start();
00511           QDPIO::cout << "Extract matrix elements" << endl;
00512 
00513           for(int spin_source=0; spin_source < Ns; ++spin_source)
00514           {
00515             QDPIO::cout << "spin_source = " << spin_source << endl; 
00516           
00517             for(int spin_sink=0; spin_sink < Ns; ++spin_sink)
00518             {
00519               QDPIO::cout << "spin_sink = " << spin_sink << endl; 
00520 
00521               // Construct the keys and values. Have to hold the buffers here given that
00522               // the source and sink spin indices are pulled out as well as the time
00523               multi1d<KeyValPropElementalOperator_t> buf(t_sources.size());
00524               for(int nn=0; nn < t_sources.size(); ++nn)
00525               {
00526                 int t = t_sources[nn];
00527 
00528                 buf[nn].key.key().t_slice      = t;
00529                 buf[nn].key.key().t_source     = t;
00530                 buf[nn].key.key().spin_src     = spin_source;
00531                 buf[nn].key.key().spin_snk     = spin_sink;
00532                 buf[nn].key.key().mass_label   = params.param.contract.mass_label;
00533                 buf[nn].val.data().mat.resize(num_vecs,num_vecs);
00534               }
00535 
00536               for(int colorvec_source=0; colorvec_source < num_vecs; ++colorvec_source)
00537               {
00538                 KeyPropColorVec_t key;
00539                 key.t_source     = t_source_start;
00540                 key.colorvec_src = colorvec_source;
00541                 key.spin_src     = spin_source;
00542                   
00543                 LatticeColorVector vec_source(peekSpin(map_obj[key], spin_sink));
00544 
00545                 for(int colorvec_sink=0; colorvec_sink < num_vecs; ++colorvec_sink)
00546                 {
00547                   const LatticeColorVector& vec_sink = eigen_source.getEvectors()[colorvec_sink];
00548 
00549                   multi1d<ComplexD> hsum(sumMulti(localInnerProduct(vec_sink, vec_source), phases.getSet()));
00550                 
00551                   for(int nn=0; nn < t_sources.size(); ++nn)
00552                   {
00553                     int t = t_sources[nn];
00554 
00555                     buf[nn].val.data().mat(colorvec_sink,colorvec_source) = conj(rng_map_obj[t]) * hsum[t];
00556                   }
00557 
00558                 } // for colorvec_sink
00559               } // for colorvec_source
00560               
00561               QDPIO::cout << "insert: spin_source= " << spin_source << " spin_sink= " << spin_sink << endl; 
00562               for(int nn=0; nn < t_sources.size(); ++nn)
00563               {
00564                 int t = t_sources[nn];
00565                 qdp_db.insert(buf[nn].key, buf[nn].val);
00566               }
00567 
00568             } // for spin_sink
00569           } // for spin_source
00570 
00571           swatch.stop();
00572           QDPIO::cout << "Matrix elements computed: time= " 
00573                       << swatch.getTimeInSeconds() 
00574                       << " secs" << endl;
00575         } // tt
00576       }
00577       catch (const std::string& e) 
00578       {
00579         QDPIO::cout << name << ": caught exception around qprop: " << e << endl;
00580         QDP_abort(1);
00581       }
00582 
00583       push(xml_out,"Relaxation_Iterations");
00584       write(xml_out, "ncg_had", ncg_had);
00585       pop(xml_out);
00586 
00587       pop(xml_out);  // prop_colorvec
00588 
00589       // Reset the seed
00590       QDP::RNG::setrn(ran_seed);
00591 
00592       snoop.stop();
00593       QDPIO::cout << name << ": total time = "
00594                   << snoop.getTimeInSeconds() 
00595                   << " secs" << endl;
00596 
00597       QDPIO::cout << name << ": ran successfully" << endl;
00598 
00599       END_CODE();
00600     } 
00601 
00602   }
00603 
00604 } // namespace Chroma

Generated on Sun Nov 22 04:31:30 2009 for CHROMA by  doxygen 1.4.7