votca 2024-dev
Loading...
Searching...
No Matches
aobasis.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 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
21#include "votca/xtp/aobasis.h"
22#include "votca/xtp/basisset.h"
25
26namespace votca {
27namespace xtp {
28
30 Index n = 0;
31 for (const auto& shell : aoshells_) {
32 n = std::max(static_cast<Index>(shell.getL()), n);
33 }
34 return n;
35}
36
38 Index n = 0;
39 for (const auto& shell : aoshells_) {
40 n = std::max(shell.getSize(), n);
41 }
42 return n;
43}
44
45std::vector<Index> AOBasis::getMapToBasisFunctions() const {
46 std::vector<Index> result;
47 result.reserve(aoshells_.size());
48
49 Index n = 0;
50 for (const auto& shell : aoshells_) {
51 result.push_back(n);
52 n += shell.getNumFunc();
53 }
54 return result;
55}
56
57AOShell& AOBasis::addShell(const Shell& shell, const QMAtom& atom,
58 Index startIndex) {
59 aoshells_.push_back(AOShell(shell, atom, startIndex));
60 return aoshells_.back();
61}
62
63const std::vector<const AOShell*> AOBasis::getShellsofAtom(Index AtomId) const {
64 std::vector<const AOShell*> result;
65 for (const auto& aoshell : aoshells_) {
66 if (aoshell.getAtomIndex() == AtomId) {
67 result.push_back(&aoshell);
68 }
69 }
70 return result;
71}
72
73void AOBasis::add(const AOBasis& other) {
74 Index atomindex_offset = Index(FuncperAtom_.size());
75 for (AOShell shell : other) {
76 shell.atomindex_ += atomindex_offset;
77 shell.startIndex_ = AOBasisSize_;
78 AOBasisSize_ += shell.getNumFunc();
79 aoshells_.push_back(shell);
80 }
81
83}
84
85void AOBasis::Fill(const BasisSet& bs, const QMMolecule& atoms) {
86 clear();
87 name_ = bs.Name();
88 // loop over atoms
89 for (const QMAtom& atom : atoms) {
90 const std::string& name = atom.getElement();
91 const Element& element = bs.getElement(name);
92 for (const Shell& shell : element) {
93 Index numfuncshell = NumFuncShell(shell.getL());
94 AOShell& aoshell = addShell(shell, atom, AOBasisSize_);
95 AOBasisSize_ += numfuncshell;
96 for (const GaussianPrimitive& gaussian : shell) {
97 aoshell.addGaussian(gaussian);
98 }
99 aoshell.CalcMinDecay();
100 aoshell.normalizeContraction();
101 }
102 }
104 return;
105}
106
108 FuncperAtom_.clear();
109 Index currentIndex = -1;
110 for (const auto& shell : aoshells_) {
111 if (shell.getAtomIndex() == currentIndex) {
112 FuncperAtom_[shell.getAtomIndex()] += shell.getNumFunc();
113 } else {
114 currentIndex = shell.getAtomIndex();
115 FuncperAtom_.push_back(shell.getNumFunc());
116 }
117 }
118}
119std::vector<libint2::Shell> AOBasis::GenerateLibintBasis() const {
120 std::vector<libint2::Shell> libintshells;
121 libintshells.reserve(aoshells_.size());
122 for (const auto& shell : aoshells_) {
123 libintshells.push_back(shell.LibintShell());
124 }
125 return libintshells;
126}
127
129 for (AOShell& shell : aoshells_) {
130 shell.pos_ = mol[shell.getAtomIndex()].getPos();
131 }
132}
133
135 name_ = "";
136 aoshells_.clear();
137 FuncperAtom_.clear();
138 AOBasisSize_ = 0;
139}
140
142 w(name_, "name");
143 w(AOBasisSize_, "basissize");
144 Index numofprimitives = 0;
145 for (const auto& shell : aoshells_) {
146 numofprimitives += shell.getSize();
147 }
148
149 CptTable table =
150 w.openTable<AOGaussianPrimitive>("Contractions", numofprimitives);
151
152 std::vector<AOGaussianPrimitive::data> dataVec(numofprimitives);
153 Index i = 0;
154 for (const auto& shell : aoshells_) {
155 for (const auto& gaussian : shell) {
156 gaussian.WriteData(dataVec[i], shell);
157 i++;
158 }
159 }
160 table.write(dataVec);
161}
162
164 clear();
165 r(name_, "name");
166 r(AOBasisSize_, "basissize");
167 if (AOBasisSize_ > 0) {
168
169 CptTable table = r.openTable<AOGaussianPrimitive>("Contractions");
170 std::vector<AOGaussianPrimitive::data> dataVec(table.numRows());
171 table.read(dataVec);
172 Index laststartindex = -1;
173 for (std::size_t i = 0; i < table.numRows(); ++i) {
174 if (dataVec[i].startindex != laststartindex) {
175 aoshells_.push_back(AOShell(dataVec[i]));
176 laststartindex = dataVec[i].startindex;
177 } else {
178 aoshells_.back().gaussians_.push_back(AOGaussianPrimitive(dataVec[i]));
179 }
180 }
181
183 }
184 for (auto& shell : aoshells_) {
185 shell.CalcMinDecay();
186 }
187}
188
189std::ostream& operator<<(std::ostream& out, const AOBasis& aobasis) {
190 out << "Name:" << aobasis.Name() << "\n";
191 out << " Functions:" << aobasis.AOBasisSize()
192 << " Shells:" << aobasis.getNumofShells() << "\n";
193 for (const auto& shell : aobasis) {
194 out << shell;
195 }
196 return out;
197}
198
199} // namespace xtp
200} // namespace votca
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Index getMaxL() const
Definition aobasis.cc:29
void UpdateShellPositions(const QMMolecule &mol)
Definition aobasis.cc:128
const std::vector< const AOShell * > getShellsofAtom(Index AtomId) const
Definition aobasis.cc:63
Index AOBasisSize() const
Definition aobasis.h:46
void Fill(const BasisSet &bs, const QMMolecule &atoms)
Definition aobasis.cc:85
void add(const AOBasis &other)
Definition aobasis.cc:73
std::vector< Index > getMapToBasisFunctions() const
Definition aobasis.cc:45
const std::string & Name() const
Definition aobasis.h:81
Index getNumofShells() const
Definition aobasis.h:56
Index getMaxNprim() const
Definition aobasis.cc:37
void WriteToCpt(CheckpointWriter &w) const
Definition aobasis.cc:141
void ReadFromCpt(CheckpointReader &r)
Definition aobasis.cc:163
std::string name_
Definition aobasis.h:92
AOShell & addShell(const Shell &shell, const QMAtom &atom, Index startIndex)
Definition aobasis.cc:57
std::vector< Index > FuncperAtom_
Definition aobasis.h:101
std::vector< AOShell > aoshells_
Definition aobasis.h:99
std::vector< libint2::Shell > GenerateLibintBasis() const
Definition aobasis.cc:119
void normalizeContraction()
Definition aoshell.cc:81
void addGaussian(const GaussianPrimitive &gaussian)
Definition aoshell.h:142
const Eigen::Vector3d & getPos() const
const Element & getElement(std::string element_type) const
Definition basisset.cc:212
const std::string & Name() const
Definition basisset.h:152
CptTable openTable(const std::string &name)
CptTable openTable(const std::string &name, std::size_t nRows, bool compact=false)
void write(void *buffer, const std::size_t &startIdx, const std::size_t &endIdx)
void read(void *buffer, const std::size_t &startIdx, const std::size_t &endIdx)
container for QM atoms
Definition qmatom.h:37
std::ostream & operator<<(std::ostream &out, const Correlate &c)
Definition correlate.h:53
Index NumFuncShell(L l)
Definition basisset.cc:100
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26