Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members | Related Pages

f_wilson.C

Go to the documentation of this file.
00001 #include<config.h>
00002 CPS_START_NAMESPACE
00008 //--------------------------------------------------------------------
00009 //  CVS keywords
00010 //
00011 //  $Source: /space/cvs/cps/cps++/src/util/lattice/f_wilson/f_wilson.C,v $
00012 //  $State: Exp $
00013 //
00014 //--------------------------------------------------------------------
00015 //------------------------------------------------------------------
00016 //
00017 // f_wilson.C
00018 //
00019 // Fwilson is derived from FwilsonTypes and is relevant to
00020 // wilson fermions
00021 //
00022 //------------------------------------------------------------------
00023 
00024 CPS_END_NAMESPACE
00025 #include <util/lattice.h>
00026 #include <util/dirac_op.h>
00027 #include <util/wilson.h>
00028 #include <util/verbose.h>
00029 #include <util/gjp.h>
00030 #include <util/error.h>
00031 #include <comms/scu.h>
00032 #include <comms/glb.h>
00033 
00034 #define BENCHMARK
00035 #ifdef BENCHMARK
00036 #include <util/qcdio.h>
00037 #include <sys/time.h>
00038 unsigned long WfmFlops;
00039 #ifndef timersub
00040 #define timersub(a, b, result)                                                \
00041   do {                                                                        \
00042     (result)->tv_sec = (a)->tv_sec - (b)->tv_sec;                             \
00043     (result)->tv_usec = (a)->tv_usec - (b)->tv_usec;                          \
00044     if ((result)->tv_usec < 0) {                                              \
00045       --(result)->tv_sec;                                                     \
00046       (result)->tv_usec += 1000000;                                           \
00047     }                                                                         \
00048   } while (0)
00049 #endif
00050 #endif
00051 
00052 CPS_START_NAMESPACE
00053 
00054 //------------------------------------------------------------------
00055 // Initialize static variables.
00056 //------------------------------------------------------------------
00057 
00058 //------------------------------------------------------------------
00059 // Constructor
00060 //------------------------------------------------------------------
00061 Fwilson::Fwilson()
00062 : FwilsonTypes()
00063 {
00064   cname = "Fwilson";
00065   char *fname = "Fwilson()";
00066   VRB.Func(cname,fname);
00067 
00068   //----------------------------------------------------------------
00069   // Check if anisotropy is present and exit since Fwilson has
00070   // not been tested for anisotropic lattices.
00071   //----------------------------------------------------------------
00072   if(GJP.XiBare() != 1 ||
00073      GJP.XiV()    != 1 ||
00074      GJP.XiVXi()  != 1   ){
00075     ERR.General(cname,fname,
00076     "XiBare=%g, XiV=%g, XiVXi=%g : Fwilson has not been tested with anisotropy\n",
00077                 GJP.XiBare(), GJP.XiV(), GJP.XiVXi());
00078   }
00079 
00080   //----------------------------------------------------------------
00081   // Do initializations before the wilson library can be used
00082   // Initialization involve memory allocation.
00083   //----------------------------------------------------------------
00084 
00085   //  static Wilson wilson_struct;
00086   //  f_dirac_op_init_ptr = &wilson_struct;
00087   //  wilson_init((Wilson *) f_dirac_op_init_ptr);
00088 
00089   static Wilson wilson_struct;
00090   f_dirac_op_init_ptr = &wilson_struct;
00091   wilson_init((Wilson*) f_dirac_op_init_ptr);
00092 }
00093 
00094 
00095 //------------------------------------------------------------------
00096 // Destructor
00097 //------------------------------------------------------------------
00098 Fwilson::~Fwilson()
00099 {
00100   char *fname = "~Fwilson()";
00101   VRB.Func(cname,fname);
00102 
00103   //----------------------------------------------------------------
00104   // Un-initialize the wilson library. Memory is set free here.
00105   //----------------------------------------------------------------
00106   wilson_end((Wilson *) f_dirac_op_init_ptr);
00107 }
00108 
00109 
00110 //------------------------------------------------------------------
00111 // FclassType Fclass():
00112 // It returns the type of fermion class.
00113 //------------------------------------------------------------------
00114 FclassType Fwilson::Fclass() const{
00115   return F_CLASS_WILSON;
00116 }
00117 
00118 int Fwilson::FsiteSize() const
00119 {
00120   return 2 * Colors() * SpinComponents();  
00121   // re/im * colors * spin_components
00122 }
00123 
00124 
00125 
00126 //------------------------------------------------------------------
00127 // int FchkbEvl() :
00128 // returns 1 => The fermion fields in the evolution
00129 //      or the CG that inverts the evolution matrix
00130 //      are defined on a single checkerboard (half the 
00131 //      lattice).
00132 //------------------------------------------------------------------
00133 int Fwilson::FchkbEvl() const
00134 {
00135   return 1;
00136 }
00137 
00138 
00139 //------------------------------------------------------------------
00140 // int FmatEvlInv(Vector *f_out, Vector *f_in, 
00141 //                CgArg *cg_arg, 
00142 //                Float *true_res,
00143 //                CnvFrmType cnv_frm = CNV_FRM_YES):
00144 // It calculates f_out where A * f_out = f_in and
00145 // A is the preconditioned fermion matrix that appears
00146 // in the HMC evolution (even/odd  preconditioning 
00147 // of [Dirac^dag Dirac]. The inversion is done
00148 // with the conjugate gradient. cg_arg is the structure
00149 // that contains all the control parameters, f_in is the
00150 // fermion field source vector, f_out should be set to be
00151 // the initial guess and on return is the solution.
00152 // f_in and f_out are defined on a checkerboard.
00153 // If true_res !=0 the value of the true residual is returned
00154 // in true_res.
00155 // *true_res = |src - MatPcDagMatPc * sol| / |src|
00156 // The function returns the total number of CG iterations.
00157 //------------------------------------------------------------------
00158 int Fwilson::FmatEvlInv(Vector *f_out, Vector *f_in, 
00159                         CgArg *cg_arg, 
00160                         Float *true_res,
00161                         CnvFrmType cnv_frm)
00162 {
00163   int iter;
00164   char *fname = "FmatEvlInv(CgArg*,V*,V*,F*,CnvFrmType)";
00165   VRB.Func(cname,fname);
00166 
00167   DiracOpWilson wilson(*this, f_out, f_in, cg_arg, cnv_frm);
00168 
00169   WfmFlops = 0;
00170   struct timeval t_start, t_stop;
00171   gettimeofday(&t_start,NULL);
00172   
00173   iter = wilson.InvCg(&(cg_arg->true_rsd));
00174   if (true_res) *true_res = cg_arg ->true_rsd;
00175 
00176   gettimeofday(&t_stop,NULL);
00177   timersub(&t_stop,&t_start,&t_start);
00178   double flops= (double)WfmFlops;
00179   double secs = t_start.tv_sec + 1.E-6 *t_start.tv_usec;
00180 //  printf("Wilson solve: %d iteratations %d flops %f Mflops per node\n",
00181 //       iter,WfmFlops,flops/(secs*1000000) );
00182 
00183   
00184   // Return the number of iterations
00185   return iter;
00186 }
00187 
00188 //------------------------------------------------------------------
00189 // int FmatEvlMInv(Vector *f_out, Vector *f_in, 
00190 //                Float shift[], int Nshift, 
00191 //                CgArg **cg_arg, Float *true_res,
00192 //                CnvFrmType cnv_frm = CNV_FRM_YES):
00193 // It calculates f_out where (A + shift)* f_out = f_in and
00194 // A is the fermion matrix that appears in the HMC 
00195 // evolution ([Dirac^dag Dirac]) and shift is a real shift of the 
00196 // fermion matrix, with Nshift such shifts. The inversion is done 
00197 // with the multishift conjugate gradient. cg_arg is the structure
00198 // that contains all the control parameters, f_in is the
00199 // fermion field source vector, f_out is the array of solution 
00200 // vectors, f_in and f_out are defined on a checkerboard.
00201 // The function returns the total number of CG iterations.
00202 //------------------------------------------------------------------
00203 int Fwilson::FmatEvlMInv(Vector **f_out, Vector *f_in, Float *shift, 
00204                          int Nshift, int isz, CgArg **cg_arg,
00205                          CnvFrmType cnv_frm, MultiShiftSolveType type, 
00206                          Float *alpha, Vector **f_out_d)
00207 {
00208   char *fname = "FmatMInv(V*, V*, .....)";
00209   VRB.Func(cname,fname);
00210 
00211   int f_size = GJP.VolNodeSites() * FsiteSize() / (FchkbEvl()+1);
00212   Float dot = f_in -> NormSqGlbSum4D(f_size);
00213 
00214   Float *RsdCG = new Float[Nshift];
00215   for (int s=0; s<Nshift; s++) RsdCG[s] = cg_arg[s]->stop_rsd;
00216 
00217   //Fake the constructor
00218   DiracOpWilson wilson(*this, f_out[0], f_in, cg_arg[0], cnv_frm);
00219 
00220   int return_value = wilson.MInvCG(f_out,f_in,dot,shift,Nshift,isz,RsdCG,type,alpha);  
00221 
00222   for (int s=0; s<Nshift; s++) cg_arg[s]->true_rsd = RsdCG[s];
00223   delete[] RsdCG;
00224   return return_value;
00225 }
00226 
00227 //------------------------------------------------------------------
00228 // Lattice class api to the chronological inverter
00229 //------------------------------------------------------------------
00230 void Fwilson::FminResExt(Vector *sol, Vector *source, Vector **sol_old, 
00231                          Vector **vm, int degree, CgArg *cg_arg, CnvFrmType cnv_frm)
00232 {
00233 
00234   char *fname = "FminResExt(V*, V*, V**, V**, int, CgArg *, CnvFrmType)";
00235   VRB.Func(cname,fname);
00236   
00237   DiracOpWilson wilson(*this, sol, source, cg_arg, cnv_frm);
00238   wilson.MinResExt(sol,source,sol_old,vm,degree);
00239   
00240 }
00241 
00242 //------------------------------------------------------------------
00243 // int FmatInv(Vector *f_out, Vector *f_in, 
00244 //             CgArg *cg_arg, 
00245 //             Float *true_res,
00246 //             CnvFrmType cnv_frm = CNV_FRM_YES,
00247 //             PreserveType prs_f_in = PRESERVE_YES):
00248 // It calculates f_out where A * f_out = f_in and
00249 // A is the fermion matrix (Dirac operator). The inversion
00250 // is done with the conjugate gradient. cg_arg is the 
00251 // structure that contains all the control parameters, f_in 
00252 // is the fermion field source vector, f_out should be set 
00253 // to be the initial guess and on return is the solution.
00254 // f_in and f_out are defined on the whole lattice.
00255 // If true_res !=0 the value of the true residual is returned
00256 // in true_res.
00257 // *true_res = |src - MatPcDagMatPc * sol| / |src|
00258 // cnv_frm is used to specify if f_in should be converted 
00259 // from canonical to fermion order and f_out from fermion 
00260 // to canonical. 
00261 // prs_f_in is used to specify if the source
00262 // f_in should be preserved or not. If not the memory usage
00263 // is less by half the size of a fermion vector.
00264 // The function returns the total number of CG iterations.
00265 //------------------------------------------------------------------
00266 int Fwilson::FmatInv(Vector *f_out, Vector *f_in, 
00267                      CgArg *cg_arg, 
00268                      Float *true_res,
00269                      CnvFrmType cnv_frm,
00270                      PreserveType prs_f_in)
00271 {
00272   int iter;
00273   char *fname = "FmatInv(CgArg*,V*,V*,F*,CnvFrmType)";
00274   VRB.Func(cname,fname);
00275 
00276   DiracOpWilson wilson(*this, f_out, f_in, cg_arg, cnv_frm);
00277   
00278   iter = wilson.MatInv(true_res, prs_f_in);
00279   
00280   // Return the number of iterations
00281   return iter;
00282 }
00283 
00284 
00285 
00286 //------------------------------------------------------------------
00287 // int FeigSolv(Vector **f_eigenv, Float *lambda, int *valid_eig,
00288 //              EigArg *eig_arg, 
00289 //              CnvFrmType cnv_frm = CNV_FRM_YES):
00290 // It solve  A * f_eigenv = lambda * f_eigenv where
00291 // A is the fermion matrix (Dirac operator). The solution
00292 // is done with the Ritz algorithm. eig_arg is the 
00293 // structure that contains all the control parameters, f_eigenv
00294 // is the fermion field eigenvectors, lambda are the
00295 // returned eigenvalues.
00296 // f_eigenv is defined on the whole lattice.
00297 // The function returns the total number of Ritz iterations.
00298 //------------------------------------------------------------------
00299 int Fwilson::FeigSolv(Vector **f_eigenv, Float *lambda,
00300                       Float *chirality, int *valid_eig,
00301                       Float **hsum,
00302                       EigArg *eig_arg, 
00303                       CnvFrmType cnv_frm)
00304 {
00305   int iter;
00306   char *fname = "FeigSolv(EigArg*,V*,F*,CnvFrmType)";
00307   VRB.Func(cname,fname);
00308   CgArg cg_arg;
00309   cg_arg.mass = eig_arg->mass;
00310   cg_arg.RitzMatOper = eig_arg->RitzMatOper;
00311   int N_eig = eig_arg->N_eig;
00312   int i;
00313 
00314   //=========================
00315   // convert fermion field
00316   //=========================
00317 
00318 
00319   for(i=0; i < N_eig; ++i)
00320     Fconvert(f_eigenv[i], WILSON, CANONICAL);
00321 
00322 
00323 
00324   //------------------------------------------------------------------
00325   //  we want both the eigenvalues of D_{hermitian} and
00326   //  D^{+}D.  To not change the arguments passed to RitzEig,
00327   //  we pass a float pointer which points to 2 * N_eig values
00328   //  and return both lambda and lambda^2 from RitzEig
00329   //------------------------------------------------------------------
00330 
00331   Float * lambda2 = (Float * ) smalloc (N_eig*2*sizeof(Float));
00332   if ( lambda2 == 0 ) ERR.Pointer(cname,fname, "lambda2");
00333   
00334   {
00335     DiracOpWilson wilson(*this, (Vector*) 0 , (Vector*) 0, &cg_arg, CNV_FRM_NO);
00336     iter = wilson.RitzEig(f_eigenv, lambda2, valid_eig, eig_arg);
00337   }
00338 
00339 
00340   for(i=0; i < N_eig; ++i)
00341     {
00342       Fconvert(f_eigenv[i], CANONICAL, WILSON);
00343     }
00344 
00345   /*
00346     the call to RitzEig returns a negative number if either the KS or CG maxes
00347     out, we wish to cope with this in alg_eig, so "pass it up". Clean up the
00348     storage order first in case we still want to use the eigenvectors as a
00349     guess.
00350   */
00351   if ( iter < 0 ) { return iter ; }
00352 
00353 
00354   // Compute chirality
00355   int f_size = (GJP.VolNodeSites() * FsiteSize());
00356 
00357   Vector* v1 = (Vector *)smalloc(f_size*sizeof(Float));
00358   if (v1 == 0)
00359     ERR.Pointer(cname, fname, "v1");
00360   VRB.Smalloc(cname, fname, "v1", v1, f_size*sizeof(Float));
00361 
00362   for(i=0; i < N_eig; ++i)
00363   {
00364     Gamma5(v1, f_eigenv[i], GJP.VolNodeSites());
00365     chirality[i] = f_eigenv[i]->ReDotProductGlbSum4D(v1, f_size);
00366   }
00367 
00368   VRB.Sfree(cname, fname, "v1", v1);
00369   sfree(v1);
00370 
00371 
00372   // rescale wilson eigenvalues to the convention  m + Dslash(U)
00373   Float factor = 4.0 + eig_arg->mass;
00374     
00375   
00376   FILE* fp=Fopen(eig_arg->fname,"a");
00377   for(i=0; i<N_eig; ++i)
00378     {
00379       lambda2[i] *= factor;                      //rescale eigenvalue
00380       lambda2[N_eig + i] *= ( factor * factor ); //rescale squared evalue
00381       lambda[i]=lambda2[i];                      //copy back
00382       
00383       //print out eigenvalue, eigenvalue^2, chirality 
00384       Fprintf(fp,"%d %g %g %g %d\n",i,
00385               (float)lambda2[i],
00386               (float)lambda2[N_eig + i],
00387               (float)chirality[i],valid_eig[i]);
00388     }
00389   Fclose(fp);
00390   sfree(lambda2); 
00391 
00392 
00393   // Slice-sum the eigenvector density to make a 1D vector
00394   if (eig_arg->print_hsum)
00395     for(i=0; i < N_eig; ++i)
00396       f_eigenv[i]->NormSqArraySliceSum(hsum[i], FsiteSize(), eig_arg->hsum_dir);
00397 
00398 
00399   // The remaining part in QCDSP version are all about "downloading
00400   // eigenvectors", supposedly not applicable here.
00401 
00402   // Return the number of iterations
00403   return iter;
00404 }
00405 
00406 
00407 //------------------------------------------------------------------
00408 // SetPhi(Vector *phi, Vector *frm1, Vector *frm2, Float mass,
00409 //        DagType dag):
00410 // It sets the pseudofermion field phi from frm1, frm2.
00411 // Note that frm2 is not used.
00412 // Modified - now returns the (trivial) value of the action
00413 //------------------------------------------------------------------
00414 Float Fwilson::SetPhi(Vector *phi, Vector *frm1, Vector *frm2,
00415                       Float mass, DagType dag){
00416   char *fname = "SetPhi(V*,V*,V*,F)";
00417   VRB.Func(cname,fname);
00418   CgArg cg_arg;
00419   cg_arg.mass = mass;
00420 
00421   if (phi == 0)
00422     ERR.Pointer(cname,fname,"phi") ;
00423 
00424   if (frm1 == 0)
00425     ERR.Pointer(cname,fname,"frm1") ;
00426 
00427   DiracOpWilson wilson(*this, frm1, frm2, &cg_arg, CNV_FRM_NO) ;
00428   
00429   if (dag == DAG_YES) wilson.MatPcDag(phi, frm1) ;
00430   else wilson.MatPc(phi, frm1) ;
00431 
00432   return FhamiltonNode(frm1, frm1);
00433 }
00434 
00435 
00436 //------------------------------------------------------------------
00437 // EvolveMomFforce(Matrix *mom, Vector *frm, Float mass, 
00438 //                 Float dt):
00439 // It evolves the canonical momentum mom by dt
00440 // using the fermion force.
00441 //------------------------------------------------------------------
00442 ForceArg Fwilson::EvolveMomFforce(Matrix *mom, Vector *chi, 
00443                               Float mass, Float dt)
00444 {
00445   char *fname = "EvolveMomFforce(M*,V*,F,F,F)";
00446   VRB.Func(cname,fname);
00447 
00448   Matrix *gauge = GaugeField() ;
00449 
00450   if (Colors() != 3)
00451     ERR.General(cname,fname,"Wrong nbr of colors.") ;
00452 
00453   if (SpinComponents() != 4)
00454     ERR.General(cname,fname,"Wrong nbr of spin comp.") ;
00455 
00456   if (mom == 0)
00457     ERR.Pointer(cname,fname,"mom") ;
00458 
00459   if (chi == 0)
00460     ERR.Pointer(cname,fname,"chi") ;
00461 
00462 //------------------------------------------------------------------
00463 // allocate space for two CANONICAL fermion fields.
00464 //------------------------------------------------------------------
00465   int f_size = FsiteSize() * GJP.VolNodeSites() ;
00466 
00467   char *str_v1 = "v1" ;
00468   Vector *v1 = (Vector *)smalloc(f_size*sizeof(Float)) ;
00469   if (v1 == 0)
00470     ERR.Pointer(cname, fname, str_v1) ;
00471   VRB.Smalloc(cname, fname, str_v1, v1, f_size*sizeof(Float)) ;
00472 
00473   char *str_v2 = "v2" ;
00474   Vector *v2 = (Vector *)smalloc(f_size*sizeof(Float)) ;
00475   if (v2 == 0)
00476     ERR.Pointer(cname, fname, str_v2) ;
00477   VRB.Smalloc(cname, fname, str_v2, v2, f_size*sizeof(Float)) ;
00478 
00479 //------------------------------------------------------------------
00480 // allocate space for two CANONICAL fermion field on a site.
00481 //------------------------------------------------------------------
00482 
00483   char *str_site_v1 = "site_v1";
00484   Float *site_v1 = (Float *)smalloc(FsiteSize()*sizeof(Float));
00485   if (site_v1 == 0)
00486     ERR.Pointer(cname, fname, str_site_v1) ;
00487   VRB.Smalloc(cname, fname, str_site_v1, site_v1,
00488     FsiteSize()*sizeof(Float)) ;
00489 
00490   char *str_site_v2 = "site_v2";
00491   Float *site_v2 = (Float *)smalloc(FsiteSize()*sizeof(Float));
00492   if (site_v2 == 0)
00493     ERR.Pointer(cname, fname, str_site_v2) ;
00494   VRB.Smalloc(cname, fname, str_site_v2, site_v2,
00495     FsiteSize()*sizeof(Float)) ;
00496 
00497   {
00498     CgArg cg_arg ;
00499     cg_arg.mass = mass ;
00500 
00501     DiracOpWilson wilson(*this, v1, v2, &cg_arg, CNV_FRM_YES) ;
00502     wilson.CalcHmdForceVecs(chi) ;
00503   }
00504 
00505   int x, y, z, t, lx, ly, lz, lt ;
00506 
00507   lx = GJP.XnodeSites() ;
00508   ly = GJP.YnodeSites() ;
00509   lz = GJP.ZnodeSites() ;
00510   lt = GJP.TnodeSites() ;
00511 
00512 //------------------------------------------------------------------
00513 // start by summing first over direction (mu) and then over site
00514 // to allow SCU transfers to happen face-by-face in the outermost
00515 // loop.
00516 //------------------------------------------------------------------
00517 
00518   int mu ;
00519 
00520   Matrix tmp, f ;
00521 
00522   Float L1 = 0.0;
00523   Float L2 = 0.0;
00524   Float Linf = 0.0;
00525 
00526   for (mu=0; mu<4; mu++) {
00527     for (t=0; t<lt; t++)
00528     for (z=0; z<lz; z++)
00529     for (y=0; y<ly; y++)
00530     for (x=0; x<lx; x++) {
00531       int gauge_offset = x+lx*(y+ly*(z+lz*t)) ;
00532       int vec_offset = FsiteSize()*gauge_offset ;
00533       gauge_offset = mu+4*gauge_offset ;
00534 
00535       Float *v1_plus_mu ;
00536       Float *v2_plus_mu ;
00537       int vec_plus_mu_offset = FsiteSize() ;
00538 
00539       Float coeff = -2.0 * dt ;
00540 
00541       switch (mu) {
00542         case 0 :
00543           vec_plus_mu_offset *= (x+1)%lx+lx*(y+ly*(z+lz*t)) ;
00544           if ((x+1) == lx) {
00545             getPlusData( (IFloat *)site_v1,
00546                          (IFloat *)v1+vec_plus_mu_offset, FsiteSize(), mu) ;
00547             getPlusData( (IFloat *)site_v2,
00548                          (IFloat *)v2+vec_plus_mu_offset, FsiteSize(), mu) ;
00549             v1_plus_mu = site_v1 ;                        
00550             v2_plus_mu = site_v2 ;                        
00551             if (GJP.XnodeBc()==BND_CND_APRD) coeff = -coeff ;
00552           } else {
00553             v1_plus_mu = (Float *)v1+vec_plus_mu_offset ;
00554             v2_plus_mu = (Float *)v2+vec_plus_mu_offset ;
00555           }
00556           break ;
00557         case 1 :
00558           vec_plus_mu_offset *= x+lx*((y+1)%ly+ly*(z+lz*t)) ;
00559           if ((y+1) == ly) {
00560             getPlusData( (IFloat *)site_v1,
00561                          (IFloat *)v1+vec_plus_mu_offset, FsiteSize(), mu) ;
00562             getPlusData( (IFloat *)site_v2,
00563                          (IFloat *)v2+vec_plus_mu_offset, FsiteSize(), mu) ;
00564             v1_plus_mu = site_v1 ;                        
00565             v2_plus_mu = site_v2 ;                        
00566             if (GJP.YnodeBc()==BND_CND_APRD) coeff = -coeff ;
00567           } else {
00568             v1_plus_mu = (Float *)v1+vec_plus_mu_offset ;
00569             v2_plus_mu = (Float *)v2+vec_plus_mu_offset ;
00570           }
00571           break ;
00572         case 2 :
00573           vec_plus_mu_offset *= x+lx*(y+ly*((z+1)%lz+lz*t)) ;
00574           if ((z+1) == lz) {
00575             getPlusData( (IFloat *)site_v1,
00576                          (IFloat *)v1+vec_plus_mu_offset, FsiteSize(), mu) ;
00577             getPlusData( (IFloat *)site_v2,
00578                          (IFloat *)v2+vec_plus_mu_offset, FsiteSize(), mu) ;
00579             v1_plus_mu = site_v1 ;
00580             v2_plus_mu = site_v2 ;
00581             if (GJP.ZnodeBc()==BND_CND_APRD) coeff = -coeff ;
00582           } else {
00583             v1_plus_mu = (Float *)v1+vec_plus_mu_offset ;
00584             v2_plus_mu = (Float *)v2+vec_plus_mu_offset ;
00585           }
00586           break ;
00587         case 3 :
00588           vec_plus_mu_offset *= x+lx*(y+ly*(z+lz*((t+1)%lt))) ;
00589           if ((t+1) == lt) {
00590             getPlusData( (IFloat *)site_v1,
00591                          (IFloat *)v1+vec_plus_mu_offset, FsiteSize(), mu) ;
00592             getPlusData( (IFloat *)site_v2,
00593                          (IFloat *)v2+vec_plus_mu_offset, FsiteSize(), mu) ;
00594             v1_plus_mu = site_v1 ;
00595             v2_plus_mu = site_v2 ;
00596             if (GJP.TnodeBc()==BND_CND_APRD) coeff = -coeff ;
00597           } else {
00598             v1_plus_mu = (Float *)v1+vec_plus_mu_offset ;
00599             v2_plus_mu = (Float *)v2+vec_plus_mu_offset ;
00600           }
00601       } // end switch mu
00602 
00603       sproj_tr[mu](   (IFloat *)&tmp,
00604                       (IFloat *)v1_plus_mu,
00605                       (IFloat *)v2+vec_offset, 1, 0, 0);
00606 
00607       sproj_tr[mu+4]( (IFloat *)&f,
00608                       (IFloat *)v2_plus_mu,
00609                       (IFloat *)v1+vec_offset, 1, 0, 0);
00610 
00611       tmp += f ;
00612 
00613       f.DotMEqual(*(gauge+gauge_offset), tmp) ;
00614 
00615       tmp.Dagger(f) ;
00616 
00617       f.TrLessAntiHermMatrix(tmp) ;
00618 
00619       f *= coeff ;
00620 
00621       *(mom+gauge_offset) += f ;
00622       Float norm = f.norm();
00623       Float tmp = sqrt(norm);
00624       L1 += tmp;
00625       L2 += norm;
00626       Linf = (tmp>Linf ? tmp : Linf);
00627     }
00628   }
00629 
00630 //------------------------------------------------------------------
00631 // deallocate space for two CANONICAL fermion fields on a site.
00632 //------------------------------------------------------------------
00633   VRB.Sfree(cname, fname, str_site_v2, site_v2) ;
00634   sfree(site_v2) ;
00635 
00636   VRB.Sfree(cname, fname, str_site_v1, site_v1) ;
00637   sfree(site_v1) ;
00638 
00639 //------------------------------------------------------------------
00640 // deallocate space for two CANONICAL fermion fields.
00641 //------------------------------------------------------------------
00642   VRB.Sfree(cname, fname, str_v2, v2) ;
00643   sfree(v2) ;
00644 
00645   VRB.Sfree(cname, fname, str_v1, v1) ;
00646   sfree(v1) ;
00647 
00648   glb_sum(&L1);
00649   glb_sum(&L2);
00650   glb_max(&Linf);
00651 
00652   L1 /= 4.0*GJP.VolSites();
00653   L2 /= 4.0*GJP.VolSites();
00654 
00655   return ForceArg(L1, sqrt(L2), Linf);
00656 }
00657 
00658 ForceArg Fwilson::RHMC_EvolveMomFforce(Matrix *mom, Vector **sol, int degree,
00659                                     int isz, Float *alpha, Float mass, 
00660                                     Float dt, Vector **sol_d, 
00661                                     ForceMeasure force_measure) {
00662   char *fname = "RHMC_EvolveMomFforce";
00663   char *force_label;
00664 
00665   ForceArg Fdt;
00666   Float L1 = 0.0;
00667   Float L2 = 0.0;
00668   Float Linf = 0.0;
00669 
00670   int g_size = GJP.VolNodeSites() * GsiteSize();
00671 
00672   Matrix *mom_tmp;
00673 
00674   if (force_measure == FORCE_MEASURE_YES) {
00675     mom_tmp = (Matrix*)smalloc(g_size*sizeof(Float),cname, fname, "mom_tmp");
00676     ((Vector*)mom_tmp) -> VecZero(g_size);
00677     force_label = new char[100];
00678   } else {
00679     mom_tmp = mom;
00680   }
00681 
00682   for (int i=0; i<degree; i++) {
00683     ForceArg Fdt = EvolveMomFforce(mom_tmp,sol[i],mass,dt*alpha[i]);
00684     if (force_measure == FORCE_MEASURE_YES) {
00685       sprintf(force_label, "Rational, mass = %e, pole = %d:", mass, i+isz);
00686       Fdt.print(dt, force_label);
00687     }
00688   }
00689 
00690   // If measuring the force, need to measure and then sum to mom
00691   if (force_measure == FORCE_MEASURE_YES) {
00692     for (int i=0; i<g_size/18; i++) {
00693       Float norm = (mom_tmp+i)->norm();
00694       Float tmp = sqrt(norm);
00695       L1 += tmp;
00696       L2 += norm;
00697       Linf = (tmp>Linf ? tmp : Linf);
00698     }
00699     glb_sum(&L1);
00700     glb_sum(&L2);
00701     glb_max(&Linf);
00702 
00703     L1 /= 4.0*GJP.VolSites();
00704     L2 /= 4.0*GJP.VolSites();
00705 
00706     fTimesV1PlusV2((IFloat*)mom, 1.0, (IFloat*)mom_tmp, (IFloat*)mom, g_size);
00707 
00708     delete[] force_label;
00709     sfree(mom_tmp, cname, fname, "mom_tmp");
00710   }
00711 
00712   return ForceArg(L1, sqrt(L2), Linf);
00713 }
00714 
00715 //------------------------------------------------------------------
00716 // Float BhamiltonNode(Vector *boson, Float mass):
00717 // The boson Hamiltonian of the node sublattice.
00718 //------------------------------------------------------------------
00719 Float Fwilson::BhamiltonNode(Vector *boson, Float mass){
00720   char *fname = "BhamiltonNode(V*,F)";
00721   VRB.Func(cname,fname);
00722   CgArg cg_arg;
00723   cg_arg.mass = mass;
00724 
00725   if (boson == 0)
00726     ERR.Pointer(cname,fname,"boson");
00727 
00728   int f_size = (GJP.VolNodeSites() * FsiteSize()) >> 1 ;
00729 
00730   Vector *bsn_tmp = (Vector *)
00731     smalloc(f_size*sizeof(Float));
00732 
00733   char *str_tmp = "bsn_tmp" ;
00734 
00735   if (bsn_tmp == 0)
00736     ERR.Pointer(cname,fname,str_tmp) ;
00737 
00738   VRB.Smalloc(cname,fname,str_tmp,bsn_tmp,f_size*sizeof(Float));
00739 
00740   DiracOpWilson wilson(*this, boson, bsn_tmp, &cg_arg, CNV_FRM_NO) ;
00741 
00742   wilson.MatPc(bsn_tmp,boson);
00743 
00744   Float ret_val = bsn_tmp->NormSqNode(f_size) ;
00745 
00746   VRB.Sfree(cname,fname,str_tmp,bsn_tmp);
00747 
00748   sfree(bsn_tmp) ;
00749 
00750   return ret_val;
00751 }
00752 
00753 ForceArg Fwilson::EvolveMomFforce(Matrix *mom, Vector *phi, Vector *eta,
00754                       Float mass, Float dt) {
00755   char *fname = "EvolveMomFforce(M*,V*,V*,F,F)";
00756   ERR.General(cname,fname,"Not Implemented\n");
00757   return ForceArg(0.0,0.0,0.0);
00758 }
00759 
00760 CPS_END_NAMESPACE

Generated on Sat Oct 10 14:11:31 2009 for Columbia Physics System by  doxygen 1.3.9.1