Esempio n. 1
0
void FieldOperator::compute(const boost::mpi::communicator& comm)
{
    if (Status < Prepared) throw (exStatusMismatch());
    if (Status >= Computed) return;

    if (!comm.rank()) INFO_NONEWLINE("Computing " << *O << " in eigenbasis of the Hamiltonian: ");
/*

    pMPI::mpi_skel<pMPI::ComputeWrap<FieldOperatorPart>> skel;
    skel.parts.resize(parts.size());
    for (size_t i=0; i<parts.size(); i++) { skel.parts[i] = pMPI::ComputeWrap<FieldOperatorPart>(*parts[i]);};
    std::map<pMPI::JobId, pMPI::WorkerId> job_map = skel.run(comm, true); // actual running - very costly

    int rank = comm.rank();
    int comm_size = comm.size(); 

    // Start distributing data
    comm.barrier();
         
    for (size_t p = 0; p<parts.size(); p++) {

        int nelem = 0;
        if (rank == job_map[p]) nelem = parts[p]->elementsColMajor.nonZeros(); 
        boost::mpi::broadcast(comm, nelem, job_map[p]);
        auto data = parts[p]->elementsColMajor.data(); // Eigen::internal::CompressedStorage 
        auto& value_start = data.value(0);
        auto& index_start = data.index(0);
        if (comm.rank()) parts[p]->elementsColMajor.resize(nelem);
        exit(0);

        //boost::mpi::broadcast(comm, parts[p]->elementsColMajor.data(), nelem, job_map[p]);
        if (rank == job_map[p]) { 
                parts[p]->Status = FieldOperatorPart::Computed;
             };
        };
    comm.barrier();
*/
    size_t Size = parts.size();
    for (size_t BlockIn = 0; BlockIn < Size; BlockIn++){
        INFO_NONEWLINE( (int) ((1.0*BlockIn/Size) * 100 ) << "  " << std::flush);
        parts[BlockIn]->compute();
    };
    INFO("");
    Status = Computed;
}
Esempio n. 2
0
void TwoParticleGFContainer::computeAll_split(bool clearTerms, const boost::mpi::communicator & comm)
{
    // split communicator
    size_t ncomponents = NonTrivialElements.size();
    size_t ncolors = std::min(int(comm.size()), int(NonTrivialElements.size()));
    RealType color_size = 1.0*comm.size()/ncolors;
    std::map<int,int> proc_colors;
    std::map<int,int> elem_colors;
    std::map<int,int> color_roots;
    bool calc = false;
    for (size_t p=0; p<comm.size(); p++) {
        int color = int (1.0*p / color_size);
        proc_colors[p] = color;
        color_roots[color]=p;
    }
    for (size_t i=0; i<ncomponents; i++) {
        int color = i*ncolors/ncomponents;
        elem_colors[i] = color;
    };

    if (!comm.rank()) {
        INFO("Splitting " << ncomponents << " components in " << ncolors << " communicators");
        for (size_t i=0; i<ncomponents; i++)
            INFO("2pgf " << i << " color: " << elem_colors[i] << " color_root: " << color_roots[elem_colors[i]]);
    };
    comm.barrier();
    int comp = 0;

    boost::mpi::communicator comm_split = comm.split(proc_colors[comm.rank()]);

    for(std::map<IndexCombination4, boost::shared_ptr<TwoParticleGF> >::iterator iter = NonTrivialElements.begin(); iter != NonTrivialElements.end(); iter++, comp++) {
        bool calc = (elem_colors[comp] == proc_colors[comm.rank()]);
        if (calc) {
            INFO("C" << elem_colors[comp] << "p" << comm.rank() << ": computing 2PGF for " << iter->first);
            if (calc) static_cast<TwoParticleGF&>(*(iter->second)).compute(clearTerms, comm_split);
        };
    };
    comm.barrier();
    // distribute data
    if (!comm.rank()) INFO_NONEWLINE("Distributing 2PGF container...");
    comp = 0;
    for(std::map<IndexCombination4, boost::shared_ptr<TwoParticleGF> >::iterator iter = NonTrivialElements.begin(); iter != NonTrivialElements.end(); iter++, comp++) {
        int sender = color_roots[elem_colors[comp]];
        TwoParticleGF& chi = *((iter)->second);
        for (size_t p = 0; p<chi.parts.size(); p++) {
            //    if (comm.rank() == sender) INFO("P" << comm.rank() << " 2pgf " << p << " " << chi.parts[p]->NonResonantTerms.size());
            boost::mpi::broadcast(comm, chi.parts[p]->NonResonantTerms, sender);
            boost::mpi::broadcast(comm, chi.parts[p]->ResonantTerms, sender);
            if (comm.rank() != sender) {
                chi.setStatus(TwoParticleGF::Computed);
            };
        };
    }
    comm.barrier();
    if (!comm.rank()) INFO("done.");
}
Esempio n. 3
0
template <bool F> matrix_type BetheSalpeter<matrix_type, F>::solve_iterations(size_t n_iter, real_type mix, bool evaluate_only_order_n)
{
    matrix_type V4     = vertex_;
    matrix_type V4_old = vertex_;
    matrix_type V4Chi  = vertex_*bubble_;
    if (verbosity_) INFO_NONEWLINE("\tEvaluating BS self-consistently. Making " << n_iter << " iterations.");
    real_type diffBS = 1.0;
    for (size_t n=0; n<n_iter && diffBS > 1e-8 * double(!evaluate_only_order_n); ++n) { 
        if (verbosity_ > 1) INFO_NONEWLINE("\t\t" << n+1 << "/" << n_iter<< ". ")
        if (fwd)
            V4 = (double(n==n_iter - 1 || !evaluate_only_order_n) * vertex_ + V4Chi*V4_old)*mix + (1.0-mix)*V4_old;
        else 
            V4 = (double(n==n_iter - 1 || !evaluate_only_order_n) * vertex_ - V4_old*V4Chi)*mix + (1.0-mix)*V4_old;
        diffBS = (V4-V4_old).norm();
        if (verbosity_ > 1) INFO("vertex diff = " << diffBS);
        V4_old = V4;
        }
    return std::move(V4);
}
Esempio n. 4
0
template <bool F> matrix_type BetheSalpeter<matrix_type, F>::solve_inversion()
{
    if (verbosity_ > 0) INFO_NONEWLINE("Running" << ((!fwd)?" inverse ":" ") << "matrix BS equation...");
    size_t size = vertex_.rows(); 
    matrix_type V4Chi = fwd ? matrix_type(matrix_type::Identity(size,size) - vertex_*bubble_) : matrix_type(matrix_type::Identity(size,size) + bubble_*vertex_);

    Eigen::PartialPivLU<matrix_type> Solver(V4Chi);
    det_ = Solver.determinant(); 

    assert(is_float_equal(det_, V4Chi.determinant(), 1e-6));

    if (std::abs(std::imag(det_))>1e-2 * std::abs(std::real(det_))) { ERROR("Determinant : " << det_); throw (std::logic_error("Complex determinant in BS. Exiting.")); };
    if ((std::real(det_))<1e-2) INFO3("Determinant : " << det_);

    if (std::real(det_) < std::numeric_limits<real_type>::epsilon()) {
        ERROR("Can't solve Bethe-Salpeter equation by inversion");
        return std::move(V4Chi);
    }

    V4Chi = fwd ? matrix_type(Solver.solve(vertex_)) : matrix_type(vertex_*Solver.inverse());
                //V4Chi=vertex_*Solver.inverse();
    if (verbosity_ > 0) INFO("done.");
    return std::move(V4Chi);
}
int main(int argc, char const *argv[]) {
  Eigen::setNbThreads(NumCores);
#ifdef MKL
  mkl_set_num_threads(NumCores);
#endif
  INFO("Eigen3 uses " << Eigen::nbThreads() << " threads.");
  int L;
  RealType J12ratio;
  int OBC;
  int N;
  RealType Uin, phi;
  std::vector<RealType> Vin;
  LoadParameters( "conf.h5", L, J12ratio, OBC, N, Uin, Vin, phi);
  HDF5IO file("BSSH.h5");
  // const int L = 5;
  // const bool OBC = true;
  // const RealType J12ratio = 0.010e0;
  INFO("Build Lattice - ");
  std::vector<ComplexType> J;
  if ( OBC ){
    J = std::vector<ComplexType>(L - 1, ComplexType(1.0, 0.0));
    for (size_t cnt = 0; cnt < L-1; cnt+=2) {
      J.at(cnt) *= J12ratio;
    }
  } else{
    J = std::vector<ComplexType>(L, ComplexType(1.0, 0.0));
    for (size_t cnt = 0; cnt < L; cnt+=2) {
      J.at(cnt) *= J12ratio;
    }
    if ( std::abs(phi) > 1.0e-10 ){
      J.at(L-1) *= exp( ComplexType(0.0e0, 1.0e0) * phi );
      // INFO(exp( ComplexType(0.0e0, 1.0e0) * phi ));
    }
  }
  for ( auto &val : J ){
    INFO_NONEWLINE(val << " ");
  }
  INFO("");
  const std::vector< Node<ComplexType>* > lattice = NN_1D_Chain(L, J, OBC);
  file.saveNumber("1DChain", "L", L);
  file.saveNumber("1DChain", "U", Uin);
  file.saveStdVector("1DChain", "J", J);
  for ( auto &lt : lattice ){
    if ( !(lt->VerifySite()) ) RUNTIME_ERROR("Wrong lattice setup!");
  }
  INFO("DONE!");
  INFO("Build Basis - ");
  // int N1 = (L+1)/2;
  Basis B1(L, N);
  B1.Boson();
  // std::vector< std::vector<int> > st = B1.getBStates();
  // std::vector< RealType > tg = B1.getBTags();
  // for (size_t cnt = 0; cnt < tg.size(); cnt++) {
  //   INFO_NONEWLINE( std::setw(3) << cnt << " - ");
  //   for (auto &j : st.at(cnt)){
  //     INFO_NONEWLINE(j << " ");
  //   }
  //   INFO("- " << tg.at(cnt));
  // }
  file.saveNumber("1DChain", "N", N);
  // file.saveStdVector("Basis", "States", st);
  // file.saveStdVector("Basis", "Tags", tg);
  INFO("DONE!");
  INFO_NONEWLINE("Build Hamiltonian - ");
  std::vector<Basis> Bases;
  Bases.push_back(B1);
  Hamiltonian<ComplexType> ham( Bases );
  std::vector< std::vector<ComplexType> > Vloc;
  std::vector<ComplexType> Vtmp;//(L, 1.0);
  for ( RealType &val : Vin ){
    Vtmp.push_back((ComplexType)val);
  }
  Vloc.push_back(Vtmp);
  std::vector< std::vector<ComplexType> > Uloc;
  // std::vector<ComplexType> Utmp(L, ComplexType(10.0e0, 0.0e0) );
  std::vector<ComplexType> Utmp(L, (ComplexType)Uin);
  Uloc.push_back(Utmp);
  ham.BuildLocalHamiltonian(Vloc, Uloc, Bases);
  ham.BuildHoppingHamiltonian(Bases, lattice);
  ham.BuildTotalHamiltonian();
  INFO("DONE!");
  INFO_NONEWLINE("Diagonalize Hamiltonian - ");
  std::vector<RealType> Val;
  Hamiltonian<ComplexType>::VectorType Vec;
  ham.eigh(Val, Vec);
  INFO("GS energy = " << Val.at(0));
  file.saveVector("GS", "EVec", Vec);
  file.saveStdVector("GS", "EVal", Val);
  INFO("DONE!");
  std::vector<ComplexType> Nbi = Ni( Bases, Vec );
  for (auto &n : Nbi ){
    INFO( n << " " );
  }
  ComplexMatrixType Nij = NiNj( Bases, Vec );
  INFO(Nij);
  INFO(Nij.diagonal());
  file.saveStdVector("Obs", "Nb", Nbi);
  file.saveMatrix("Obs", "Nij", Nij);
  return 0;
}
Esempio n. 6
0
void TwoParticleGFPart::compute()
{
    NonResonantTerms.clear();
    ResonantTerms.clear();

    RealType beta = DMpart1.beta;
    // I don't have any pen now, so I'm writing here:
    // <1 | O1 | 2> <2 | O2 | 3> <3 | O3 |4> <4| CX4 |1>
    // Iterate over all values of |1><1| and |3><3|
    // Chase indices |2> and <2|, |4> and <4|.
    const RowMajorMatrixType& O1matrix = O1.getRowMajorValue();
    const ColMajorMatrixType& O2matrix = O2.getColMajorValue();    
    const RowMajorMatrixType& O3matrix = O3.getRowMajorValue();
    const ColMajorMatrixType& CX4matrix = CX4.getColMajorValue();

    InnerQuantumState index1;
    InnerQuantumState index1Max = CX4matrix.outerSize(); // One can not make a cutoff in external index for evaluating 2PGF 

    InnerQuantumState index3;
    InnerQuantumState index3Max = O2matrix.outerSize();

    std::list<InnerQuantumState> Index4List;

    unsigned long ResonantTermsUnreducedSize=0;
    unsigned long NonResonantTermsUnreducedSize=0;
    unsigned long ResonantTermsPreviousSize=0;
    unsigned long NonResonantTermsPreviousSize=0;

    for(index1=0; index1<index1Max; ++index1){
    for(index3=0; index3<index3Max; ++index3){
        ColMajorMatrixType::InnerIterator index4bra_iter(CX4matrix,index1);       
        RowMajorMatrixType::InnerIterator index4ket_iter(O3matrix,index3);
        Index4List.clear();
        while (index4bra_iter && index4ket_iter){
            if(chaseIndices(index4ket_iter,index4bra_iter)){
                Index4List.push_back(index4bra_iter.index());
                ++index4bra_iter;
                ++index4ket_iter;
            }
        };

        if (!Index4List.empty())
        {
            RealType E1 = Hpart1.getEigenValue(index1);
            RealType E3 = Hpart3.getEigenValue(index3);
            RealType weight1 = DMpart1.getWeight(index1);
            RealType weight3 = DMpart3.getWeight(index3);

            ColMajorMatrixType::InnerIterator index2bra_iter(O2matrix,index3);
            RowMajorMatrixType::InnerIterator index2ket_iter(O1matrix,index1);       
            while (index2bra_iter && index2ket_iter){
                if (chaseIndices(index2ket_iter,index2bra_iter)){

                    InnerQuantumState index2 = index2ket_iter.index();
                    RealType E2 = Hpart2.getEigenValue(index2);
                    RealType weight2 = DMpart2.getWeight(index2);

                    for (std::list<InnerQuantumState>::iterator pIndex4 = Index4List.begin(); pIndex4!=Index4List.end(); ++pIndex4) 
                    {
                        InnerQuantumState index4 = *pIndex4;
                        RealType E4 = Hpart4.getEigenValue(index4);                       
                        RealType weight4 = DMpart4.getWeight(index4);
                        if (weight1 + weight2 + weight3 + weight4 < CoefficientTolerance) break;  // optimization 

                        ComplexType MatrixElement = index2ket_iter.value()*
                                                    index2bra_iter.value()*
                                                    O3matrix.coeff(index3,index4)*
                                                    CX4matrix.coeff(index4,index1);

                        MatrixElement *= Permutation.sign;

                        addMultiterm(MatrixElement,beta,E1,E2,E3,E4,weight1,weight2,weight3,weight4);
                    }
                    ++index2bra_iter;
                    ++index2ket_iter;
                };
            }
        };
    }
    if (ResonantTerms.size()-ResonantTermsPreviousSize + NonResonantTerms.size() - NonResonantTermsPreviousSize > ReduceInvocationThreshold ){ 
        INFO_NONEWLINE(NonResonantTerms.size()-NonResonantTermsPreviousSize << "+" << ResonantTerms.size() - ResonantTermsPreviousSize << "=");
        INFO_NONEWLINE(ResonantTerms.size()-ResonantTermsPreviousSize + NonResonantTerms.size() - NonResonantTermsPreviousSize);
        INFO_NONEWLINE(" terms -> ");

        NonResonantTermsUnreducedSize+=NonResonantTerms.size();
        ResonantTermsUnreducedSize+= ResonantTerms.size();

        RealType NonResonantTolerance = MultiTermCoefficientTolerance*(index1+1)/NonResonantTermsUnreducedSize/(index1Max+1);
        RealType ResonantTolerance = MultiTermCoefficientTolerance*(index1+1)/ResonantTermsUnreducedSize/(index1Max+1);

        this->reduceTerms(NonResonantTolerance, ResonantTolerance, NonResonantTerms, ResonantTerms); 
        NonResonantTermsPreviousSize = NonResonantTerms.size();
        ResonantTermsPreviousSize = ResonantTerms.size(); 

        INFO_NONEWLINE(NonResonantTermsPreviousSize << "+" << ResonantTermsPreviousSize << "=");
        INFO(NonResonantTermsPreviousSize + ResonantTermsPreviousSize << " , tol = " << std::max(NonResonantTolerance,ResonantTolerance));
        };
    };

    NonResonantTermsUnreducedSize=(NonResonantTermsUnreducedSize>0)?NonResonantTermsUnreducedSize:NonResonantTerms.size();
    ResonantTermsUnreducedSize=(ResonantTermsUnreducedSize>0)?ResonantTermsUnreducedSize:ResonantTerms.size();
    if (ResonantTermsUnreducedSize + NonResonantTermsUnreducedSize > 0){
        INFO_NONEWLINE("Total " << NonResonantTermsUnreducedSize << "+" << ResonantTermsUnreducedSize << "=");
        INFO_NONEWLINE(NonResonantTermsUnreducedSize+ResonantTermsUnreducedSize << " terms -> ");
        this->reduceTerms(MultiTermCoefficientTolerance/(NonResonantTermsUnreducedSize+1), MultiTermCoefficientTolerance/(ResonantTermsUnreducedSize+1), NonResonantTerms, ResonantTerms);
        INFO_NONEWLINE(NonResonantTerms.size() << "+" << ResonantTerms.size() << "=");
        INFO(NonResonantTerms.size() + ResonantTerms.size()  << ", \ttols = " << std::setw(4)  
             << std::max(MultiTermCoefficientTolerance/(NonResonantTermsUnreducedSize+1), MultiTermCoefficientTolerance/(ResonantTermsUnreducedSize+1))
             << " (coeff), " << ReduceResonanceTolerance << " (res)" );
    };

    Status = Computed;
}
int main(int argc, char const *argv[]) {
  Eigen::setNbThreads(NumCores);
#ifdef MKL
  mkl_set_num_threads(NumCores);
#endif
  INFO("Eigen3 uses " << Eigen::nbThreads() << " threads.");
  int L;
  RealType J12ratio;
  int OBC;
  int N1, N2;
  int dynamics, Tsteps;
  RealType Uin, phi, dt;
  std::vector<RealType> Vin;
  LoadParameters( "conf.h5", L, J12ratio, OBC, N1, N2, Uin, Vin, phi, dynamics, Tsteps, dt);
  HDF5IO *file = new HDF5IO("FSSH.h5");
  // const int L = 5;
  // const bool OBC = true;
  // const RealType J12ratio = 0.010e0;
  INFO("Build Lattice - ");
  std::vector<DT> J;
  if ( OBC ){
    // J = std::vector<DT>(L - 1, 0.0);// NOTE: Atomic limit testing
    J = std::vector<DT>(L - 1, 1.0);
    if ( J12ratio > 1.0e0 ) {
      for (size_t cnt = 1; cnt < L-1; cnt+=2) {
        J.at(cnt) /= J12ratio;
      }
    } else{
      for (size_t cnt = 0; cnt < L-1; cnt+=2) {
        J.at(cnt) *= J12ratio;
      }
    }
  } else{
    J = std::vector<DT>(L, 1.0);
    if ( J12ratio > 1.0e0 ) {
      for (size_t cnt = 1; cnt < L; cnt+=2) {
        J.at(cnt) /= J12ratio;
      }
    } else{
      for (size_t cnt = 0; cnt < L; cnt+=2) {
        J.at(cnt) *= J12ratio;
      }
    }
#ifndef DTYPE
    if ( std::abs(phi) > 1.0e-10 ){
      J.at(L-1) *= exp( DT(0.0, 1.0) * phi );
    }
#endif
  }
  for ( auto &val : J ){
    INFO_NONEWLINE(val << " ");
  }
  INFO("");
  const std::vector< Node<DT>* > lattice = NN_1D_Chain(L, J, OBC);
  file->saveNumber("1DChain", "L", L);
  file->saveStdVector("1DChain", "J", J);
  for ( auto &lt : lattice ){
    if ( !(lt->VerifySite()) ) RUNTIME_ERROR("Wrong lattice setup!");
  }
  INFO("DONE!");
  INFO("Build Basis - ");
  // int N1 = (L+1)/2;
  Basis F1(L, N1, true);
  F1.Fermion();
  std::vector<int> st1 = F1.getFStates();
  std::vector<size_t> tg1 = F1.getFTags();
  // for (size_t cnt = 0; cnt < st1.size(); cnt++) {
  //   INFO_NONEWLINE( std::setw(3) << st1.at(cnt) << " - ");
  //   F1.printFermionBasis(st1.at(cnt));
  //   INFO("- " << tg1.at(st1.at(cnt)));
  // }
  // int N2 = (L-1)/2;
  Basis F2(L, N2, true);
  F2.Fermion();
  std::vector<int> st2 = F2.getFStates();
  std::vector<size_t> tg2 = F2.getFTags();
  // for (size_t cnt = 0; cnt < st2.size(); cnt++) {
  //   INFO_NONEWLINE( std::setw(3) << st2.at(cnt) << " - ");
  //   F2.printFermionBasis(st2.at(cnt));
  //   INFO("- " << tg2.at(st2.at(cnt)));
  // }
  file->saveNumber("Basis", "N1", N1);
  file->saveStdVector("Basis", "F1States", st1);
  file->saveStdVector("Basis", "F1Tags", tg1);
  file->saveNumber("Basis", "N2", N2);
  file->saveStdVector("Basis", "F2States", st2);
  file->saveStdVector("Basis", "F2Tags", tg2);
  INFO("DONE!");
  INFO_NONEWLINE("Build Hamiltonian - ");
  std::vector<Basis> Bases;
  Bases.push_back(F1);
  Bases.push_back(F2);
  Hamiltonian<DT> ham( Bases );
  std::vector< std::vector<DT> > Vloc;
  std::vector<DT> Vtmp;//(L, 1.0);
  for ( RealType &val : Vin ){
    Vtmp.push_back((DT)val);
  }
  Vloc.push_back(Vtmp);
  Vloc.push_back(Vtmp);
  std::vector< std::vector<DT> > Uloc;
  // std::vector<DT> Utmp(L, DT(10.0e0, 0.0e0) );
  std::vector<DT> Utmp(L, (DT)Uin);
  Uloc.push_back(Utmp);
  Uloc.push_back(Utmp);
  ham.BuildLocalHamiltonian(Vloc, Uloc, Bases);
  INFO(" - BuildLocalHamiltonian DONE!");
  ham.BuildHoppingHamiltonian(Bases, lattice);
  INFO(" - BuildHoppingHamiltonian DONE!");
  ham.BuildTotalHamiltonian();
  INFO("DONE!");
  INFO_NONEWLINE("Diagonalize Hamiltonian - ");
  std::vector<RealType> Val;
  Hamiltonian<DT>::VectorType Vec;
  ham.eigh(Val, Vec);
  INFO("GS energy = " << Val.at(0));
  file->saveVector("GS", "EVec", Vec);
  file->saveStdVector("GS", "EVal", Val);
  INFO("DONE!");
  std::vector< DTV > Nfi = Ni( Bases, Vec, ham );
  INFO(" Up Spin - ");
  INFO(Nfi.at(0));
  INFO(" Down Spin - ");
  INFO(Nfi.at(1));
  INFO(" N_i - ");
  DTV Niall = Nfi.at(0) + Nfi.at(1);
  INFO(Niall);
  DTM Nud = NupNdn( Bases, Vec, ham );
  INFO(" Correlation NupNdn");
  INFO(Nud);
  DTM Nuu = NupNup( Bases, Vec, ham );
  INFO(" Correlation NupNup");
  INFO(Nuu);
  DTM Ndd = NdnNdn( Bases, Vec, ham );
  INFO(" Correlation NdnNdn");
  INFO(Ndd);
  INFO(" N_i^2 - ");
  DTM Ni2 = Nuu.diagonal() + Ndd.diagonal() + 2.0e0 * Nud.diagonal();
  INFO(Ni2);
  file->saveVector("Obs", "Nup", Nfi.at(0));
  file->saveVector("Obs", "Ndn", Nfi.at(1));
  file->saveMatrix("Obs", "NupNdn", Nud);
  file->saveMatrix("Obs", "NupNup", Nuu);
  file->saveMatrix("Obs", "NdnNdn", Ndd);
  delete file;
  if ( dynamics ){
    ComplexType Prefactor = ComplexType(0.0, -1.0e0*dt);/* NOTE: hbar = 1 */
    std::cout << "Begin dynamics......" << std::endl;
    std::cout << "Cut the boundary." << std::endl;
    J.pop_back();
    std::vector< Node<DT>* > lattice2 = NN_1D_Chain(L, J, true);// cut to open
    ham.BuildHoppingHamiltonian(Bases, lattice2);
    INFO(" - Update Hopping Hamiltonian DONE!");
    ham.BuildTotalHamiltonian();
    INFO(" - Update Total Hamiltonian DONE!");
    for (size_t cntT = 1; cntT <= Tsteps; cntT++) {
      ham.expH(Prefactor, Vec);
      if ( cntT % 2 == 0 ){
        HDF5IO file2("DYN.h5");
        std::string gname = "Obs-";
        gname.append(std::to_string((unsigned long long)cntT));
        gname.append("/");
        Nfi = Ni( Bases, Vec, ham );
        Nud = NupNdn( Bases, Vec, ham );
        Nuu = NupNup( Bases, Vec, ham );
        Ndd = NdnNdn( Bases, Vec, ham );
        file2.saveVector(gname, "Nup", Nfi.at(0));
        file2.saveVector(gname, "Ndn", Nfi.at(1));
        file2.saveMatrix(gname, "NupNdn", Nud);
        file2.saveMatrix(gname, "NupNup", Nuu);
        file2.saveMatrix(gname, "NdnNdn", Ndd);
      }
    }
  }
  return 0;
}