votca 2024.1-dev
Loading...
Searching...
No Matches
aopotential.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
22
23namespace votca {
24namespace xtp {
25
26template <class T>
27Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> AOPotential<T>::Fill(
28 const AOBasis& aobasis) const {
29 typedef Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> MatrixXcdd;
30
31 MatrixXcdd result =
32 MatrixXcdd::Zero(aobasis.AOBasisSize(), aobasis.AOBasisSize());
33 // AOMatrix is symmetric, restrict explicit calculation of lower triangular
34 // matrix
35#pragma omp parallel for schedule(guided)
36 for (Index col = 0; col < aobasis.getNumofShells(); col++) {
37 const AOShell& shell_col = aobasis.getShell(col);
38 Index col_start = shell_col.getStartIndex();
39 for (Index row = col; row < aobasis.getNumofShells(); row++) {
40 const AOShell& shell_row = aobasis.getShell(row);
41 Index row_start = shell_row.getStartIndex();
42 // figure out the submatrix
43 Eigen::Block<MatrixXcdd> block = result.block(
44 row_start, col_start, shell_row.getNumFunc(), shell_col.getNumFunc());
45 // Fill block
46 FillBlock(block, shell_row, shell_col);
47 }
48 }
49 // Fill whole matrix by copying
50 return result.template selfadjointView<Eigen::Lower>();
51}
52
53template class AOPotential<double>;
54template class AOPotential<std::complex<double> >;
55
56} // namespace xtp
57} // namespace votca
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Index AOBasisSize() const
Definition aobasis.h:46
const AOShell & getShell(Index idx) const
Definition aobasis.h:52
Index getNumofShells() const
Definition aobasis.h:56
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > Fill(const AOBasis &aobasis) const
Index getStartIndex() const
Definition aoshell.h:105
Index getNumFunc() const
Definition aoshell.h:103
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26