votca 2026-dev
Loading...
Searching...
No Matches
orbitals.h
Go to the documentation of this file.
1
2/*
3 * Copyright 2009-2023 The VOTCA Development Team
4 * (http://www.votca.org)
5 *
6 * Licensed under the Apache License, Version 2.0 (the "License")
7 *
8 * You may not use this file except in compliance with the License.
9 * You may obtain a copy of the License at
10 *
11 * http://www.apache.org/licenses/LICENSE-2.0
12 *
13 * Unless required by applicable law or agreed to in writing, software
14 * distributed under the License is distributed on an "AS IS" BASIS,
15 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
16 * See the License for the specific language governing permissions and
17 * limitations under the License.
18 *
19 */
20
21#pragma once
22#ifndef VOTCA_XTP_ORBITALS_H
23#define VOTCA_XTP_ORBITALS_H
24
25// VOTCA includes
26#include <votca/tools/globals.h>
28
29// Local VOTCA includes
30#include "aobasis.h"
31#include "checkpoint.h"
32#include "classicalsegment.h"
33#include "eigen.h"
34#include "qmmolecule.h"
35#include "qmstate.h"
36
37namespace votca {
38namespace xtp {
39
47class Orbitals {
48 public:
50 Orbitals();
51
53 bool hasBasisSetSize() const {
54 return (dftbasis_.AOBasisSize() > 0) ? true : false;
55 }
56
59
62
64 void setTruncMOsFullBasis(const Eigen::MatrixXd &expandedMOs) {
65 expandedMOs_ = expandedMOs;
66 }
67
69 const Eigen::MatrixXd getTruncMOsFullBasis() const { return expandedMOs_; }
70
72 Index getBasisSetSize() const { return dftbasis_.AOBasisSize(); }
73
78 }
79 return occupied_levels_;
80 }
81
83 Index getHomoAlpha() const { return getLumoAlpha() - 1; }
84
88 if (number_beta_electrons_ > 0) {
90 }
91
92 // Legacy restricted closed-shell fallback
93 if (!hasUnrestrictedOrbitals() && total_spin_ == 1) {
94 return occupied_levels_;
95 }
96
98 }
99
101 Index getHomoBeta() const { return getLumoBeta() - 1; }
102
103 // Current generic convention: use alpha frontier orbitals
105 Index getLumo() const { return getLumoAlpha(); }
107 Index getHomo() const { return getHomoAlpha(); }
108
109 // access to DFT number of levels, new, tested
111 bool hasNumberOfLevels() const {
112 return ((occupied_levels_ > 0) ? true : false);
113 }
114
117 return ((occupied_levels_beta_ > 0) ? true : false);
118 }
119
122 void setNumberOfOccupiedLevels(Index occupied_levels) {
123 occupied_levels_ = occupied_levels;
124
125 // Backward compatibility for legacy restricted/singlet workflows:
126 // many callers only set occupied_levels_ and expect a closed-shell density.
127 if (!hasUnrestrictedOrbitals() && total_spin_ == 1) {
128 number_alpha_electrons_ = occupied_levels;
129 number_beta_electrons_ = occupied_levels;
130 }
131 }
132
134 void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta) {
135 occupied_levels_beta_ = occupied_levels_beta;
136 }
137
140 number_alpha_electrons_ = electrons;
141 }
142
145 number_beta_electrons_ = electrons;
146 }
147
148 // access to DFT number of electrons, new, tested
151 return (number_alpha_electrons_ > 0) ? true : false;
152 }
153
155 return (number_beta_electrons_ > 0) ? true : false;
156 }
157
162
164 bool hasECPName() const { return (ECP_ != "") ? true : false; }
165
167 const std::string &getECPName() const { return ECP_; };
168
170 void setECPName(const std::string &ECP) { ECP_ = ECP; };
171
172 // access to QM package name, new, tested
173
175 bool hasQMpackage() const { return (!qm_package_.empty()); }
176
178 const std::string &getQMpackage() const { return qm_package_; }
179
181 void setQMpackage(const std::string &qmpackage) { qm_package_ = qmpackage; }
182
183 // access to DFT molecular orbital energies, new, tested
185 bool hasMOs() const { return (mos_.eigenvalues().size() > 0) ? true : false; }
187 bool hasBetaMOs() const {
188 return (mos_beta_.eigenvalues().size() > 0) ? true : false;
189 }
190
192 const tools::EigenSystem &MOs() const { return mos_; }
195
197 const Eigen::MatrixXd &Occupations() const { return occupations_; }
199 Eigen::MatrixXd &Occupations() { return occupations_; }
200
202 const tools::EigenSystem &MOs_beta() const { return mos_beta_; }
205
206 // determine (pseudo-)degeneracy of a DFT molecular orbital
209 std::vector<Index> CheckDegeneracy(Index level,
210 double energy_difference) const;
211
214 switch (type.Type()) {
216 return Index(BSE_singlet_.eigenvalues().size());
217 break;
219 return Index(BSE_triplet_.eigenvalues().size());
220 break;
222 return Index(mos_.eigenvalues().size());
223 break;
225 return Index(QPpert_energies_.size());
226 break;
228 return Index(QPdiag_.eigenvalues().size());
229 break;
231 return Index(BSE_uks_.eigenvalues().size());
232 break;
233 default:
234 return 1;
235 break;
236 }
237 }
238
240 void setCalculationType(std::string CalcType) { CalcType_ = CalcType; }
242 std::string getCalculationType() const { return CalcType_; }
243
246 void setChargeAndSpin(Index charge, Index spin) {
247 total_charge_ = charge;
248 total_spin_ = spin;
249 }
250
252 Index getSpin() const { return total_spin_; }
254 Index getCharge() const { return total_charge_; }
256 bool isOpenShell() const { return (total_spin_ > 1) ? true : false; }
257
259 bool hasQMAtoms() const { return (atoms_.size() > 0) ? true : false; }
260
262 const QMMolecule &QMAtoms() const { return atoms_; }
263
265 QMMolecule &QMAtoms() { return atoms_; }
266
268 void updateAtomPostion(Index atom_index, Eigen::Vector3d new_position) {
269 this->QMAtoms()[atom_index].setPos(new_position);
270 dftbasis_.UpdateShellPositions(this->QMAtoms());
271 auxbasis_.UpdateShellPositions(this->QMAtoms());
272 }
273
276 void setXCFunctionalName(std::string functionalname) {
277 functionalname_ = functionalname;
278 }
279
280 const std::string &getXCFunctionalName() const { return functionalname_; }
281
283 void setXCGrid(std::string grid) { grid_quality_ = grid; }
285 const std::string &getXCGrid() const { return grid_quality_; }
286
287 // access to QM total energy, new, tested
289 bool hasQMEnergy() const { return (qm_energy_ != 0.0) ? true : false; }
290
292 double getDFTTotalEnergy() const { return qm_energy_; }
293
295 void setQMEnergy(double qmenergy) { qm_energy_ = qmenergy; }
296
297 // access to DFT basis set name
299 bool hasDFTbasisName() const {
300 return (!dftbasis_.Name().empty()) ? true : false;
301 }
302
304 const std::string &getDFTbasisName() const { return dftbasis_.Name(); }
305
307 void SetupDftBasis(std::string basis_name);
308
311 void SetupAuxBasis(std::string aux_basis_name);
312
314 const AOBasis &getDftBasis() const {
315 if (dftbasis_.AOBasisSize() == 0) {
316 throw std::runtime_error(
317 "Requested the DFT basis, but no basis is present. Make sure "
318 "SetupDftBasis is called.");
319 } else {
320 return dftbasis_;
321 }
322 }
323
325 const AOBasis &getAuxBasis() const {
326 if (auxbasis_.AOBasisSize() == 0) {
327 throw std::runtime_error(
328 "Requested the Aux basis, but no basis is present. Make sure "
329 "SetupAuxBasis is called.");
330 } else {
331 return auxbasis_;
332 }
333 }
334
335 /*
336 * ======= GW-BSE related functions =======
337 */
338
339 // access to auxiliary basis set name
340
342 bool hasAuxbasisName() const {
343 return (!auxbasis_.Name().empty()) ? true : false;
344 }
345
347 const std::string getAuxbasisName() const { return auxbasis_.Name(); }
348
349 // access to list of indices used in GWA
350
352 bool hasGWAindices() const { return (qpmax_ > 0) ? true : false; }
353
355 void setGWindices(Index qpmin, Index qpmax) {
356 qpmin_ = qpmin;
357 qpmax_ = qpmax;
358 }
359
361 Index getGWAmin() const { return qpmin_; }
362
364 Index getGWAmax() const { return qpmax_; }
365
366 // access to list of indices used in RPA
367
369 bool hasRPAindices() const { return (rpamax_ > 0) ? true : false; }
370
372 void setRPAindices(Index rpamin, Index rpamax) {
373 rpamin_ = rpamin;
374 rpamax_ = rpamax;
375 }
376
378 Index getRPAmin() const { return rpamin_; }
379
381 Index getRPAmax() const { return rpamax_; }
382
383 // access to list of indices used in BSE
384
386 void setTDAApprox(bool usedTDA) { useTDA_ = usedTDA; }
388 bool getTDAApprox() const { return useTDA_; }
389
391 bool hasBSEindices() const { return (bse_cmax_ > 0) ? true : false; }
392
395 void setBSEindices(Index vmin, Index cmax) {
396 bse_vmin_ = vmin;
397 bse_vmax_ = getHomo();
398 bse_cmin_ = getLumo();
399 bse_cmax_ = cmax;
403 return;
404 }
405
407 Index getBSEvmin() const { return bse_vmin_; }
408
410 Index getBSEvmax() const { return bse_vmax_; }
411
413 Index getBSEcmin() const { return bse_cmin_; }
414
416 Index getBSEcmax() const { return bse_cmax_; }
417
419 double getScaHFX() const { return ScaHFX_; }
420
422 void setScaHFX(double ScaHFX) { ScaHFX_ = ScaHFX; }
423
424 // access to perturbative QP energies
426 bool hasRPAInputEnergies() const { return (rpa_inputenergies_.size() > 0); }
427
429 const Eigen::VectorXd &RPAInputEnergies() const { return rpa_inputenergies_; }
430
432 Eigen::VectorXd &RPAInputEnergies() { return rpa_inputenergies_; }
433
434 // access to RPA input energies energies
436 bool hasQPpert() const {
437 return (QPpert_energies_.size() > 0) ? true : false;
438 }
439
441 const Eigen::VectorXd &QPpertEnergies() const { return QPpert_energies_; }
442
444 Eigen::VectorXd &QPpertEnergies() { return QPpert_energies_; }
445
446 // access to diagonalized QP energies and wavefunctions
447
450 bool hasQPdiag() const {
451 return (QPdiag_.eigenvalues().size() > 0) ? true : false;
452 }
453
454 const tools::EigenSystem &QPdiag() const { return QPdiag_; }
457
459 bool hasBSETriplets() const {
460 return (BSE_triplet_.eigenvectors().cols() > 0) ? true : false;
461 }
462
464 const tools::EigenSystem &BSETriplets() const { return BSE_triplet_; }
465
468
469 // access to singlet energies and wave function coefficients
470
472 bool hasBSESinglets() const {
473 return (BSE_singlet_.eigenvectors().cols() > 0) ? true : false;
474 }
475
477 const tools::EigenSystem &BSESinglets() const { return BSE_singlet_; }
478
481
482 // access to BSE energies with dynamical screening
485 return (BSE_singlet_energies_dynamic_.size() > 0) ? true : false;
486 }
487
489 const Eigen::VectorXd &BSESinglets_dynamic() const {
491 }
492
494 Eigen::VectorXd &BSESinglets_dynamic() {
496 }
497
500 return (BSE_triplet_energies_dynamic_.size() > 0) ? true : false;
501 }
502
504 const Eigen::VectorXd &BSETriplets_dynamic() const {
506 }
507
509 Eigen::VectorXd &BSETriplets_dynamic() {
511 }
512
513 // access to transition dipole moments
514
516 bool hasTransitionDipoles() const {
517 return (transition_dipoles_.size() > 0) ? true : false;
518 }
519
521 const std::vector<Eigen::Vector3d> &TransitionDipoles() const {
522 return transition_dipoles_;
523 }
524
527 Eigen::VectorXd Oscillatorstrengths() const;
528
530 Eigen::VectorXd Oscillatorstrengths(const QMStateType &type) const;
531
533 Eigen::Vector3d CalcElDipole(const QMState &state) const;
534
535 // Calculates full electron density for state or transition density, if you
536 // want to calculate only the density contribution of hole or electron use
537 // DensityMatrixExcitedState
540 Eigen::MatrixXd DensityMatrixFull(const QMState &state) const;
543 Eigen::MatrixXd DensityMatrixWithoutGS(const QMState &state) const;
544
545 // functions for calculating density matrices
547 Eigen::MatrixXd DensityMatrixGroundState() const;
549 std::array<Eigen::MatrixXd, 2> DensityMatrixExcitedState(
550 const QMState &state) const;
552 Eigen::MatrixXd DensityMatrixQuasiParticle(const QMState &state) const;
554 Eigen::MatrixXd DensityMatrixKSstate(const QMState &state) const;
556 Eigen::MatrixXd CalculateQParticleAORepresentation() const;
558 double getTotalStateEnergy(const QMState &state) const; // Hartree
560 double getExcitedStateEnergy(const QMState &state) const; // Hartree
561
564 void OrderMOsbyEnergy();
565
568 void PrepareDimerGuess(const Orbitals &orbitalsA, const Orbitals &orbitalsB);
569
573
576
578 void WriteToCpt(const std::string &filename) const;
579
581 void ReadFromCpt(const std::string &filename);
582
584 void WriteToCpt(CheckpointWriter &w) const;
592
595 bool GetFlagUseHqpOffdiag() const { return use_Hqp_offdiag_; };
598 void SetFlagUseHqpOffdiag(bool flag) { use_Hqp_offdiag_ = flag; };
599
601 const Eigen::MatrixXd &getLMOs() const { return lmos_; };
603 void setLMOs(const Eigen::MatrixXd &matrix) { lmos_ = matrix; }
604
606 const Eigen::VectorXd &getLMOs_energies() const { return lmos_energies_; };
608 void setLMOs_energies(const Eigen::VectorXd &energies) {
609 lmos_energies_ = energies;
610 }
611
615 void setNumofActiveElectrons(const Index active_electrons) {
616 active_electrons_ = active_electrons;
617 }
618
620 const Eigen::MatrixXd &getInactiveDensity() const { return inactivedensity_; }
622 void setInactiveDensity(Eigen::MatrixXd inactivedensity) {
623 inactivedensity_ = inactivedensity;
624 }
625
626 /************************************
627 * Extensions spin DFT
628 *************************************/
630 bool hasUnrestrictedOrbitals() const { return hasBetaMOs(); }
631
634 return total_spin_ > 1 && !hasUnrestrictedOrbitals();
635 }
636
640 return total_spin_ > 1 && hasUnrestrictedOrbitals();
641 }
642
644 std::array<Eigen::MatrixXd, 2> DensityMatrixGroundStateSpinResolved() const;
645
646 /************************************
647 * Extensions spin GW
648 *************************************/
650 return (rpa_inputenergies_alpha_.size() > 0);
651 }
653 return (rpa_inputenergies_beta_.size() > 0);
654 }
655 const Eigen::VectorXd &RPAInputEnergiesAlpha() const {
657 }
658 const Eigen::VectorXd &RPAInputEnergiesBeta() const {
660 }
661 Eigen::VectorXd &RPAInputEnergiesAlpha() { return rpa_inputenergies_alpha_; }
662 Eigen::VectorXd &RPAInputEnergiesBeta() { return rpa_inputenergies_beta_; }
663
664 bool hasQPpertAlpha() const { return (QPpert_energies_alpha_.size() > 0); }
665 bool hasQPpertBeta() const { return (QPpert_energies_beta_.size() > 0); }
666 const Eigen::MatrixXd &QPpertEnergiesAlpha() const {
668 }
669 const Eigen::MatrixXd &QPpertEnergiesBeta() const {
671 }
672 Eigen::MatrixXd &QPpertEnergiesAlpha() { return QPpert_energies_alpha_; }
673 Eigen::MatrixXd &QPpertEnergiesBeta() { return QPpert_energies_beta_; }
674
675 bool hasQPdiagAlpha() const {
676 return (QPdiag_alpha_.eigenvalues().size() > 0);
677 }
678 bool hasQPdiagBeta() const { return (QPdiag_beta_.eigenvalues().size() > 0); }
680 const tools::EigenSystem &QPdiagBeta() const { return QPdiag_beta_; }
683
687 bool hasBSEUKS() const {
688 return (BSE_uks_.eigenvectors().cols() > 0) ? true : false;
689 }
690 const tools::EigenSystem &BSEUKS() const { return BSE_uks_; }
692
693 bool hasBSEUKS_dynamic() const {
694 return (BSE_uks_energies_dynamic_.size() > 0) ? true : false;
695 }
696 const Eigen::VectorXd &BSEUKS_dynamic() const {
698 }
699 Eigen::VectorXd &BSEUKS_dynamic() { return BSE_uks_energies_dynamic_; }
700
701 private:
702 std::array<Eigen::MatrixXd, 3> CalcFreeTransition_Dipoles() const;
703
704 // returns indeces of a re-sorted vector of energies from lowest to highest
705 std::vector<Index> SortEnergies();
706
707 void WriteToCpt(CheckpointFile &f) const;
708
710 Eigen::MatrixXd TransitionDensityMatrix(const QMState &state) const;
711 std::array<Eigen::MatrixXd, 2> DensityMatrixExcitedState_R(
712 const QMState &state) const;
713 std::array<Eigen::MatrixXd, 2> DensityMatrixExcitedState_AR(
714 const QMState &state) const;
715 Eigen::MatrixXd CalcAuxMat_cc(const Eigen::VectorXd &coeffs) const;
716 Eigen::MatrixXd CalcAuxMat_vv(const Eigen::VectorXd &coeffs) const;
717
722 std::string ECP_ = "";
724
725 std::string CalcType_ = "NoEmbedding";
726
729 Eigen::MatrixXd occupations_;
730
732
733 Eigen::MatrixXd lmos_;
734 Eigen::VectorXd lmos_energies_;
736 Eigen::MatrixXd inactivedensity_;
737 Eigen::MatrixXd expandedMOs_;
738
740
743
744 double qm_energy_ = 0;
745
748
749 // new variables for GW-BSE storage
752
762
763 double ScaHFX_ = 0;
764
765 std::string functionalname_ = "";
766 std::string grid_quality_ = "";
767 std::string qm_package_ = "";
768
769 Eigen::VectorXd rpa_inputenergies_;
770 // perturbative quasiparticle energies
771 Eigen::VectorXd QPpert_energies_;
772
773 // quasiparticle energies and coefficients after diagonalization
775
777 std::vector<Eigen::Vector3d> transition_dipoles_;
779
780 // singlet and triplet energies after perturbative dynamical screening
783
784 bool use_Hqp_offdiag_ = false;
785
786 // Spin-GW additions
788 Eigen::VectorXd rpa_inputenergies_beta_;
789
790 Eigen::MatrixXd QPpert_energies_alpha_;
791 Eigen::MatrixXd QPpert_energies_beta_;
792
795
798
799 // Version 2: adds BSE energies after perturbative dynamical screening
800 // Version 3: changed shell ordering
801 // Version 4: added vxc grid quality
802 // Version 5: added the dft and aux basisset
803 // Version 6: added spin in dft
804 static constexpr int orbitals_version() { return 8; }
805};
806
807} // namespace xtp
808} // namespace votca
809
810#endif // VOTCA_XTP_ORBITALS_H
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
tools::EigenSystem QPdiag_alpha_
Definition orbitals.h:793
Index getBSEvmax() const
Return the highest valence orbital included in BSE.
Definition orbitals.h:410
const Eigen::MatrixXd getTruncMOsFullBasis() const
Return truncated active-region orbitals represented in the full AO basis.
Definition orbitals.h:69
void setScaHFX(double ScaHFX)
Store the fraction of exact exchange associated with the functional.
Definition orbitals.h:422
const Eigen::VectorXd & BSEUKS_dynamic() const
Definition orbitals.h:696
void setCalculationType(std::string CalcType)
Store the calculation-type tag used by downstream workflows.
Definition orbitals.h:240
Eigen::VectorXd & RPAInputEnergiesBeta()
Definition orbitals.h:662
const tools::EigenSystem & BSETriplets() const
Return read-only access to triplet BSE eigenpairs.
Definition orbitals.h:464
bool GetFlagUseHqpOffdiag() const
Definition orbitals.h:595
Index getNumOfActiveElectrons()
Return the number of electrons assigned to the active region.
Definition orbitals.h:613
Index number_beta_electrons_
Definition orbitals.h:721
void CalcCoupledTransition_Dipoles()
Definition orbitals.cc:760
bool hasQPdiag() const
Definition orbitals.h:450
Index getCharge() const
Return the stored total charge.
Definition orbitals.h:254
Eigen::MatrixXd occupations_
Definition orbitals.h:729
bool hasQPpertAlpha() const
Definition orbitals.h:664
Index getBSEcmin() const
Return the lowest conduction orbital included in BSE.
Definition orbitals.h:413
bool hasRPAindices() const
Report whether the RPA window has been defined.
Definition orbitals.h:369
Eigen::VectorXd rpa_inputenergies_
Definition orbitals.h:769
const tools::EigenSystem & QPdiag() const
Return read-only access to the diagonalized quasiparticle representation.
Definition orbitals.h:454
Index getHomo() const
Return the default HOMO index used by spin-agnostic callers.
Definition orbitals.h:107
Eigen::MatrixXd & QPpertEnergiesBeta()
Definition orbitals.h:673
Eigen::VectorXd BSE_uks_energies_dynamic_
Definition orbitals.h:797
Eigen::MatrixXd TransitionDensityMatrix(const QMState &state) const
Definition orbitals.cc:316
double getDFTTotalEnergy() const
Return the stored total DFT energy.
Definition orbitals.h:292
bool hasBSESinglets_dynamic() const
Report whether dynamically screened singlet BSE energies are available.
Definition orbitals.h:484
void setInactiveDensity(Eigen::MatrixXd inactivedensity)
Store the inactive-region density matrix used in embedding workflows.
Definition orbitals.h:622
Eigen::VectorXd & BSEUKS_dynamic()
Definition orbitals.h:699
bool hasBSEindices() const
Report whether the BSE valence/conduction window has been defined.
Definition orbitals.h:391
bool hasAuxbasisName() const
Report whether an auxiliary basis-set name has been stored.
Definition orbitals.h:342
Index getSpin() const
Return the stored spin multiplicity.
Definition orbitals.h:252
Eigen::MatrixXd CalculateQParticleAORepresentation() const
Transform quasiparticle amplitudes into the AO representation.
Definition orbitals.cc:250
Eigen::MatrixXd & QPpertEnergiesAlpha()
Definition orbitals.h:672
std::array< Eigen::MatrixXd, 2 > DensityMatrixExcitedState_R(const QMState &state) const
Definition orbitals.cc:535
const tools::EigenSystem & MOs_beta() const
Return read-only access to beta-spin molecular orbitals.
Definition orbitals.h:202
void setNumofActiveElectrons(const Index active_electrons)
Store the number of electrons assigned to the active region.
Definition orbitals.h:615
bool hasUnrestrictedOrbitals() const
Report whether separate beta-spin orbitals are present.
Definition orbitals.h:630
Eigen::VectorXd lmos_energies_
Definition orbitals.h:734
const tools::EigenSystem & BSESinglets() const
Return read-only access to singlet BSE eigenpairs.
Definition orbitals.h:477
Index getBSEcmax() const
Return the highest conduction orbital included in BSE.
Definition orbitals.h:416
Eigen::VectorXd & BSESinglets_dynamic()
Return writable access to dynamically screened singlet BSE energies.
Definition orbitals.h:494
const Eigen::VectorXd & QPpertEnergies() const
Return read-only access to perturbative quasiparticle energies.
Definition orbitals.h:441
bool getTDAApprox() const
Return whether the Tamm-Dancoff approximation is enabled.
Definition orbitals.h:388
tools::EigenSystem BSE_singlet_
Definition orbitals.h:776
tools::EigenSystem QPdiag_beta_
Definition orbitals.h:794
void PrepareDimerGuess(const Orbitals &orbitalsA, const Orbitals &orbitalsB)
Guess for a dimer based on monomer orbitals.
Definition orbitals.cc:906
bool hasNumberOfLevelsBeta() const
Definition orbitals.h:116
tools::EigenSystem mos_beta_
Definition orbitals.h:728
tools::EigenSystem & QPdiagAlpha()
Definition orbitals.h:681
std::string grid_quality_
Definition orbitals.h:766
bool hasMOs() const
Report whether alpha/restricted molecular orbitals are available.
Definition orbitals.h:185
void setTruncMOsFullBasis(const Eigen::MatrixXd &expandedMOs)
Store truncated active-region orbitals expanded back to the full AO basis.
Definition orbitals.h:64
bool hasBSESinglets() const
Report whether singlet BSE eigenpairs are available.
Definition orbitals.h:472
double getTotalStateEnergy(const QMState &state) const
Return the absolute total energy of the requested state in Hartree.
Definition orbitals.cc:676
Eigen::VectorXd BSE_triplet_energies_dynamic_
Definition orbitals.h:782
bool hasQPpertBeta() const
Definition orbitals.h:665
const Eigen::MatrixXd & Occupations() const
Return the stored fractional occupation matrix.
Definition orbitals.h:197
void ReadBasisSetsFromCpt(CheckpointReader &r)
Deserialize attached AO basis sets from an already-open checkpoint reader.
Definition orbitals.cc:1076
std::array< Eigen::MatrixXd, 2 > DensityMatrixGroundStateSpinResolved() const
Definition orbitals.cc:1231
Eigen::MatrixXd QPpert_energies_beta_
Definition orbitals.h:791
void WriteBasisSetsToCpt(CheckpointWriter &w) const
Definition orbitals.cc:983
tools::EigenSystem & QPdiagBeta()
Definition orbitals.h:682
bool hasRPAInputEnergies() const
Report whether RPA input energies are available.
Definition orbitals.h:426
bool isRestrictedOpenShell() const
Report whether the orbitals represent a restricted open-shell reference.
Definition orbitals.h:633
std::vector< Index > CheckDegeneracy(Index level, double energy_difference) const
Definition orbitals.cc:56
tools::EigenSystem & BSETriplets()
Return writable access to triplet BSE eigenpairs.
Definition orbitals.h:467
const AOBasis & getAuxBasis() const
Return the auxiliary AO basis, throwing if it has not been initialized.
Definition orbitals.h:325
Eigen::MatrixXd inactivedensity_
Definition orbitals.h:736
Index getNumberOfAlphaElectrons() const
Return the stored number of alpha electrons.
Definition orbitals.h:159
void setLMOs_energies(const Eigen::VectorXd &energies)
Store the energies associated with localized molecular orbitals.
Definition orbitals.h:608
std::string getCalculationType() const
Return the stored calculation-type tag.
Definition orbitals.h:242
bool hasQMAtoms() const
Report whether a molecular geometry has been stored.
Definition orbitals.h:259
std::array< Eigen::MatrixXd, 3 > CalcFreeTransition_Dipoles() const
Definition orbitals.cc:742
void setLMOs(const Eigen::MatrixXd &matrix)
Store localized molecular orbitals.
Definition orbitals.h:603
bool hasNumberOfBetaElectrons() const
Report whether the beta-electron count has been set explicitly.
Definition orbitals.h:154
const std::string & getQMpackage() const
Return the stored QM package name.
Definition orbitals.h:178
QMMolecule & QMAtoms()
Return writable access to the molecular geometry.
Definition orbitals.h:265
const Eigen::MatrixXd & QPpertEnergiesBeta() const
Definition orbitals.h:669
Eigen::MatrixXd CalcAuxMat_cc(const Eigen::VectorXd &coeffs) const
Definition orbitals.cc:584
void setRPAindices(Index rpamin, Index rpamax)
Store the orbital window used in RPA calculations.
Definition orbitals.h:372
bool hasNumberOfLevels() const
Report whether the number of occupied spatial orbitals has been set.
Definition orbitals.h:111
Eigen::MatrixXd expandedMOs_
Definition orbitals.h:737
void setNumberOfAlphaElectrons(Index electrons)
Store the total number of alpha electrons.
Definition orbitals.h:139
Index getBSEvmin() const
Return the lowest valence orbital included in BSE.
Definition orbitals.h:407
void updateAtomPostion(Index atom_index, Eigen::Vector3d new_position)
Update one atomic position and keep attached AO basis shells synchronized.
Definition orbitals.h:268
const tools::EigenSystem & BSEUKS() const
Definition orbitals.h:690
Eigen::VectorXd & RPAInputEnergies()
Return writable access to the RPA input energies.
Definition orbitals.h:432
Eigen::VectorXd Oscillatorstrengths() const
Definition orbitals.cc:643
std::vector< Eigen::Vector3d > transition_dipoles_
Definition orbitals.h:777
const tools::EigenSystem & QPdiagAlpha() const
Definition orbitals.h:679
static constexpr int orbitals_version()
Definition orbitals.h:804
std::string functionalname_
Definition orbitals.h:765
Eigen::MatrixXd lmos_
Definition orbitals.h:733
Eigen::MatrixXd & Occupations()
Return writable access to the fractional occupation matrix.
Definition orbitals.h:199
Eigen::MatrixXd DensityMatrixWithoutGS(const QMState &state) const
Definition orbitals.cc:121
bool hasGWAindices() const
Report whether the GW quasiparticle window has been defined.
Definition orbitals.h:352
tools::EigenSystem & MOs()
Return writable access to alpha/restricted molecular orbitals.
Definition orbitals.h:194
bool hasQMEnergy() const
Report whether a total electronic energy has been stored.
Definition orbitals.h:289
Eigen::MatrixXd DensityMatrixFull(const QMState &state) const
Definition orbitals.cc:150
bool hasTransitionDipoles() const
Report whether transition dipole moments have been computed.
Definition orbitals.h:516
const Eigen::MatrixXd & getLMOs() const
Return localized molecular orbitals, if available.
Definition orbitals.h:601
Eigen::MatrixXd CalcAuxMat_vv(const Eigen::VectorXd &coeffs) const
Definition orbitals.cc:578
tools::EigenSystem & MOs_beta()
Return writable access to beta-spin molecular orbitals.
Definition orbitals.h:204
tools::EigenSystem & QPdiag()
Return writable access to the diagonalized quasiparticle representation.
Definition orbitals.h:456
void SetupAuxBasis(std::string aux_basis_name)
Definition orbitals.cc:108
std::string qm_package_
Definition orbitals.h:767
Index getLumoAlpha() const
Return the alpha-spin LUMO index inferred from the stored electron counts.
Definition orbitals.h:75
void setNumberOfBetaElectrons(Index electrons)
Store the total number of beta electrons.
Definition orbitals.h:144
tools::EigenSystem BSE_triplet_
Definition orbitals.h:778
Index getRPAmax() const
Return the upper RPA orbital index.
Definition orbitals.h:381
bool hasBSEUKS() const
Definition orbitals.h:687
void setBSEindices(Index vmin, Index cmax)
Definition orbitals.h:395
bool hasBasisSetSize() const
Report whether a DFT AO basis has already been attached.
Definition orbitals.h:53
const Eigen::MatrixXd & getInactiveDensity() const
Return the inactive-region density matrix used in embedding workflows.
Definition orbitals.h:620
void setGWindices(Index qpmin, Index qpmax)
Store the orbital window used in GW calculations.
Definition orbitals.h:355
std::vector< Index > SortEnergies()
Definition orbitals.cc:81
Index getHomoBeta() const
Return the beta-spin HOMO index.
Definition orbitals.h:101
std::array< Eigen::MatrixXd, 2 > DensityMatrixExcitedState_AR(const QMState &state) const
Definition orbitals.cc:590
bool hasQMpackage() const
Report whether the originating QM package name has been stored.
Definition orbitals.h:175
double getExcitedStateEnergy(const QMState &state) const
Return the excitation energy of the requested state in Hartree.
Definition orbitals.cc:685
bool isUnrestrictedOpenShell() const
Definition orbitals.h:639
tools::EigenSystem & BSESinglets()
Return writable access to singlet BSE eigenpairs.
Definition orbitals.h:480
bool hasRPAInputEnergiesBeta() const
Definition orbitals.h:652
void setECPName(const std::string &ECP)
Store the effective core potential label.
Definition orbitals.h:170
void setXCGrid(std::string grid)
Store the numerical XC grid quality label.
Definition orbitals.h:283
void setNumberOfOccupiedLevels(Index occupied_levels)
Definition orbitals.h:122
const std::vector< Eigen::Vector3d > & TransitionDipoles() const
Return the stored transition dipole moments.
Definition orbitals.h:521
Index occupied_levels_beta_
Definition orbitals.h:719
Orbitals()
Construct an empty orbital container with default metadata.
Definition orbitals.cc:46
tools::EigenSystem mos_
Definition orbitals.h:727
tools::EigenSystem QPdiag_
Definition orbitals.h:774
Eigen::VectorXd rpa_inputenergies_beta_
Definition orbitals.h:788
const Eigen::MatrixXd & QPpertEnergiesAlpha() const
Definition orbitals.h:666
bool hasBSETriplets() const
Report whether triplet BSE eigenpairs are available.
Definition orbitals.h:459
Eigen::MatrixXd DensityMatrixKSstate(const QMState &state) const
Build the AO density matrix for a single KS orbital state.
Definition orbitals.cc:232
bool hasBSETriplets_dynamic() const
Report whether dynamically screened triplet BSE energies are available.
Definition orbitals.h:499
const std::string getAuxbasisName() const
Return the auxiliary basis-set name.
Definition orbitals.h:347
Eigen::MatrixXd DensityMatrixQuasiParticle(const QMState &state) const
Build the AO density matrix for a quasiparticle state.
Definition orbitals.cc:265
bool hasQPpert() const
Report whether perturbative quasiparticle energies are available.
Definition orbitals.h:436
double getScaHFX() const
Return the fraction of exact exchange associated with the functional.
Definition orbitals.h:419
Index getBasisSetSize() const
Return the number of AO basis functions in the DFT basis.
Definition orbitals.h:72
void setQMEnergy(double qmenergy)
Store the total DFT energy.
Definition orbitals.h:295
const Eigen::VectorXd & RPAInputEnergies() const
Return read-only access to the RPA input energies.
Definition orbitals.h:429
Eigen::VectorXd & RPAInputEnergiesAlpha()
Definition orbitals.h:661
Eigen::MatrixXd QPpert_energies_alpha_
Definition orbitals.h:790
const std::string & getXCGrid() const
Return the numerical XC grid quality label.
Definition orbitals.h:285
Eigen::VectorXd BSE_singlet_energies_dynamic_
Definition orbitals.h:781
Index getGWAmin() const
Return the lower GW orbital index.
Definition orbitals.h:361
bool hasRPAInputEnergiesAlpha() const
Definition orbitals.h:649
Eigen::VectorXd & QPpertEnergies()
Return writable access to perturbative quasiparticle energies.
Definition orbitals.h:444
bool hasECPName() const
Report whether an effective core potential label has been stored.
Definition orbitals.h:164
const Eigen::VectorXd & BSETriplets_dynamic() const
Return dynamically screened triplet BSE energies.
Definition orbitals.h:504
const tools::EigenSystem & MOs() const
Return read-only access to alpha/restricted molecular orbitals.
Definition orbitals.h:192
void setQMpackage(const std::string &qmpackage)
Store the name of the QM package that produced these orbitals.
Definition orbitals.h:181
Eigen::VectorXd & BSETriplets_dynamic()
Return writable access to dynamically screened triplet BSE energies.
Definition orbitals.h:509
Eigen::MatrixXd DensityMatrixGroundState() const
Build the ground-state AO density matrix from the stored orbitals.
Definition orbitals.cc:220
std::string CalcType_
Definition orbitals.h:725
Index getLumoBeta() const
Definition orbitals.h:87
const QMMolecule & QMAtoms() const
Return read-only access to the molecular geometry.
Definition orbitals.h:262
tools::EigenSystem mos_embedding_
Definition orbitals.h:731
const std::string & getXCFunctionalName() const
Return the exchange-correlation functional label.
Definition orbitals.h:280
Eigen::VectorXd rpa_inputenergies_alpha_
Definition orbitals.h:787
bool hasQPdiagBeta() const
Definition orbitals.h:678
Index getGWAmax() const
Return the upper GW orbital index.
Definition orbitals.h:364
void SetFlagUseHqpOffdiag(bool flag)
Definition orbitals.h:598
void ReadFromCpt(const std::string &filename)
Read the orbital container from a checkpoint file on disk.
Definition orbitals.cc:1063
void setEmbeddedMOs(tools::EigenSystem &system)
Store molecular orbitals obtained from an embedding calculation.
Definition orbitals.h:58
tools::EigenSystem BSE_uks_
Definition orbitals.h:796
void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta)
Store the number of occupied beta-spin orbitals.
Definition orbitals.h:134
const tools::EigenSystem & QPdiagBeta() const
Definition orbitals.h:680
void setChargeAndSpin(Index charge, Index spin)
Definition orbitals.h:246
const Eigen::VectorXd & RPAInputEnergiesBeta() const
Definition orbitals.h:658
Eigen::VectorXd QPpert_energies_
Definition orbitals.h:771
void WriteToCpt(const std::string &filename) const
Write the orbital container to a checkpoint file on disk.
Definition orbitals.cc:970
bool hasDFTbasisName() const
Report whether a DFT basis-set name has been stored.
Definition orbitals.h:299
const std::string & getECPName() const
Return the effective core potential label.
Definition orbitals.h:167
const std::string & getDFTbasisName() const
Return the DFT basis-set name.
Definition orbitals.h:304
bool hasBetaMOs() const
Report whether beta-spin molecular orbitals are available.
Definition orbitals.h:187
void SetupDftBasis(std::string basis_name)
Build and attach the DFT AO basis from the stored molecular geometry.
Definition orbitals.cc:99
void setTDAApprox(bool usedTDA)
Enable or disable the Tamm-Dancoff approximation flag.
Definition orbitals.h:386
const Eigen::VectorXd & RPAInputEnergiesAlpha() const
Definition orbitals.h:655
tools::EigenSystem & BSEUKS()
Definition orbitals.h:691
std::array< Eigen::MatrixXd, 2 > DensityMatrixExcitedState(const QMState &state) const
Build separate hole and electron AO densities for an excited state.
Definition orbitals.cc:413
bool isOpenShell() const
Report whether the stored state corresponds to an open-shell system.
Definition orbitals.h:256
const Eigen::VectorXd & BSESinglets_dynamic() const
Return dynamically screened singlet BSE energies.
Definition orbitals.h:489
std::string ECP_
Definition orbitals.h:722
Index number_alpha_electrons_
Definition orbitals.h:720
const tools::EigenSystem & getEmbeddedMOs() const
Return molecular orbitals obtained from an embedding calculation.
Definition orbitals.h:61
Index NumberofStates(QMStateType type) const
Return the number of states available for the requested state family.
Definition orbitals.h:213
Index getRPAmin() const
Return the lower RPA orbital index.
Definition orbitals.h:378
Eigen::Vector3d CalcElDipole(const QMState &state) const
Compute the electronic dipole moment associated with a state density.
Definition orbitals.cc:288
Index getNumberOfBetaElectrons() const
Return the stored number of beta electrons.
Definition orbitals.h:161
Index getHomoAlpha() const
Return the alpha-spin HOMO index.
Definition orbitals.h:83
const AOBasis & getDftBasis() const
Return the DFT AO basis, throwing if it has not been initialized.
Definition orbitals.h:314
const Eigen::VectorXd & getLMOs_energies() const
Return the energies associated with localized molecular orbitals.
Definition orbitals.h:606
QMMolecule atoms_
Definition orbitals.h:739
bool hasQPdiagAlpha() const
Definition orbitals.h:675
bool hasNumberOfAlphaElectrons() const
Report whether the alpha-electron count has been set explicitly.
Definition orbitals.h:150
Index getLumo() const
Return the default LUMO index used by spin-agnostic callers.
Definition orbitals.h:105
bool hasBSEUKS_dynamic() const
Definition orbitals.h:693
void setXCFunctionalName(std::string functionalname)
Definition orbitals.h:276
statetype Type() const
Definition qmstate.h:53
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:135
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26