3
0
mirror of https://github.com/triqs/dft_tools synced 2024-11-01 11:43:47 +01:00
dft_tools/triqs/gfs/local/fourier_matsubara.hpp
Olivier Parcollet 39edb2f846 [API change] gf : factories -> constructors
- Make more general constructors for the gf.
  gf( mesh, target_shape_t)
- remove the old make_gf for the basic gf.
- 2 var non generic gf removed.
- clean evaluator
- add tensor_valued
- add a simple vertex test.
- clean specialisation
- Fix bug introduced in 1906dc3
- forgot to resize the gf in new version of operator =
- Fix make_singularity in gf.hpp

- clean resize in operator =

- update h5 read/write for block gf
  - changed a bit the general trait to save *all* the gf.
  - allows a more general specialization, then a correct for blocks

- NOT FINISHED : need to save the block indice for python.
  How to reread ?
  Currently it read the blocks names and reconstitute the mesh from it.
  Is it sufficient ?

- clean block constructors

 - block constructors simplest possible : an int for the number of blocks
 - rest in free factories.
 - fixed the generic constructor from GfType for the regular type :
   only enable iif GfType is ImmutableGreenFunction

- multivar. fix linear index in C, and h5 format

  - linear index now correctly flatten in C mode
    (was in fortran mode), using a simple reverse of the tuple in the folding.
  - fix the h5 read write of the multivar fonctions
   in order to write an array on dimension # variables + dim_target
   i.e. without flattening the indices of the meshes.
   Easier for later data analysis, e.g. in Python.

- merge matrix/tensor_valued. improve factories

  - matrix_valued now = tensor_valued<2>
    (simplifies generic code for h5).
  - factories_one_var -> factories : this is the generic case ...
    only a few specialization, code is simpler.

- clef expression call with rvalue for *this
- generalize matrix_proxy to tensor and clean

 - clean exception catch in tests

  - exception catching catch in need in test
    because the silly OS X does not print anything, just "exception occurred".
    Very convenient for the developer...
  - BUT, one MUST add return 1, or the make test will *pass* !!
  - --> systematically replace the catch by a macro TRIQS_CATCH_AND_ABORT
    which return a non zero error code.
   - exception : curry_and_fourier which does not work at this stage
   (mesh incompatible).

- gf: clean draft of gf 2 times
  - comment the python interface for the moment.
  - rm useless tests
2013-10-21 15:11:44 +02:00

85 lines
4.1 KiB
C++

/*******************************************************************************
*
* 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 <http://www.gnu.org/licenses/>.
*
******************************************************************************/
#ifndef TRIQS_GF_LOCAL_FOURIER_MATSU_H
#define TRIQS_GF_LOCAL_FOURIER_MATSU_H
#include "fourier_base.hpp"
#include <triqs/gfs/imfreq.hpp>
#include <triqs/gfs/imtime.hpp>
namespace triqs { namespace gfs {
// First the implementation of the fourier transform
void fourier_impl (gf_view<imfreq,scalar_valued> gw , gf_view<imtime,scalar_valued> const gt, scalar_valued);
void fourier_impl (gf_view<imfreq,matrix_valued> gw , gf_view<imtime,matrix_valued> const gt, matrix_valued);
void inverse_fourier_impl (gf_view<imtime,scalar_valued> gt, gf_view<imfreq,scalar_valued> const gw, scalar_valued);
void inverse_fourier_impl (gf_view<imtime,matrix_valued> gt, gf_view<imfreq,matrix_valued> const gw, matrix_valued);
inline gf_view<imfreq,matrix_valued> fourier (gf_view<imtime,matrix_valued> const gt) {
int L = (gt.mesh().kind() == full_bins ? gt.mesh().size()-1 : gt.mesh().size() );
auto gw = gf<imfreq,matrix_valued>{ {gt.domain(),L}, gt.data().shape().front_pop() };
auto V = gw();
fourier_impl(V, gt, matrix_valued());
return gw;
}
inline gf_view<imfreq,scalar_valued> fourier (gf_view<imtime,scalar_valued> const gt) {
int L = (gt.mesh().kind() == full_bins ? gt.mesh().size()-1 : gt.mesh().size() );
auto gw = gf<imfreq,scalar_valued>{ {gt.domain(),L} };
auto V = gw();
fourier_impl(V, gt, scalar_valued());
return gw;
}
inline gf_view<imtime, matrix_valued> inverse_fourier (gf_view<imfreq, matrix_valued> const gw, mesh_kind mk = half_bins) {
double pi = std::acos(-1);
int L = (mk == full_bins ? gw.mesh().size()+1 : gw.mesh().size() );
auto gt = gf<imtime,matrix_valued>{ {gw.domain(),L}, gw.data().shape().front_pop()};
auto V = gt();
inverse_fourier_impl(V, gw, matrix_valued());
return gt;
}
inline gf_view<imtime,scalar_valued> inverse_fourier (gf_view<imfreq,scalar_valued> const gw, mesh_kind mk = half_bins) {
double pi = std::acos(-1);
int L = (mk == full_bins ? gw.mesh().size()+1 : gw.mesh().size() );
auto gt = gf<imtime,scalar_valued>{ {gw.domain(),L} };
auto V = gt();
inverse_fourier_impl(V, gw,scalar_valued());
return gt;
}
inline gf_keeper<tags::fourier,imtime,scalar_valued> lazy_fourier (gf_view<imtime,scalar_valued> const & g) { return g;}
inline gf_keeper<tags::fourier,imfreq,scalar_valued> lazy_inverse_fourier (gf_view<imfreq,scalar_valued> const & g) { return g;}
inline gf_keeper<tags::fourier,imtime,matrix_valued> lazy_fourier (gf_view<imtime,matrix_valued> const & g) { return g;}
inline gf_keeper<tags::fourier,imfreq,matrix_valued> lazy_inverse_fourier (gf_view<imfreq,matrix_valued> const & g) { return g;}
void triqs_gf_view_assign_delegation( gf_view<imfreq,scalar_valued> g, gf_keeper<tags::fourier,imtime,scalar_valued> const & L);
void triqs_gf_view_assign_delegation( gf_view<imfreq,matrix_valued> g, gf_keeper<tags::fourier,imtime,matrix_valued> const & L);
void triqs_gf_view_assign_delegation( gf_view<imtime,scalar_valued> g, gf_keeper<tags::fourier,imfreq,scalar_valued> const & L);
void triqs_gf_view_assign_delegation( gf_view<imtime,matrix_valued> g, gf_keeper<tags::fourier,imfreq,matrix_valued> const & L);
}}
namespace triqs { namespace clef {
TRIQS_CLEF_MAKE_FNT_LAZY (lazy_fourier);
TRIQS_CLEF_MAKE_FNT_LAZY (lazy_inverse_fourier);
}}
#endif