diff --git a/pyblock2/driver/core.py b/pyblock2/driver/core.py index d5c6f037..1d575757 100644 --- a/pyblock2/driver/core.py +++ b/pyblock2/driver/core.py @@ -478,10 +478,12 @@ def su2_to_sgf(self): [self.orb_sym[i // 2] for i in range(self.n_sites)] ) - def get_conventional_qc_mpo(self, fcidump, algo_type=None): + def get_conventional_qc_mpo(self, fcidump, algo_type=None, iprint=1): """This method cannot take care of parallelization!""" bw = self.bw hamil = bw.bs.HamiltonianQC(self.vacuum, self.n_sites, self.orb_sym, fcidump) + import time + tt = time.perf_counter() if algo_type is not None and MPOAlgorithmTypes.NC in algo_type: mpo = bw.bs.MPOQC(hamil, bw.b.QCTypes.NC, "HQC") elif algo_type is not None and MPOAlgorithmTypes.CN in algo_type: @@ -497,6 +499,12 @@ def get_conventional_qc_mpo(self, fcidump, algo_type=None): use_r_intermediates = ( algo_type is None or MPOAlgorithmTypes.NoRIntermed not in algo_type ) + + if iprint: + nnz, sz, bdim = mpo.get_summary() + print('Ttotal = %10.3f MPO method = %s bond dimension = %7d NNZ = %12d SIZE = %12d SPT = %6.4f' + % (time.perf_counter() - tt, algo_type.name, bdim, nnz, sz, (1.0 * sz - nnz) / sz)) + mpo = bw.bs.SimplifiedMPO( mpo, rule, @@ -694,7 +702,7 @@ def get_qc_mpo( if algo_type is not None and MPOAlgorithmTypes.Conventional in algo_type: fd = self.write_fcidump(h1e, g2e, ecore=ecore) - return self.get_conventional_qc_mpo(fd, algo_type=algo_type) + return self.get_conventional_qc_mpo(fd, algo_type=algo_type, iprint=iprint) # build Hamiltonian expression b = self.expr_builder() diff --git a/src/dmrg/mpo.hpp b/src/dmrg/mpo.hpp index 1006d937..30c9cd4f 100644 --- a/src/dmrg/mpo.hpp +++ b/src/dmrg/mpo.hpp @@ -282,6 +282,18 @@ template struct MPO { if (tensors[m] != nullptr) tensors[m]->deallocate(); } + // nnz, size, bond dimension + tuple get_summary() const { + size_t lnnz = 0, lsz = 0, rnnz = 0, rsz = 0; + int bdim = 0; + for (int ii = 0; ii < n_sites; ii++) { + shared_ptr> opt = tensors[ii]; + lnnz += opt->lmat->nnz(), rnnz += opt->rmat->nnz(); + lsz += opt->lmat->size(), rsz += opt->rmat->size(); + bdim = max(max(bdim, (int)opt->lmat->n), (int)opt->lmat->m); + } + return make_tuple(max(lnnz, rnnz), max(lsz, rsz), bdim); + } string get_filename(int i, int ixtag, const string &dir = "") const { const static string xtag[] = {"TENSOR", "LEFT.OP", "RIGHT.OP", "MIDDLE.OP"}; diff --git a/src/pybind/pybind_dmrg.hpp b/src/pybind/pybind_dmrg.hpp index 01497275..0c02acfd 100644 --- a/src/pybind/pybind_dmrg.hpp +++ b/src/pybind/pybind_dmrg.hpp @@ -1576,6 +1576,7 @@ template void bind_fl_mpo(py::module &m) { .def_readwrite("tag", &MPO::tag) .def_readwrite("tread", &MPO::tread) .def_readwrite("twrite", &MPO::twrite) + .def("get_summary", &MPO::get_summary) .def("get_filename", &MPO::get_filename) .def("load_left_operators", &MPO::load_left_operators) .def("save_left_operators", &MPO::save_left_operators)