예제 #1
0
//fix phase of eigensystem and store phase of first entry of each eigenvector
void fix_phase(Eigen::MatrixXcd& V, Eigen::MatrixXcd& V_fix, std::vector<double>& phase) {
   const int V3 = pars -> get_int("V3");
  //helper variables:
  //Number of eigenvectors
  int n_ev;
  //negative imaginary
  std::complex<double> i_neg (0,-1);
  //tmp factor and phase
  std::complex<double> fac (1.,1.);
  double tmp_phase = 0;
  //get sizes right, resize if necessary
  n_ev = V.cols();
  if (phase.size() != n_ev) phase.resize(n_ev);
  if (V_fix.cols() != n_ev) V_fix.resize(3*V3,n_ev);
  //loop over all eigenvectors of system
  for (int n = 0; n < n_ev; ++n) {

    tmp_phase = std::arg(V(0,n));
    phase.at(n) = tmp_phase;
    fac = std::exp(i_neg*tmp_phase);
    //Fix phase of eigenvector with negative polar angle of first entry
    V_fix.col(n) = fac * V.col(n); 

  }
}
예제 #2
0
//Reads in Eigenvectors from one Timeslice in binary format to V
void read_evectors_bin_ts(const char * path, const char* prefix, const int config_i, const int t,
    const int nb_ev, Eigen::MatrixXcd& V) {
  int V3 = pars -> get_int("V3");
  //bool thorough = pars -> get_int("strict");
  const int dim_row = 3 * V3;
  //TODO: Change path getting to something keyword independent
  //buffer for read in
  std::complex<double>* eigen_vec = new std::complex<double>[dim_row];
  //setting up file
  char filename[200];
  sprintf(filename, "%s/%s.%04d.%03d", path, prefix, config_i, t);
  std::cout << "Reading file: " << filename << std::endl;
  std::ifstream infile(filename, std::ifstream::binary);
  for (int nev = 0; nev < nb_ev; ++nev) {
    infile.read( reinterpret_cast<char*> (eigen_vec), 2*dim_row*sizeof(double));
    V.col(nev) = Eigen::Map<Eigen::VectorXcd, 0 >(eigen_vec, dim_row);
    eof_check(t,nev,nb_ev,infile.eof());
  }
  if(check_trace(V, nb_ev) != true){
    std::cout << "Timeslice: " << t << ": Eigenvectors damaged, exiting" << std::endl;
    exit(0);
  }

  //clean up
  delete[] eigen_vec;
  infile.close();
}
예제 #3
0
//TODO complex values as matrix entries too
//TODO support endomorphisms over Rn too
Expression EigenvectorsCommand::operator()(const QList< Analitza::Expression >& args)
{
    Expression ret;
    
    QStringList errors;
    const Eigen::MatrixXcd eigeninfo = executeEigenSolver(args, true, errors);
    
    if (!errors.isEmpty()) {
        ret.addError(errors.first());
        
        return ret;
    }
    const int neigenvectors = eigeninfo.rows();
    
    QScopedPointer<Analitza::List> list(new Analitza::List);
    
    for (int j = 0; j < neigenvectors; ++j) {
        const Eigen::VectorXcd col = eigeninfo.col(j);
        QScopedPointer<Analitza::Vector> eigenvector(new Analitza::Vector(neigenvectors));
        
        for (int i = 0; i < neigenvectors; ++i) {
            const std::complex<double> eigenvalue = col(i);
            const double realpart = eigenvalue.real();
            const double imagpart = eigenvalue.imag();
            
            if (std::isnan(realpart) || std::isnan(imagpart)) {
                ret.addError(QCoreApplication::tr("Returned eigenvalue is NaN", "NaN means Not a Number, is an invalid float number"));
                return ret;
            } else if (std::isinf(realpart) || std::isinf(imagpart)) {
                ret.addError(QCoreApplication::tr("Returned eigenvalue is too big"));
                return ret;
            } else {
                bool isonlyreal = true;
                
                if (std::isnormal(imagpart)) {
                    isonlyreal = false;
                }
                
                Analitza::Cn * eigenvalueobj = 0;
                
                if (isonlyreal) {
                    eigenvalueobj = new Analitza::Cn(realpart);
                } else {
                    eigenvalueobj = new Analitza::Cn(realpart, imagpart);
                }
                
                eigenvector->appendBranch(eigenvalueobj);
            }
        }
        
        list->appendBranch(eigenvector.take());
    }
    
    ret.setTree(list.take());
    
    return ret;
}
예제 #4
0
// TODO: work on interface with eigenvector class
// transform matrix of eigenvectors with gauge array
Eigen::MatrixXcd GaugeField::trafo_ev(const Eigen::MatrixXcd &eig_sys) {
  const ssize_t dim_row = eig_sys.rows();
  const ssize_t dim_col = eig_sys.cols();
  Eigen::MatrixXcd ret(dim_row, dim_col);
  if (omega.shape()[0] == 0)
    build_gauge_array(1);
  // write_gauge_matrices("ev_trafo_log.bin",Omega);

  for (ssize_t nev = 0; nev < dim_col; ++nev) {
    for (ssize_t vol = 0; vol < dim_row; ++vol) {
      int ind_c = vol % 3;
      Eigen::Vector3cd tmp =
          omega[0][ind_c].adjoint() * (eig_sys.col(nev)).segment(ind_c, 3);
      (ret.col(nev)).segment(ind_c, 3) = tmp;
    }  // end loop nev
  }    // end loop vol
  return ret;
}
예제 #5
0
// -----------------------------------------------------------------------------
// -----------------------------------------------------------------------------
void LapH::OperatorsForMesons::build_rvdaggervr(
                                            const LapH::RandomVector& rnd_vec) {

  // check of vdaggerv is already build
  if(not is_vdaggerv_set){
    std::cout << "\n\n\tCaution: vdaggerv is not set and rvdaggervr cannot be" 
              << " computed\n\n" << std::endl;
    exit(0);
  }

  clock_t t2 = clock();
  std::cout << "\tbuild rvdaggervr:";

  for(auto& rvdvr_level1 : rvdaggervr)
    for(auto& rvdvr_level2 : rvdvr_level1)
      for(auto& rvdvr_level3 : rvdvr_level2)
        rvdvr_level3 = Eigen::MatrixXcd::Zero(4*dilE, 4*dilE);

#pragma omp parallel for schedule(dynamic)
  for(size_t t = 0; t < Lt; t++){

  // rvdaggervr is calculated by multiplying vdaggerv with the same quantum
  // numbers with random vectors from right and left.
  for(const auto& op : operator_lookuptable.rvdaggervr_lookuptable){

    Eigen::MatrixXcd vdv;
    if(op.need_vdaggerv_daggering == false)
      vdv = vdaggerv[op.id_vdaggerv][t];
    else
      vdv = vdaggerv[op.id_vdaggerv][t].adjoint();

    size_t rid = 0;
    int check = -1;
    Eigen::MatrixXcd M; // Intermediate memory
    for(const auto& rnd_id : 
              operator_lookuptable.ricQ2_lookup[op.id_ricQ_lookup].rnd_vec_ids){

      if(check != rnd_id.first){ // this avoids recomputation
        M = Eigen::MatrixXcd::Zero(nb_ev, 4*dilE);
        for(size_t block = 0; block < 4; block++){
        for(size_t vec_i = 0; vec_i < nb_ev; vec_i++) {
          size_t blk =  block + (vec_i + nb_ev * t) * 4;
          M.block(0, vec_i%dilE + dilE*block, nb_ev, 1) += 
               vdv.col(vec_i) * rnd_vec(rnd_id.first, blk);
        }}
      }
      for(size_t block_x = 0; block_x < 4; block_x++){
      for(size_t block_y = 0; block_y < 4; block_y++){
      for(size_t vec_y = 0; vec_y < nb_ev; ++vec_y) {
        size_t blk =  block_y + (vec_y + nb_ev * t) * 4;
        rvdaggervr[op.id][t][rid].block(
                            dilE*block_y + vec_y%dilE, dilE*block_x, 1, dilE) +=
                M.block(vec_y, dilE*block_x, 1, dilE) * 
                std::conj(rnd_vec(rnd_id.second, blk));
      }}} 
      check = rnd_id.first;
      rid++;
    }
  }}// time and operator loops end here

  t2 = clock() - t2;
  std::cout << std::setprecision(1) << "\t\tSUCCESS - " << std::fixed 
    << ((float) t2)/CLOCKS_PER_SEC << " seconds" << std::endl;
}
예제 #6
0
// -----------------------------------------------------------------------------
// -----------------------------------------------------------------------------
void LapH::distillery::create_source(const size_t dil_t, const size_t dil_e,
                                     std::complex<double>* source) {
  
  if(dil_t >= param.dilution_size_so[0] ||
     dil_e >= param.dilution_size_so[1] ){
        std::cout << "dilution is out of bounds in \"create_source\"" <<
        std::endl;
        exit(0);
  }

  const size_t Lt = param.Lt;
  const size_t Ls = param.Ls;
  const size_t number_of_eigen_vec = param.nb_ev;
  const size_t Vs = Ls*Ls*Ls;
  const size_t dim_row = Vs*3;

  int numprocs = 0, myid = 0;
  MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
  MPI_Comm_rank(MPI_COMM_WORLD, &myid);

  // setting source to zero 
  for(size_t dil_d = 0; dil_d < param.dilution_size_so[2]; ++dil_d)
    for(size_t i = 0; i < 12*Vs*Lt/numprocs; ++i)
      source[dil_d*12*Vs*Lt/numprocs + i] = {0.0, 0.0};

  // indices of timeslices with non-zero entries
  size_t nb_of_nonzero_t = Lt/param.dilution_size_so[0];  
  // TODO: think about smart pointer here!
  size_t* t_index = new size_t[nb_of_nonzero_t];
  create_dilution_lookup(nb_of_nonzero_t, param.dilution_size_so[0], 
                         dil_t, param.dilution_type_so[0], t_index);

  // indices of eigenvectors to be combined
  size_t nb_of_ev_combined = number_of_eigen_vec/param.dilution_size_so[1];
  size_t* ev_index = new size_t[nb_of_ev_combined];
  create_dilution_lookup(nb_of_ev_combined, param.dilution_size_so[1], 
                         dil_e, param.dilution_type_so[1], ev_index);

  // indices of Dirac components to be combined 
  // TODO: This is needed only once - could be part of class
  size_t nb_of_dirac_combined = 4/param.dilution_size_so[2];
  size_t** d_index = new size_t*[param.dilution_size_so[2]];
  for(size_t i = 0; i < param.dilution_size_so[2]; ++i)
    d_index[i] = new size_t[nb_of_dirac_combined];
  for(size_t dil_d = 0; dil_d < param.dilution_size_so[2]; ++dil_d)
    create_dilution_lookup(nb_of_dirac_combined, param.dilution_size_so[2],
                           dil_d, param.dilution_type_so[2], d_index[dil_d]);

  // creating the source 
  // running over nonzero timeslices
  for(size_t t = 0; t < nb_of_nonzero_t; ++t){ 
    // intermidiate memory
    Eigen::MatrixXcd S = Eigen::MatrixXcd::Zero(dim_row, 4);
    size_t time = t_index[t];                 // helper index
    if((int) (time / (Lt/numprocs)) == myid) {
      size_t time_id = time % (Lt/numprocs);  // helper index
      // building source on one timeslice
      for(size_t ev = 0; ev < nb_of_ev_combined; ++ev){
        size_t ev_h = ev_index[ev] * 4; // helper index 
        for(size_t d = 0; d < 4; ++d){
          S.col(d) += random_vector[time](ev_h+d) *
                          (V[time_id]).col(ev_index[ev]); 
        }
      }
      // copy the created source into output array
      size_t t_h = time_id*Vs; // helper index
      for(size_t x = 0; x < Vs; ++x){
        size_t x_h  = (t_h + x)*12;    // helper index
        size_t x_h2 = x*3;            // helper index
        for(size_t d2 = 0; d2 < param.dilution_size_so[2]; ++d2){
          for(size_t d3 = 0; d3 < nb_of_dirac_combined; ++d3){
            size_t d_h = x_h + d_index[d2][d3]*3; // helper index
            for(size_t c = 0; c < 3; ++c){
              source[d2*12*Vs*Lt/numprocs + d_h + c] = 
                                          S(x_h2 + c, d_index[d2][d3]);
            }
          }
        }    
      }
    }
  } // end of loop over nonzero timeslices
  MPI_Barrier(MPI_COMM_WORLD);
  // freeing memory
  delete[] t_index;
  delete[] ev_index;
  for(size_t i = 0; i < param.dilution_size_so[2]; ++i)
    delete[] d_index[i];
  delete[] d_index;
  t_index = NULL;
  ev_index = NULL;
  d_index = NULL;
}