mirror of
				https://github.com/paboyle/Grid.git
				synced 2025-11-04 05:54:32 +00:00 
			
		
		
		
	Stout smearing compiles (untested)
This commit is contained in:
		@@ -5,146 +5,120 @@
 | 
			
		||||
#ifndef APE_SMEAR_
 | 
			
		||||
#define APE_SMEAR_
 | 
			
		||||
 | 
			
		||||
/*!  @brief APE type smearing of link variables. */
 | 
			
		||||
namespace Grid {
 | 
			
		||||
  namespace QCD {
 | 
			
		||||
 | 
			
		||||
template <class Gimpl> 
 | 
			
		||||
class Smear_APE: public Smear<Gimpl>{
 | 
			
		||||
private:
 | 
			
		||||
  const std::vector<double> rho;/*!< Array of weights */
 | 
			
		||||
 | 
			
		||||
  //This member must be private - we do not want to control from outside 
 | 
			
		||||
  std::vector<double> set_rho(const double)const {
 | 
			
		||||
    std::vector<double> res;
 | 
			
		||||
    /*!  @brief APE type smearing of link variables. */
 | 
			
		||||
    template <class Gimpl> 
 | 
			
		||||
    class Smear_APE: public Smear<Gimpl>{
 | 
			
		||||
    private:
 | 
			
		||||
      const std::vector<double> rho;/*!< Array of weights */
 | 
			
		||||
 | 
			
		||||
      //This member must be private - we do not want to control from outside 
 | 
			
		||||
      std::vector<double> set_rho(const double common_rho)const {
 | 
			
		||||
	std::vector<double> res;
 | 
			
		||||
    
 | 
			
		||||
    for(int mn=0; mn<Nd*Nd; ++mn) res.push_back(common_rho);
 | 
			
		||||
    for(int mu=0; mu<Nd; ++mu) res[mu + mu*Nd] = 0.0;
 | 
			
		||||
    return res;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
public:
 | 
			
		||||
  INHERIT_GIMPL_TYPES(Gimpl)
 | 
			
		||||
 | 
			
		||||
  Smear_APE(const std::vector<double>& rho_):rho(rho_){}
 | 
			
		||||
  Smear_APE(double rho_val):rho(set_rho(rho_val)){}
 | 
			
		||||
  Smear_APE():rho(set_rho(1.0)){}
 | 
			
		||||
  ~Smear_APE(){}
 | 
			
		||||
 | 
			
		||||
  void smear(GaugeField& u_smr, const GaugeField& U)const{
 | 
			
		||||
    double d_rho;
 | 
			
		||||
    GaugeLinkField Cup, tmp_stpl;
 | 
			
		||||
    WilsonLoops<Gimpl> WL;
 | 
			
		||||
    u_smr = zero;
 | 
			
		||||
 | 
			
		||||
    for(int mu=0; mu<Nd; ++mu){
 | 
			
		||||
      Cup = zero;
 | 
			
		||||
      for(int nu=0; nu<Nd; ++nu){
 | 
			
		||||
	d_rho = rho[mu + Nd * nu];
 | 
			
		||||
	WL.Staple(tmp_stpl, U, mu, nu);
 | 
			
		||||
	Cup += tmp_stpl*d_rho;
 | 
			
		||||
	for(int mn=0; mn<Nd*Nd; ++mn) res.push_back(common_rho);
 | 
			
		||||
	for(int mu=0; mu<Nd; ++mu) res[mu + mu*Nd] = 0.0;
 | 
			
		||||
	return res;
 | 
			
		||||
      }
 | 
			
		||||
      pokeLorentz(u_smr, Cup, mu);
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void derivative(GaugeField& SigmaTerm,
 | 
			
		||||
		  const GaugeField& iLambda,
 | 
			
		||||
		  const GaugeField& U)const{
 | 
			
		||||
    public:
 | 
			
		||||
      // Defines the gauge field types
 | 
			
		||||
      INHERIT_GIMPL_TYPES(Gimpl)
 | 
			
		||||
 | 
			
		||||
    /*
 | 
			
		||||
  
 | 
			
		||||
      // Constructors and destructors
 | 
			
		||||
      Smear_APE(const std::vector<double>& rho_):rho(rho_){}
 | 
			
		||||
      Smear_APE(double rho_val):rho(set_rho(rho_val)){}
 | 
			
		||||
      Smear_APE():rho(set_rho(1.0)){}
 | 
			
		||||
      ~Smear_APE(){}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    WilsonLoops<Gimpl> WL;
 | 
			
		||||
    GaugeLinkField staple, u_tmp, iLambda_mu, iLambda_nu;
 | 
			
		||||
    GaugeLinkField U_mu, U_nu;
 | 
			
		||||
    GaugeLinkField sh_field ;
 | 
			
		||||
    GaugeLinkField temp_Sigma;
 | 
			
		||||
      void smear(GaugeField& u_smr, const GaugeField& U)const{
 | 
			
		||||
	GridBase *grid = U._grid;
 | 
			
		||||
	double d_rho;
 | 
			
		||||
	GaugeLinkField Cup(grid), tmp_stpl(grid);
 | 
			
		||||
	WilsonLoops<Gimpl> WL;
 | 
			
		||||
	u_smr = zero; // probably unecessary
 | 
			
		||||
 | 
			
		||||
    SU<N>::Matrix temp_mat, temp_mat2;
 | 
			
		||||
    Real rho_munu, rho_numu;
 | 
			
		||||
    
 | 
			
		||||
    // to be completed 
 | 
			
		||||
    int Nvol = CommonPrms::instance()->Nvol();
 | 
			
		||||
    
 | 
			
		||||
    for(int mu = 0; mu < Nd; ++mu){
 | 
			
		||||
      U_mu       = PeekIndex<LorentzIndex>(      U, mu);
 | 
			
		||||
      iLambda_mu = PeekIndex<LorentzIndex>(iLambda, mu);
 | 
			
		||||
      
 | 
			
		||||
      
 | 
			
		||||
      for(int nu = 0; nu < Nd; ++nu){
 | 
			
		||||
	if(nu==mu) continue;
 | 
			
		||||
	U_nu       = PeekIndex<LorentzIndex>(      U, nu);
 | 
			
		||||
	iLambda_nu = PeekIndex<LorentzIndex>(iLambda, nu);
 | 
			
		||||
	
 | 
			
		||||
	rho_munu = rho[mu + Nd * nu];
 | 
			
		||||
	rho_numu = rho[nu + Nd * mu];
 | 
			
		||||
	
 | 
			
		||||
	WL.StapleUpper(staple, U, mu, nu);
 | 
			
		||||
	
 | 
			
		||||
	temp_Sigma = adj(staple)*iLambda_nu;
 | 
			
		||||
	temp_Sigma *= - rho_numu;
 | 
			
		||||
	//-r_numu*U_nu(x+mu)*Udag_mu(x+nu)*Udag_nu(x)*Lambda_nu(x)
 | 
			
		||||
	SigmaTerm .................
 | 
			
		||||
 | 
			
		||||
	for (int site = 0; site < Nvol; ++site){
 | 
			
		||||
	  temp_mat = mat_dag(staple,site) * mat(iLambda_nu,site);
 | 
			
		||||
	  temp_mat *= - rho_numu;
 | 
			
		||||
	  AddMat(SigmaTerm, temp_mat, site, mu);
 | 
			
		||||
	for(int mu=0; mu<Nd; ++mu){
 | 
			
		||||
	  Cup = zero;
 | 
			
		||||
	  for(int nu=0; nu<Nd; ++nu){
 | 
			
		||||
	    d_rho = rho[mu + Nd * nu];
 | 
			
		||||
	    WL.Staple(tmp_stpl, U, mu, nu);  //nb staple conventions of IroIro and Grid differ by a dag
 | 
			
		||||
	    Cup += tmp_stpl*d_rho;
 | 
			
		||||
	  }
 | 
			
		||||
	  pokeLorentz(u_smr, adj(Cup), mu); // u_smr[mu] = Cup^dag
 | 
			
		||||
	}
 | 
			
		||||
	sh_field = shiftField(iLambda_nu, mu, Forward());
 | 
			
		||||
	
 | 
			
		||||
	for (int site = 0; site < Nvol; ++site){
 | 
			
		||||
	  temp_mat = mat(sh_field,site) * mat_dag(staple,site);
 | 
			
		||||
	  temp_mat *= rho_numu;
 | 
			
		||||
	  AddMat(SigmaTerm, temp_mat, site, mu);
 | 
			
		||||
	}//r_numu*Lambda_nu(mu)*U_nu(x+mu)*Udag_mu(x+nu)*Udag_nu(x)
 | 
			
		||||
	
 | 
			
		||||
	sh_field = shiftField(iLambda_mu, nu, Forward());
 | 
			
		||||
	
 | 
			
		||||
	for (int site = 0; site < Nvol; ++site){
 | 
			
		||||
	  temp_mat = mat(U_nu,site) * mat(sh_field,site) * mat_dag(U_nu,site);
 | 
			
		||||
	  temp_mat = mat_dag(staple,site) * temp_mat;
 | 
			
		||||
	  temp_mat *= - rho_munu;
 | 
			
		||||
	  AddMat(SigmaTerm, temp_mat, site, mu);
 | 
			
		||||
	}//-r_munu*U_nu(x+mu)*Udag_mu(x+nu)*Lambda_mu(x+nu)*Udag_nu(x)
 | 
			
		||||
	
 | 
			
		||||
	staple = 0.0;
 | 
			
		||||
	sh_field = shiftField(U_nu, mu, Forward());
 | 
			
		||||
	
 | 
			
		||||
	for (int site = 0; site < Nvol; ++site){
 | 
			
		||||
	  temp_mat2 = mat_dag(sh_field,site) * mat_dag(U_mu,site);
 | 
			
		||||
	  temp_mat = temp_mat2  * mat(iLambda_mu,site) * mat(U_nu,site);
 | 
			
		||||
	  temp_mat *= - rho_munu;
 | 
			
		||||
	  AddMat(staple, temp_mat, site);
 | 
			
		||||
	  temp_mat = temp_mat2 * mat(iLambda_nu,site) * mat(U_nu,site);
 | 
			
		||||
	  temp_mat *= rho_numu;
 | 
			
		||||
	  AddMat(staple, temp_mat, site);
 | 
			
		||||
	} 
 | 
			
		||||
	
 | 
			
		||||
	for (int site = 0; site < Nvol; ++site){
 | 
			
		||||
	  temp_mat = mat_dag(U_nu,site) * mat(iLambda_nu,site);
 | 
			
		||||
	  SetMat(u_tmp, temp_mat, site);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	sh_field = shiftField(u_tmp, mu, Forward()); 
 | 
			
		||||
	
 | 
			
		||||
	for (int site = 0; site < Nvol; ++site){
 | 
			
		||||
	  temp_mat = mat(sh_field,site) * mat_dag(U_mu,site) * mat(U_nu,site);
 | 
			
		||||
	  temp_mat *= - rho_numu;
 | 
			
		||||
	  AddMat(staple, temp_mat, site);
 | 
			
		||||
	}   
 | 
			
		||||
	
 | 
			
		||||
	sh_field = shiftField(staple, nu, Backward());
 | 
			
		||||
	
 | 
			
		||||
	AddSlice(SigmaTerm, sh_field, mu);
 | 
			
		||||
      }
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    */
 | 
			
		||||
  }
 | 
			
		||||
      void derivative(GaugeField& SigmaTerm,
 | 
			
		||||
		      const GaugeField& iLambda,
 | 
			
		||||
		      const GaugeField& U)const{
 | 
			
		||||
    
 | 
			
		||||
	// Reference 
 | 
			
		||||
	// Morningstar, Peardon, Phys.Rev.D69,054501(2004)
 | 
			
		||||
	// Equation 75
 | 
			
		||||
 | 
			
		||||
    	GridBase *grid = U._grid;
 | 
			
		||||
	int vol = U._grid->gSites();
 | 
			
		||||
	
 | 
			
		||||
	WilsonLoops<Gimpl> WL;
 | 
			
		||||
	GaugeLinkField staple(grid), u_tmp(grid), iLambda_mu(grid), iLambda_nu(grid);
 | 
			
		||||
	GaugeLinkField U_mu(grid), U_nu(grid), sh_field(grid), temp_Sigma(grid);
 | 
			
		||||
	Real rho_munu, rho_numu;
 | 
			
		||||
    
 | 
			
		||||
	for(int mu = 0; mu < Nd; ++mu){
 | 
			
		||||
	  U_mu       = PeekIndex<LorentzIndex>(      U, mu);
 | 
			
		||||
	  iLambda_mu = PeekIndex<LorentzIndex>(iLambda, mu);
 | 
			
		||||
      
 | 
			
		||||
	  for(int nu = 0; nu < Nd; ++nu){
 | 
			
		||||
	    if(nu==mu) continue;
 | 
			
		||||
	    U_nu       = PeekIndex<LorentzIndex>(      U, nu);
 | 
			
		||||
	    iLambda_nu = PeekIndex<LorentzIndex>(iLambda, nu);
 | 
			
		||||
	
 | 
			
		||||
	    rho_munu = rho[mu + Nd * nu];
 | 
			
		||||
	    rho_numu = rho[nu + Nd * mu];
 | 
			
		||||
	
 | 
			
		||||
	    WL.StapleUpper(staple, U, mu, nu);
 | 
			
		||||
 | 
			
		||||
	    temp_Sigma = -rho_numu*staple*iLambda_nu;
 | 
			
		||||
	    //-r_numu*U_nu(x+mu)*Udag_mu(x+nu)*Udag_nu(x)*Lambda_nu(x)
 | 
			
		||||
	    AddGaugeLink(SigmaTerm, temp_Sigma, mu);
 | 
			
		||||
 | 
			
		||||
	    sh_field = Cshift(iLambda_nu, mu, 1);// general also for Gparity?
 | 
			
		||||
	
 | 
			
		||||
	    temp_Sigma = rho_numu*sh_field*staple;
 | 
			
		||||
	    //r_numu*Lambda_nu(mu)*U_nu(x+mu)*Udag_mu(x+nu)*Udag_nu(x)
 | 
			
		||||
	    AddGaugeLink(SigmaTerm, temp_Sigma, mu);
 | 
			
		||||
 | 
			
		||||
	    sh_field = Cshift(iLambda_mu, nu, 1);
 | 
			
		||||
 | 
			
		||||
	    temp_Sigma = -rho_munu*staple*U_nu*sh_field*adj(U_nu);
 | 
			
		||||
	    //-r_munu*U_nu(x+mu)*Udag_mu(x+nu)*Lambda_mu(x+nu)*Udag_nu(x)
 | 
			
		||||
	    AddGaugeLink(SigmaTerm, temp_Sigma, mu);
 | 
			
		||||
 | 
			
		||||
	    staple = zero;
 | 
			
		||||
	    sh_field = Cshift(U_nu, mu, 1);
 | 
			
		||||
 | 
			
		||||
	    temp_Sigma = -rho_munu*adj(sh_field)*adj(U_mu)*iLambda_mu*U_nu;
 | 
			
		||||
	    temp_Sigma += rho_numu*adj(sh_field)*adj(U_mu)*iLambda_nu*U_nu;
 | 
			
		||||
 | 
			
		||||
	    u_tmp = adj(U_nu)*iLambda_nu;
 | 
			
		||||
	    sh_field = Cshift(u_tmp, mu, 1);
 | 
			
		||||
	    temp_Sigma += -rho_numu*sh_field*adj(U_mu)*U_nu;
 | 
			
		||||
	    sh_field = Cshift(temp_Sigma, nu, -1);
 | 
			
		||||
	    AddGaugeLink(SigmaTerm, sh_field, mu);
 | 
			
		||||
	    
 | 
			
		||||
	  }
 | 
			
		||||
	}
 | 
			
		||||
      }
 | 
			
		||||
    };
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
  }// namespace QCD
 | 
			
		||||
}//namespace Grid
 | 
			
		||||
#endif  
 | 
			
		||||
 
 | 
			
		||||
@@ -22,7 +22,7 @@ namespace Grid {
 | 
			
		||||
      It stores a list of smeared configurations.
 | 
			
		||||
    */
 | 
			
		||||
    template <class Gimpl>
 | 
			
		||||
    class GaugeConfiguration {
 | 
			
		||||
    class SmearedConfiguration {
 | 
			
		||||
    public:
 | 
			
		||||
      INHERIT_GIMPL_TYPES(Gimpl)
 | 
			
		||||
      private:
 | 
			
		||||
@@ -58,7 +58,7 @@ namespace Grid {
 | 
			
		||||
				    links configuration */
 | 
			
		||||
      
 | 
			
		||||
      /*! @brief Standard constructor */
 | 
			
		||||
      GaugeConfiguration(GridCartesian * UGrid,
 | 
			
		||||
      SmearedConfiguration(GridCartesian * UGrid,
 | 
			
		||||
			 unsigned int Nsmear, 
 | 
			
		||||
			 Smear_Stout& Stout):
 | 
			
		||||
	smearingLevels(Nsmear),
 | 
			
		||||
@@ -69,7 +69,7 @@ namespace Grid {
 | 
			
		||||
      }
 | 
			
		||||
      
 | 
			
		||||
      /*! For just thin links */
 | 
			
		||||
      GaugeConfiguration(GridCartesian * UGrid):
 | 
			
		||||
      SmearedConfiguration(GridCartesian * UGrid):
 | 
			
		||||
	smearingLevels(0),
 | 
			
		||||
	StoutSmearing(),
 | 
			
		||||
	SmearedSet(0),
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										9
									
								
								lib/qcd/smearing/Smearing.h
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										9
									
								
								lib/qcd/smearing/Smearing.h
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,9 @@
 | 
			
		||||
#ifndef GRID_QCD_SMEARING_H
 | 
			
		||||
#define GRID_QCD_SMEARING_H
 | 
			
		||||
 | 
			
		||||
#include <qcd/smearing/BaseSmearing.h>
 | 
			
		||||
#include <qcd/smearing/APEsmearing.h>
 | 
			
		||||
#include <qcd/smearing/StoutSmearing.h>
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#endif
 | 
			
		||||
@@ -5,28 +5,132 @@
 | 
			
		||||
#ifndef STOUT_SMEAR_
 | 
			
		||||
#define STOUT_SMEAR_
 | 
			
		||||
 | 
			
		||||
/*!  @brief Stout smearing of link variable. */
 | 
			
		||||
template <class Gimpl> 
 | 
			
		||||
class Smear_Stout: public Smear<Gimpl> {
 | 
			
		||||
private:
 | 
			
		||||
  const std::valarray<double> d_rho;
 | 
			
		||||
  const Smear* SmearBase;
 | 
			
		||||
namespace Grid {
 | 
			
		||||
  namespace QCD {
 | 
			
		||||
 | 
			
		||||
  double func_xi0(double w) const;
 | 
			
		||||
public:
 | 
			
		||||
  INHERIT_GIMPL_TYPES(Gimpl)
 | 
			
		||||
  Smear_Stout(Smear* base):SmearBase(base){}
 | 
			
		||||
    /*!  @brief Stout smearing of link variable. */
 | 
			
		||||
    template <class Gimpl> 
 | 
			
		||||
    class Smear_Stout: public Smear<Gimpl> {
 | 
			
		||||
    private:
 | 
			
		||||
      const std::vector<double> d_rho;
 | 
			
		||||
      const Smear < Gimpl > * SmearBase;
 | 
			
		||||
      
 | 
			
		||||
      LatticeReal func_xi0(LatticeReal w) const{
 | 
			
		||||
	// Define a function to do the check
 | 
			
		||||
	//if( w < 1e-4 ) std::cout << GridLogWarning << "[Smear_stout] w too small: "<< w <<"\n";
 | 
			
		||||
	return  sin(w)/w;
 | 
			
		||||
      }
 | 
			
		||||
      
 | 
			
		||||
    public:
 | 
			
		||||
      INHERIT_GIMPL_TYPES(Gimpl)
 | 
			
		||||
      
 | 
			
		||||
      Smear_Stout(Smear < Gimpl >* base):SmearBase(base){}
 | 
			
		||||
      
 | 
			
		||||
      /*! Default constructor */
 | 
			
		||||
      Smear_Stout():SmearBase(new Smear_APE < Gimpl > ()){}
 | 
			
		||||
      
 | 
			
		||||
      ~Smear_Stout(){}
 | 
			
		||||
      
 | 
			
		||||
      void smear(GaugeField& u_smr,const GaugeField& U) const{
 | 
			
		||||
	long double timing;
 | 
			
		||||
	
 | 
			
		||||
	GaugeField u_tmp1, q_mu;
 | 
			
		||||
	
 | 
			
		||||
	std::cout<< GridLogDebug << "Stout smearing started\n";
 | 
			
		||||
	
 | 
			
		||||
	//Smear the configurations
 | 
			
		||||
	SmearBase->smear(u_tmp1, U);
 | 
			
		||||
	
 | 
			
		||||
	q_mu = Ta(u_tmp1*adj(u_tmp1)); // q_mu = Ta(Omega_mu)
 | 
			
		||||
	
 | 
			
		||||
	exponentiate_iQ(u_tmp1, q_mu);
 | 
			
		||||
 | 
			
		||||
  /*! Default constructor */
 | 
			
		||||
  Smear_Stout():SmearBase(new Smear_APE()){}
 | 
			
		||||
	u_smr = u_tmp1*U;
 | 
			
		||||
	
 | 
			
		||||
	std::cout<< GridLogDebug << "Stout smearing completed\n";
 | 
			
		||||
      }
 | 
			
		||||
      void derivative(GaugeField& SigmaTerm,
 | 
			
		||||
		      const GaugeField& iLambda,
 | 
			
		||||
		      const GaugeField& Gauge) const{
 | 
			
		||||
	SmearBase->derivative(SigmaTerm, iLambda, Gauge);
 | 
			
		||||
      }
 | 
			
		||||
      
 | 
			
		||||
      
 | 
			
		||||
      void BaseSmear(GaugeField& C,
 | 
			
		||||
		     const GaugeField& U) const{
 | 
			
		||||
	SmearBase->smear(C, U);
 | 
			
		||||
      }
 | 
			
		||||
      
 | 
			
		||||
      void exponentiate_iQ(GaugeField& e_iQ,
 | 
			
		||||
			   const GaugeField& iQ) const{
 | 
			
		||||
	// Put this outside 
 | 
			
		||||
	// only valid for SU(3) matrices
 | 
			
		||||
 | 
			
		||||
  ~Smear_Stout(){}
 | 
			
		||||
	GridBase *grid = iQ._grid;
 | 
			
		||||
	Real one_over_three = 1.0/3.0;
 | 
			
		||||
	Real one_over_two = 1.0/2.0;
 | 
			
		||||
 | 
			
		||||
  void smear(GaugeField&,const GaugeField&) const;
 | 
			
		||||
  void BaseSmear(GaugeField&, const GaugeField&) const;
 | 
			
		||||
  void derivative(GaugeField&, const GaugeField&, const GaugeField&) const;
 | 
			
		||||
  void exponentiate_iQ(GaugeField&, const GaugeField&) const;
 | 
			
		||||
	GaugeField unity;
 | 
			
		||||
	GaugeLinkField Umu(iQ._grid);
 | 
			
		||||
	Umu=1.0;
 | 
			
		||||
	for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
	  pokeLorentz(unity,Umu,mu);
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
	
 | 
			
		||||
	GaugeField iQ2, iQ3;
 | 
			
		||||
	LatticeReal c0(grid), c1(grid), c0max(grid), u_val(grid), tmp(grid);
 | 
			
		||||
	LatticeReal w(grid), theta(grid), xi0(grid), u2(grid), w2(grid), cosw(grid);
 | 
			
		||||
	LatticeComplex fden(grid);
 | 
			
		||||
	LatticeComplex f0(grid), f1(grid), f2(grid), h0(grid), h1(grid), h2(grid);
 | 
			
		||||
	LatticeComplex e2iu(grid), emiu(grid), ixi0(grid), qt(grid);
 | 
			
		||||
 | 
			
		||||
	iQ2 = iQ * iQ;
 | 
			
		||||
	iQ3 = iQ * iQ2;
 | 
			
		||||
 | 
			
		||||
	c0    = - imag(trace(iQ3)) * one_over_three;
 | 
			
		||||
	c1    = - real(trace(iQ2)) * one_over_two;
 | 
			
		||||
	tmp   = c1 * one_over_three;
 | 
			
		||||
    	c0max = 2.0 * pow(tmp, 1.5);
 | 
			
		||||
 | 
			
		||||
	theta = acos(c0/c0max);
 | 
			
		||||
 | 
			
		||||
	u_val = sqrt(tmp) * cos( theta * one_over_three);
 | 
			
		||||
	w     = sqrt(c1) * sin ( theta * one_over_three);
 | 
			
		||||
	xi0   = func_xi0(w);
 | 
			
		||||
	u2    = u_val * u_val;
 | 
			
		||||
	w2    = w * w;
 | 
			
		||||
	cosw  = cos(w);
 | 
			
		||||
 | 
			
		||||
	ixi0  = timesI(toComplex(xi0));
 | 
			
		||||
	emiu  = toComplex(cos(u_val)) - timesI(toComplex(u_val));
 | 
			
		||||
	e2iu  = toComplex(cos(2.0*u_val)) + timesI(toComplex(2.0*u_val));
 | 
			
		||||
	
 | 
			
		||||
	h0    = e2iu * toComplex(u2 - w2) + emiu *( toComplex(8.0*u2*cosw) +
 | 
			
		||||
						    toComplex(2.0*u_val*(3.0*u2 + w2))*ixi0);
 | 
			
		||||
	
 | 
			
		||||
	h1    = toComplex(2.0*u_val) * e2iu - emiu*( toComplex(2.0*u_val*cosw) -
 | 
			
		||||
						     toComplex(3.0*u2-w2)*ixi0);
 | 
			
		||||
	
 | 
			
		||||
	h2    = e2iu - emiu * (toComplex(cosw) + toComplex(3.0*u_val)*ixi0);
 | 
			
		||||
 | 
			
		||||
	
 | 
			
		||||
	tmp   = 9.0*u2 - w2;
 | 
			
		||||
	fden  = toComplex(pow(tmp, -1.0));
 | 
			
		||||
	f0    = h0 * fden;
 | 
			
		||||
	f1    = h1 * fden;
 | 
			
		||||
	f2    = h2 * fden;
 | 
			
		||||
 | 
			
		||||
	
 | 
			
		||||
	e_iQ = f0*unity + f1 * timesMinusI(iQ) - f2 * iQ2;
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
      };
 | 
			
		||||
      
 | 
			
		||||
    };
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
#endif  
 | 
			
		||||
 
 | 
			
		||||
		Reference in New Issue
	
	Block a user