votca 2024.2-dev
Loading...
Searching...
No Matches
cubefile_writer.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20// Local VOTCA includes
25
26namespace votca {
27namespace xtp {
28
29std::vector<std::vector<double>> CubeFile_Writer::CalculateValues(
30 const Orbitals& orb, QMState state, bool dostateonly,
31 const Regular_Grid& grid) const {
32
33 bool do_amplitude = (state.Type().isSingleParticleState());
34 if (do_amplitude) {
35 Eigen::VectorXd amplitude;
36 if (state.Type() == QMStateType::DQPstate) {
37 Index amplitudeindex = state.StateIdx() - orb.getGWAmin();
38 amplitude = orb.CalculateQParticleAORepresentation().col(amplitudeindex);
39 } else if (state.Type() == QMStateType::LMOstate) {
40 amplitude = orb.getLMOs().col(state.StateIdx());
41 } else {
42 amplitude = orb.MOs().eigenvectors().col(state.StateIdx());
43 }
45 return ampl.IntegrateAmplitude(amplitude);
46 } else {
47 Eigen::MatrixXd mat;
48 if (state.Type().isExciton() && dostateonly) {
49 mat = orb.DensityMatrixWithoutGS(state);
50 } else {
51 mat = orb.DensityMatrixFull(state);
52 }
54 density.IntegrateDensity(mat);
55 return density.getDensities();
56 }
57}
58
59std::vector<double> FlattenValues(
60 const std::vector<std::vector<double>>& gridvalues) {
61 Index size = 0;
62 for (const auto& box : gridvalues) {
63 size += Index(box.size());
64 }
65 std::vector<double> result;
66 result.reserve(size);
67 for (const auto& box : gridvalues) {
68 result.insert(result.end(), box.begin(), box.end());
69 }
70 return result;
71}
72
73void CubeFile_Writer::WriteFile(const std::string& filename,
74 const Orbitals& orb, QMState state,
75 bool dostateonly) const {
76
77 Regular_Grid grid;
78 Eigen::Array3d padding = Eigen::Array3d::Ones() * padding_;
79 AOBasis basis = orb.getDftBasis();
80 XTP_LOG(Log::info, log_) << " Loaded DFT Basis Set " << orb.getDFTbasisName()
81 << std::flush;
82 grid.GridSetup(steps_, padding, orb.QMAtoms(), basis);
83 XTP_LOG(Log::info, log_) << " Calculating Gridvalues " << std::flush;
84 auto temp = CalculateValues(orb, state, dostateonly, grid);
85 std::vector<double> grid_values = FlattenValues(temp);
86 XTP_LOG(Log::info, log_) << " Calculated Gridvalues " << std::flush;
87 bool do_amplitude = (state.Type().isSingleParticleState());
88 std::ofstream out(filename);
89 if (!out.is_open()) {
90 throw std::runtime_error("Bad file handle: " + filename);
91 }
92
93 // write cube header
94 if (state.isTransition()) {
95 out << boost::format("Transition state: %1$s \n") % state.ToString();
96 } else if (do_amplitude) {
97 out << boost::format("%1$s with energy %2$f eV \n") % state.ToString() %
99 } else {
100 if (dostateonly) {
101 out << boost::format(
102 "Difference electron density of excited state %1$s \n") %
103 state.ToString();
104 } else {
105 out << boost::format("Total electron density of %1$s state\n") %
106 state.ToString();
107 }
108 }
109 Eigen::Vector3d start = grid.getStartingPoint();
110 out << "Created by VOTCA-XTP \n";
111 if (do_amplitude) {
112 out << boost::format("-%1$lu %2$f %3$f %4$f \n") % orb.QMAtoms().size() %
113 start.x() % start.y() % start.z();
114 } else {
115 out << boost::format("%1$lu %2$f %3$f %4$f \n") % orb.QMAtoms().size() %
116 start.x() % start.y() % start.z();
117 }
118 Eigen::Array<Index, 3, 1> steps = grid.getSteps();
119 Eigen::Array3d stepsizes = grid.getStepSizes();
120 out << boost::format("%1$d %2$f 0.0 0.0 \n") % steps.x() % stepsizes.x();
121 out << boost::format("%1$d 0.0 %2$f 0.0 \n") % steps.y() % stepsizes.y();
122 out << boost::format("%1$d 0.0 0.0 %2$f \n") % steps.z() % stepsizes.z();
123 tools::Elements elements;
124 for (const QMAtom& atom : orb.QMAtoms()) {
125 double x = atom.getPos().x();
126 double y = atom.getPos().y();
127 double z = atom.getPos().z();
128 std::string element = atom.getElement();
129 Index atnum = elements.getEleNum(element);
130 Index crg = atom.getNuccharge();
131 out << boost::format("%1$d %2$d %3$f %4$f %5$f\n") % atnum % crg % x % y %
132 z;
133 }
134
135 if (do_amplitude) {
136 out << boost::format(" 1 %1$d \n") % (state.StateIdx() + 1);
137 }
138
139 Eigen::TensorMap<Eigen::Tensor<double, 3>> gridvalues(
140 grid_values.data(), steps.z(), steps.y(), steps.x());
141 for (Index ix = 0; ix < steps.x(); ix++) {
142 for (Index iy = 0; iy < steps.y(); iy++) {
143 Index Nrecord = 0;
144 for (Index iz = 0; iz < steps.z(); iz++) {
145 out << boost::format("%1$E ") % gridvalues(iz, iy, ix);
146 Nrecord++;
147 if (Nrecord == 6 || iz == (steps.z() - 1)) {
148 out << "\n";
149 Nrecord = 0;
150 }
151 }
152 }
153 }
154 out.close();
155}
156
157} // namespace xtp
158} // namespace votca
const Eigen::MatrixXd & eigenvectors() const
Definition eigensystem.h:33
information about an element
Definition elements.h:42
Index getEleNum(std::string name)
Definition elements.cc:49
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
std::vector< std::vector< double > > IntegrateAmplitude(const Eigen::VectorXd &amplitude)
void WriteFile(const std::string &filename, const Orbitals &orb, QMState state, bool dostateonly) const
std::vector< std::vector< double > > CalculateValues(const Orbitals &orb, QMState state, bool dostateonly, const Regular_Grid &grid) const
Eigen::Array< Index, 3, 1 > steps_
const std::vector< std::vector< double > > & getDensities() const
double IntegrateDensity(const Eigen::MatrixXd &density_matrix)
container for molecular orbitals
Definition orbitals.h:46
Eigen::MatrixXd CalculateQParticleAORepresentation() const
Definition orbitals.cc:217
Eigen::MatrixXd DensityMatrixWithoutGS(const QMState &state) const
Definition orbitals.cc:112
Eigen::MatrixXd DensityMatrixFull(const QMState &state) const
Definition orbitals.cc:140
const Eigen::MatrixXd & getLMOs() const
Definition orbitals.h:411
double getExcitedStateEnergy(const QMState &state) const
Definition orbitals.cc:451
Index getGWAmin() const
Definition orbitals.h:249
const tools::EigenSystem & MOs() const
Definition orbitals.h:122
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
const std::string & getDFTbasisName() const
Definition orbitals.h:202
const AOBasis & getDftBasis() const
Definition orbitals.h:208
container for QM atoms
Definition qmatom.h:37
bool isExciton() const
Definition qmstate.h:71
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
std::string ToString() const
Definition qmstate.cc:146
Index StateIdx() const
Definition qmstate.h:154
bool isTransition() const
Definition qmstate.h:153
const QMStateType & Type() const
Definition qmstate.h:151
Eigen::Vector3d getStartingPoint() const
Eigen::Array< Index, 3, 1 > getSteps() const
void GridSetup(const Eigen::Array< Index, 3, 1 > &steps, const Eigen::Array3d &padding, const QMMolecule &atoms, const AOBasis &basis)
Eigen::Array3d getStepSizes() const
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
std::vector< double > FlattenValues(const std::vector< std::vector< double > > &gridvalues)
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26