1
0
mirror of https://github.com/paboyle/Grid.git synced 2025-06-12 20:27:06 +01:00

Elemental force term for Wilson dslash added and tests thereof passing.

Now need to construct pseudofermion two flavour, ratio, one flavour, ratio
action fragments.
This commit is contained in:
Peter Boyle
2015-07-26 10:54:38 +09:00
parent ba4989dd45
commit d7e6b65a76
18 changed files with 758 additions and 42 deletions

View File

@ -33,9 +33,6 @@
/* Support AVX2 (Advanced Vector Extensions 2) instructions */
#undef HAVE_AVX2
/* define if the compiler supports basic C++11 syntax */
#undef HAVE_CXX11
/* Define to 1 if you have the declaration of `be64toh', and to 0 if you
don't. */
#undef HAVE_DECL_BE64TOH

View File

@ -60,6 +60,11 @@ public:
GridBase *_grid;
int checkerboard;
std::vector<vobj,alignedAllocator<vobj> > _odata;
// to pthread need a computable loop where loop induction is not required
int begin(void) { return 0;};
int end(void) { return _odata.size(); }
vobj & operator[](int i) { return _odata[i]; };
public:
typedef typename vobj::scalar_type scalar_type;

View File

@ -39,10 +39,23 @@ namespace Grid {
virtual void Dhop (const FermionField &in, FermionField &out,int dag)=0;
virtual void DhopOE(const FermionField &in, FermionField &out,int dag)=0;
virtual void DhopEO(const FermionField &in, FermionField &out,int dag)=0;
virtual void DhopDir(const FermionField &in, FermionField &out,int dir,int disp)=0; // implemented by WilsonFermion and WilsonFermion5D
// force terms; five routines; default to Dhop on diagonal
virtual void MDeriv (LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag){DhopDeriv(mat,U,V,dag);};
virtual void MoeDeriv(LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag){DhopDerivOE(mat,U,V,dag);};
virtual void MeoDeriv(LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag){DhopDerivEO(mat,U,V,dag);};
virtual void MooDeriv(LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag){mat=zero;};
virtual void MeeDeriv(LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag){mat=zero;};
virtual void DhopDeriv (LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag)=0;
virtual void DhopDerivEO(LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag)=0;
virtual void DhopDerivOE(LatticeGaugeField &mat,const FermionField &U,const FermionField &V,int dag)=0;
virtual void Mdiag (const FermionField &in, FermionField &out) { Mooee(in,out);}; // Same as Mooee applied to both CB's
virtual void Mdir (const FermionField &in, FermionField &out,int dir,int disp)=0; // case by case Wilson, Clover, Cayley, ContFrac, PartFrac
virtual void Mdiag(const FermionField &in, FermionField &out) { Mooee(in,out);}; // Same as Mooee applied to both CB's
virtual void Mdir (const FermionField &in, FermionField &out,int dir,int disp)=0; // case by case Wilson, Clover, Cayley, ContFrac, PartFrac
virtual void DhopDir(const FermionField &in, FermionField &out,int dir,int disp)=0; // implemented by WilsonFermion and WilsonFermion5D
};

View File

@ -98,7 +98,9 @@ void WilsonFermion::Mdir (const LatticeFermion &in, LatticeFermion &out,int dir,
DhopDir(in,out,dir,disp);
}
void WilsonFermion::DhopDir(const LatticeFermion &in, LatticeFermion &out,int dir,int disp){
WilsonCompressor compressor(DaggerNo);
Stencil.HaloExchange<vSpinColourVector,vHalfSpinColourVector,WilsonCompressor>(in,comm_buf,compressor);
assert( (disp==1)||(disp==-1) );
@ -109,9 +111,22 @@ void WilsonFermion::DhopDir(const LatticeFermion &in, LatticeFermion &out,int di
PARALLEL_FOR_LOOP
for(int sss=0;sss<in._grid->oSites();sss++){
DiracOptDhopDir(Stencil,Umu,comm_buf,sss,sss,in,out,dirdisp);
DiracOptDhopDir(Stencil,Umu,comm_buf,sss,sss,in,out,dirdisp,dirdisp);
}
};
void WilsonFermion::DhopDirDisp(const LatticeFermion &in, LatticeFermion &out,int dirdisp,int gamma,int dag)
{
WilsonCompressor compressor(dag);
Stencil.HaloExchange<vSpinColourVector,vHalfSpinColourVector,WilsonCompressor>(in,comm_buf,compressor);
PARALLEL_FOR_LOOP
for(int sss=0;sss<in._grid->oSites();sss++){
DiracOptDhopDir(Stencil,Umu,comm_buf,sss,sss,in,out,dirdisp,gamma);
}
};
void WilsonFermion::DhopInternal(CartesianStencil & st,LatticeDoubledGaugeField & U,
@ -177,6 +192,77 @@ void WilsonFermion::Dhop(const LatticeFermion &in, LatticeFermion &out,int dag)
DhopInternal(Stencil,Umu,in,out,dag);
}
void WilsonFermion::DerivInternal(CartesianStencil & st,LatticeDoubledGaugeField & U,
LatticeGaugeField &mat,const LatticeFermion &A,const LatticeFermion &B,int dag)
{
assert((dag==DaggerNo) ||(dag==DaggerYes));
WilsonCompressor compressor(dag);
LatticeColourMatrix tmp(B._grid);
LatticeFermion Btilde(B._grid);
st.HaloExchange<vSpinColourVector,vHalfSpinColourVector,WilsonCompressor>(B,comm_buf,compressor);
for(int mu=0;mu<Nd;mu++){
////////////////////////////////////////////////////////////////////////
// Flip gamma (1+g)<->(1-g) if dag
////////////////////////////////////////////////////////////////////////
int gamma = mu;
if ( dag ) gamma+= Nd;
////////////////////////
// Call the single hop
////////////////////////
PARALLEL_FOR_LOOP
for(int sss=0;sss<B._grid->oSites();sss++){
DiracOptDhopDir(st,U,comm_buf,sss,sss,B,Btilde,mu,gamma);
}
//////////////////////////////////////////////////
// spin trace outer product
//////////////////////////////////////////////////
tmp = TraceIndex<SpinIndex>(outerProduct(Btilde,A));
PokeIndex<LorentzIndex>(mat,tmp,mu);
}
}
void WilsonFermion::DhopDeriv(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag)
{
conformable(U._grid,_grid);
conformable(U._grid,V._grid);
conformable(U._grid,mat._grid);
mat.checkerboard = U.checkerboard;
DerivInternal(Stencil,Umu,mat,U,V,dag);
}
void WilsonFermion::DhopDerivOE(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag)
{
conformable(U._grid,_cbgrid);
conformable(U._grid,V._grid);
conformable(U._grid,mat._grid);
assert(V.checkerboard==Even);
assert(U.checkerboard==Odd);
mat.checkerboard = Odd;
DerivInternal(StencilEven,UmuOdd,mat,U,V,dag);
}
void WilsonFermion::DhopDerivEO(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag)
{
conformable(U._grid,_cbgrid);
conformable(U._grid,V._grid);
conformable(U._grid,mat._grid);
assert(V.checkerboard==Odd);
assert(U.checkerboard==Even);
mat.checkerboard = Even;
DerivInternal(StencilOdd,UmuEven,mat,U,V,dag);
}
}}

View File

@ -24,11 +24,95 @@ namespace Grid {
// half checkerboard operaions
void Meooe (const LatticeFermion &in, LatticeFermion &out);
void MeooeDag (const LatticeFermion &in, LatticeFermion &out);
virtual void Mooee (const LatticeFermion &in, LatticeFermion &out); // remain virtual so we
virtual void MooeeDag (const LatticeFermion &in, LatticeFermion &out); // can derive Clover
virtual void MooeeInv (const LatticeFermion &in, LatticeFermion &out); // from Wilson base
virtual void MooeeInvDag (const LatticeFermion &in, LatticeFermion &out);
////////////////////////
//
// Force term: d/dtau S = 0
//
// It is simplest to consider the two flavour force term
//
// S[U,phi] = phidag (MdagM)^-1 phi
//
// But simplify even this to
//
// S[U,phi] = phidag MdagM phi
//
// (other options exist depending on nature of action fragment.)
//
// Require momentum be traceless anti-hermitian to move within group manifold [ P = i P^a T^a ]
//
// Define the HMC hamiltonian
//
// H = 1/2 Tr P^2 + S(U,phi)
//
// .
// U = P U (lorentz & color indices multiplied)
//
// Hence
//
// .c c c c
// U = U P = - U P (c == dagger)
//
// So, taking some liberty with implicit indices
// . . .c c
// dH/dt = 0 = Tr P P +Tr[ U dS/dU + U dS/dU ]
//
// . c c
// = Tr P P + i Tr[ P U dS/dU - U P dS/dU ]
//
// . c c
// = Tr P (P + i ( U dS/dU - P dS/dU U ]
//
// . c c
// => P = -i [ U dS/dU - dS/dU U ] generates HMC EoM
//
// Simple case work this out using S = phi^dag MdagM phi for wilson:
// c c
// dSdt = dU_xdt dSdUx + dUxdt dSdUx
//
// = Tr i P U_x [ (\phi^\dag)_x (1+g) (M \phi)_x+\mu +(\phi^\dag M^\dag)_x (1-g) \phi_{x+\mu} ]
// c
// - i U_x P [ (\phi^\dag)_x+mu (1-g) (M \phi)_x +(\phi^\dag M^\dag)_(x+\mu) (1+g) \phi_{x} ]
//
// = i [(\phi^\dag)_x ]_j P_jk [U_x(1+g) (M \phi)_x+\mu]_k (1)
// + i [(\phi^\dagM^\dag)_x]_j P_jk [U_x(1-g) (\phi)_x+\mu]_k (2)
// - i [(\phi^\dag)_x+mu (1-g) U^dag_x]_j P_jk [(M \phi)_xk (3)
// - i [(\phi^\dagM^\dag)_x+mu (1+g) U^dag_x]_j P_jk [ \phi]_xk (4)
//
// Observe that (1)* = (4)
// (2)* = (3)
//
// Write as .
// P_{kj} = - i ( [U_x(1+g) (M \phi)_x+\mu] (x) [(\phi^\dag)_x] + [U_x(1-g) (\phi)_x+\mu] (x) [(\phi^\dagM^\dag)_x] - h.c )
//
// where (x) denotes outer product in colour and spins are traced.
//
// Need only evaluate (1) and (2) [Chroma] or (2) and (4) [IroIro] and take the
// traceless anti hermitian part (of term in brackets w/o the "i")
//
// Generalisation to S=phi^dag (MdagM)^{-1} phi is simple:
//
// For more complicated DWF etc... apply product rule in differentiation
//
////////////////////////
void DhopDeriv (LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag);
void DhopDerivEO(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag);
void DhopDerivOE(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag);
// Extra support internal
void DerivInternal(CartesianStencil & st,
LatticeDoubledGaugeField & U,
LatticeGaugeField &mat,
const LatticeFermion &A,
const LatticeFermion &B,
int dag);
// non-hermitian hopping term; half cb or both
void Dhop (const LatticeFermion &in, LatticeFermion &out,int dag);
void DhopOE(const LatticeFermion &in, LatticeFermion &out,int dag);
@ -37,6 +121,7 @@ namespace Grid {
// Multigrid assistance
void Mdir (const LatticeFermion &in, LatticeFermion &out,int dir,int disp);
void DhopDir(const LatticeFermion &in, LatticeFermion &out,int dir,int disp);
void DhopDirDisp(const LatticeFermion &in, LatticeFermion &out,int dirdisp,int gamma,int dag);
///////////////////////////////////////////////////////////////
// Extra methods added by derived
@ -59,7 +144,8 @@ namespace Grid {
static int HandOptDslash; // these are a temporary hack
static int MortonOrder;
protected:
// protected:
public:
RealD mass;

View File

@ -105,14 +105,105 @@ PARALLEL_FOR_LOOP
for(int s=0;s<Ls;s++){
int sU=ss;
int sF = s+Ls*sU;
DiracOptDhopDir(Stencil,Umu,comm_buf,sF,sU,in,out,dirdisp);
DiracOptDhopDir(Stencil,Umu,comm_buf,sF,sU,in,out,dirdisp,dirdisp);
}
}
};
void WilsonFermion5D::DerivInternal(CartesianStencil & st,
LatticeDoubledGaugeField & U,
LatticeGaugeField &mat,
const LatticeFermion &A,
const LatticeFermion &B,
int dag)
{
assert((dag==DaggerNo) ||(dag==DaggerYes));
WilsonCompressor compressor(dag);
LatticeColourMatrix tmp(B._grid);
LatticeFermion Btilde(B._grid);
st.HaloExchange<vSpinColourVector,vHalfSpinColourVector,WilsonCompressor>(B,comm_buf,compressor);
for(int mu=0;mu<Nd;mu++){
////////////////////////////////////////////////////////////////////////
// Flip gamma if dag
////////////////////////////////////////////////////////////////////////
int gamma = mu;
if ( dag ) gamma+= Nd;
////////////////////////
// Call the single hop
////////////////////////
PARALLEL_FOR_LOOP
for(int sss=0;sss<B._grid->oSites();sss++){
for(int s=0;s<Ls;s++){
int sU=sss;
int sF = s+Ls*sU;
DiracOptDhopDir(st,U,comm_buf,sF,sU,B,Btilde,mu,gamma);
}
}
////////////////////////////
// spin trace outer product
////////////////////////////
tmp = TraceIndex<SpinIndex>(outerProduct(Btilde,A)); // ordering here
PokeIndex<LorentzIndex>(mat,tmp,mu);
}
}
void WilsonFermion5D::DhopDeriv( LatticeGaugeField &mat,
const LatticeFermion &A,
const LatticeFermion &B,
int dag)
{
conformable(A._grid,FermionGrid());
conformable(A._grid,B._grid);
conformable(GaugeGrid(),mat._grid);
mat.checkerboard = A.checkerboard;
DerivInternal(Stencil,Umu,mat,A,B,dag);
}
void WilsonFermion5D::DhopDerivEO(LatticeGaugeField &mat,
const LatticeFermion &A,
const LatticeFermion &B,
int dag)
{
conformable(A._grid,FermionRedBlackGrid());
conformable(GaugeRedBlackGrid(),mat._grid);
conformable(A._grid,B._grid);
assert(B.checkerboard==Odd);
assert(A.checkerboard==Even);
mat.checkerboard = Even;
DerivInternal(StencilOdd,UmuEven,mat,A,B,dag);
}
void WilsonFermion5D::DhopDerivOE(LatticeGaugeField &mat,
const LatticeFermion &A,
const LatticeFermion &B,
int dag)
{
conformable(A._grid,FermionRedBlackGrid());
conformable(GaugeRedBlackGrid(),mat._grid);
conformable(A._grid,B._grid);
assert(B.checkerboard==Even);
assert(A.checkerboard==Odd);
mat.checkerboard = Odd;
DerivInternal(StencilEven,UmuOdd,mat,A,B,dag);
}
void WilsonFermion5D::DhopInternal(CartesianStencil & st, LebesgueOrder &lo,
LatticeDoubledGaugeField & U,
const LatticeFermion &in, LatticeFermion &out,int dag)
const LatticeFermion &in, LatticeFermion &out,int dag)
{
// assert((dag==DaggerNo) ||(dag==DaggerYes));

View File

@ -44,12 +44,18 @@ namespace Grid {
// half checkerboard operations; leave unimplemented as abstract for now
virtual void Meooe (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void MeooeDag (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void Mooee (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void MooeeDag (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void MooeeInv (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void MeooeDag (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void MooeeDag (const LatticeFermion &in, LatticeFermion &out){assert(0);};
virtual void MooeeInvDag (const LatticeFermion &in, LatticeFermion &out){assert(0);};
// These can be overridden by fancy 5d chiral actions
virtual void DhopDeriv (LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag);
virtual void DhopDerivEO(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag);
virtual void DhopDerivOE(LatticeGaugeField &mat,const LatticeFermion &U,const LatticeFermion &V,int dag);
// Implement hopping term non-hermitian hopping term; half cb or both
// Implement s-diagonal DW
void DW (const LatticeFermion &in, LatticeFermion &out,int dag);
@ -64,6 +70,14 @@ namespace Grid {
///////////////////////////////////////////////////////////////
// New methods added
///////////////////////////////////////////////////////////////
void DerivInternal(CartesianStencil & st,
LatticeDoubledGaugeField & U,
LatticeGaugeField &mat,
const LatticeFermion &A,
const LatticeFermion &B,
int dag);
void DhopInternal(CartesianStencil & st,
LebesgueOrder &lo,
LatticeDoubledGaugeField &U,

View File

@ -294,8 +294,8 @@ void DiracOptDhopSiteDag(CartesianStencil &st,LatticeDoubledGaugeField &U,
}
void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
std::vector<vHalfSpinColourVector,alignedAllocator<vHalfSpinColourVector> > &buf,
int sF,int sU,const LatticeFermion &in, LatticeFermion &out,int dirdisp)
std::vector<vHalfSpinColourVector,alignedAllocator<vHalfSpinColourVector> > &buf,
int sF,int sU,const LatticeFermion &in, LatticeFermion &out,int dir,int gamma)
{
vHalfSpinColourVector tmp;
vHalfSpinColourVector chi;
@ -304,13 +304,13 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
int offset,local,perm, ptype;
int ss=sF;
offset = st._offsets [dirdisp][ss];
local = st._is_local[dirdisp][ss];
perm = st._permute[dirdisp][ss];
ptype = st._permute_type[dirdisp];
offset = st._offsets [dir][ss];
local = st._is_local[dir][ss];
perm = st._permute[dir][ss];
ptype = st._permute_type[dir];
// Xp
if(dirdisp==Xp){
if(gamma==Xp){
if ( local && perm ) {
spProjXp(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -319,12 +319,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Xp),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconXp(result,Uchi);
}
// Yp
if ( dirdisp==Yp ){
if ( gamma==Yp ){
if ( local && perm ) {
spProjYp(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -333,12 +333,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Yp),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconYp(result,Uchi);
}
// Zp
if ( dirdisp ==Zp ){
if ( gamma ==Zp ){
if ( local && perm ) {
spProjZp(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -347,12 +347,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Zp),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconZp(result,Uchi);
}
// Tp
if ( dirdisp ==Tp ){
if ( gamma ==Tp ){
if ( local && perm ) {
spProjTp(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -361,12 +361,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Tp),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconTp(result,Uchi);
}
// Xm
if ( dirdisp==Xm ){
if ( gamma==Xm ){
if ( local && perm ) {
spProjXm(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -375,12 +375,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Xm),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconXm(result,Uchi);
}
// Ym
if ( dirdisp == Ym ){
if ( gamma == Ym ){
if ( local && perm ) {
spProjYm(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -389,12 +389,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Ym),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconYm(result,Uchi);
}
// Zm
if ( dirdisp == Zm ){
if ( gamma == Zm ){
if ( local && perm ) {
spProjZm(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -403,12 +403,12 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Zm),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconZm(result,Uchi);
}
// Tm
if ( dirdisp==Tm ) {
if ( gamma==Tm ) {
if ( local && perm ) {
spProjTm(tmp,in._odata[offset]);
permute(chi,tmp,ptype);
@ -417,7 +417,7 @@ void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
} else {
chi=buf[offset];
}
mult(&Uchi(),&U._odata[sU](Tm),&chi());
mult(&Uchi(),&U._odata[sU](dir),&chi());
spReconTm(result,Uchi);
}

View File

@ -22,7 +22,7 @@ namespace Grid {
int sF,int sU,const LatticeFermion &in, LatticeFermion &out);
void DiracOptDhopDir(CartesianStencil &st,LatticeDoubledGaugeField &U,
std::vector<vHalfSpinColourVector,alignedAllocator<vHalfSpinColourVector> > &buf,
int sF,int sU,const LatticeFermion &in, LatticeFermion &out,int dirdisp);
int sF,int sU,const LatticeFermion &in, LatticeFermion &out,int dirdisp,int gamma);
// };

View File

@ -68,14 +68,12 @@ namespace Grid{
Integrator<Algorithm>& MolDyn):
Params(Pms),MD(MolDyn){
//FIXME
// initialize RNGs also with seed
sRNG.SeedRandomDevice();
}
~HybridMonteCarlo(){};
void evolve(LatticeLorentzColourMatrix& Uin){
Real DeltaH;

View File

@ -414,7 +414,7 @@ namespace Grid {
template<class S, class V >
inline Grid_simd< S, V> outerProduct(const Grid_simd< S, V> &l, const Grid_simd< S, V> & r)
{
return l*r;
return l*conjugate(r);
}
template<class S, class V >

View File

@ -28,11 +28,13 @@ auto outerProduct (const iScalar<l>& lhs,const iScalar<r>& rhs) -> iScalar<declt
inline ComplexF outerProduct(const ComplexF &l, const ComplexF& r)
{
return l*r;
std::cout << "outer product taking conj "<<r<<" "<<conj(r)<<std::endl;
return l*conj(r);
}
inline ComplexD outerProduct(const ComplexD &l, const ComplexD& r)
{
return l*r;
std::cout << "outer product taking conj "<<r<<" "<<conj(r)<<std::endl;
return l*conj(r);
}
inline RealF outerProduct(const RealF &l, const RealF& r)
{

View File

@ -31,6 +31,16 @@ inline auto trace(const iScalar<vtype> &arg) -> iScalar<decltype(trace(arg._inte
return ret;
}
template<class vtype,int N>
inline auto trace(const iVector<vtype,N> &arg) -> iVector<decltype(trace(arg._internal[0])),N>
{
iVector<decltype(trace(arg._internal[0])),N> ret;
for(int i=0;i<N;i++){
ret._internal[i]=trace(arg._internal[i]);
}
return ret;
}
}
#endif