votca 2024.2-dev
Loading...
Searching...
No Matches
bse_operator.h
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#pragma once
21#ifndef VOTCA_XTP_BSE_OPERATOR_H
22#define VOTCA_XTP_BSE_OPERATOR_H
23
24// Local VOTCA includes
25#include "eigen.h"
26#include "matrixfreeoperator.h"
27#include "threecenter.h"
28
29namespace votca {
30namespace xtp {
31
39
40template <Index cqp, Index cx, Index cd, Index cd2>
41class BSE_OPERATOR final : public MatrixFreeOperator {
42
43 public:
44 BSE_OPERATOR(const Eigen::VectorXd& Hd_operator, const TCMatrix_gwbse& Mmn,
45 const Eigen::MatrixXd& Hqp)
46 : epsilon_0_inv_(Hd_operator), Mmn_(Mmn), Hqp_(Hqp) {};
47
49
50 // This method sets up the diagonal of the hermitian BSE hamiltonian.
51 // Otherwise see the matmul function
52 Eigen::VectorXd diagonal() const;
53 /*
54 * This is the main routine for setting up the hermitian parts of the BSE
55 * For the non-hermitian case look at bseoperator_btda.h which combines
56 * the bse_operators. The operator is never set up explicitly, instead only
57 * the product of it with an input matrix is computed. using the template
58 * arguements Different parts of the hamiltonian can be constructed. In
59 * general it is inefficient to set them up independently (unless you need it
60 * for analysis) thus the function combines all parts.
61 */
62 Eigen::MatrixXd matmul(const Eigen::MatrixXd& input) const;
63
64 private:
65 Eigen::VectorXd Hqp_row(Index v1, Index c1) const;
66
72
73 const Eigen::VectorXd& epsilon_0_inv_;
75 const Eigen::MatrixXd& Hqp_;
76};
77
78// type defs for the different operators
81
83
88
89} // namespace xtp
90} // namespace votca
91
92#endif // VOTCA_XTP_BSE_OPERATOR_H
BSE_OPERATOR(const Eigen::VectorXd &Hd_operator, const TCMatrix_gwbse &Mmn, const Eigen::MatrixXd &Hqp)
const Eigen::VectorXd & epsilon_0_inv_
const Eigen::MatrixXd & Hqp_
Eigen::VectorXd Hqp_row(Index v1, Index c1) const
void configure(BSEOperator_Options opt)
const TCMatrix_gwbse & Mmn_
Eigen::MatrixXd matmul(const Eigen::MatrixXd &input) const
BSEOperator_Options opt_
Eigen::VectorXd diagonal() const
BSE_OPERATOR< 1, 0, 0, 0 > HqpOperator
BSE_OPERATOR< 0, 0, 1, 0 > HdOperator
BSE_OPERATOR< 0, 0, 0, 1 > Hd2Operator
BSE_OPERATOR< 1, 2, 1, 0 > SingletOperator_TDA
Definition bse.h:35
BSE_OPERATOR< 1, 0, 1, 0 > TripletOperator_TDA
Definition bse.h:36
BSE_OPERATOR< 0, 2, 0, 1 > SingletOperator_BTDA_B
BSE_OPERATOR< 0, 1, 0, 0 > HxOperator
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26