/*******************************************************************************
*
* TRIQS: a Toolbox for Research in Interacting Quantum Systems
*
* Copyright (C) 2011 by M. Ferrero, O. Parcollet
*
* TRIQS is free software: you can redistribute it and/or modify it under the
* terms of the GNU General Public License as published by the Free Software
* Foundation, either version 3 of the License, or (at your option) any later
* version.
*
* TRIQS is distributed in the hope that it will be useful, but WITHOUT ANY
* WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
* FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
* details.
*
* You should have received a copy of the GNU General Public License along with
* TRIQS. If not, see .
*
******************************************************************************/
#include "tight_binding.hpp"
#include
#include
#include "grid_generator.hpp"
#include "functors.hpp"
namespace triqs { namespace lattice_tools {
using namespace std;
using namespace tqa;
//------------------------------------------------------
array_view hopping_stack (tight_binding const & TB, array_view const & k_stack) {
auto TK = Fourier(TB);
array res(TB.n_bands(), TB.n_bands(), k_stack.shape(1));
for(size_t i = 0; i energies_on_bz_path(tight_binding const & TB, K_view_type K1, K_view_type K2, size_t n_pts) {
auto TK = Fourier(TB);
const size_t norb=TB.lattice().n_orbitals();
const size_t ndim=TB.lattice().dim();
array eval(norb,n_pts);
K_type dk = (K2 - K1)/double(n_pts), k = K1;
for (size_t i =0; i energies_on_bz_grid(tight_binding const & TB, size_t n_pts) {
auto TK = Fourier(TB);
const size_t norb=TB.lattice().n_orbitals();
const size_t ndim=TB.lattice().dim();
grid_generator grid(ndim,n_pts);
array eval(norb,grid.size());
for (; grid ; ++grid) {
eval(range(),grid.index()) = linalg::eigenvalues( TK( (*grid) (range(0,ndim))), false);
}
return eval;
}
//------------------------------------------------------
std::pair, array > dos(tight_binding const & TB, size_t nkpts, size_t neps) {
// The Fourier transform of TK
// auto TK = Fourier(TB); // C++0x ....
auto TK = Fourier(TB);
// loop on the BZ
const size_t ndim=TB.lattice().dim();
const size_t norb=TB.lattice().n_orbitals();
grid_generator grid(ndim,nkpts);
array tempeval(norb);
array evec(norb,norb,grid.size());
array eval(norb,grid.size());
if (norb ==1)
for (; grid ; ++grid) {
double ee = real(TK( (*grid) (range(0,ndim)))(0,0));
eval(0,grid.index()) =ee;
evec(0,0,grid.index()) =1;
}
else
for (; grid ; ++grid) {
//cerr<<" index = "< eval_sl = eval(range(),grid.index());
array_view evec_sl = evec(range(),range(),grid.index());
std::tie (eval_sl,evec_sl) = linalg::eigenelements( TK( (*grid) (range(0,ndim)))); //, true);
//cerr<< " point "<< *grid << " value "<< eval_sl<< endl; //" "<< (*grid) (range(0,ndim)) << endl;
}
// define the epsilon mesh, etc.
array epsilon(neps);
double epsmax = tqa::max_element(eval);
double epsmin = tqa::min_element(eval);
double deps=(epsmax-epsmin)/neps;
//for (size_t i =0; i< neps; ++i) epsilon(i)= epsmin+i/(neps-1.0)*(epsmax-epsmin);
for (size_t i =0; i< neps; ++i) epsilon(i)=epsmin+(i+0.5)*deps;
// bin the eigenvalues according to their energy
// NOTE: a is defined as an integer. it is the index for the DOS.
//REPORT <<"Starting Binning ...."< rho (neps,norb);rho()=0;
for(size_t l=0;l, array > dos_patch(tight_binding const & TB, const array & triangles, size_t neps, size_t ndiv) {
// WARNING: This version only works for a single band Hamiltonian in 2 dimensions!!!!
// triangles is an array of points defining the triangles of the patch
// neps in the number of bins in energy
// ndiv in the number of divisions used to divide the triangles
//const size_t ndim=TB.lattice().dim();
//const size_t norb=TB.lattice().n_orbitals();
int ntri = triangles.shape(0)/3;
array dos(neps);
// Check consistency
const size_t ndim=TB.lattice().dim();
//const size_t norb=TB.lattice().n_orbitals();
if (ndim !=2) TRIQS_RUNTIME_ERROR<<"dos_patch : dimension 2 only !";
if (triangles.shape(1) != ndim) TRIQS_RUNTIME_ERROR<<"dos_patch : the second dimension of the 'triangle' array in not "< thop(1,1);
array energ(nk), weight(nk);
// a, b, c are the corners of the triangle
// g the center of gravity taken from a
array a(ndim), b(ndim), c(ndim), g(ndim), rv(ndim);
int pt = 0;
double s, t;
// The Fourier transform of TK
auto TK = Fourier(TB);
// loop over the triangles
for (int tri = 0; tri < ntri; tri++) {
a = triangles(pt,range());
pt++;
b = triangles(pt,range());
pt++;
c = triangles(pt,range());
pt++;
g = ((a+b+c)/3.0-a)/double(ndiv);
// the area around a k point might be different from one triangle to the other
// so I use it to weight the sum in the dos
double area = abs(0.5*((b(0)-a(0))*(c(1)-a(1))-(b(1)-a(1))*(c(0)-a(0)))/(ndiv*ndiv));
for (size_t i = 0; iepsmax) epsmax=energ(k_index);
if (energ(k_index) epsilon(neps);
double deps=(epsmax-epsmin)/neps;
for (size_t i =0; i< neps; ++i) epsilon(i)= epsmin+i/(neps-1.0)*(epsmax-epsmin);
// bin the eigenvalues according to their energy
int ind;
double totalweight(0.0);
dos() = 0.0;
for (size_t j = 0; j < nk; j++) {
ind=int((energ(j)-epsmin)/deps);
if (ind == int(neps)) ind--;
dos(ind) += weight(j);
totalweight += weight(j);
}
dos /= deps;// Normalize the DOS
return std::make_pair(epsilon, dos);
}
}}