votca 2026-dev
Loading...
Searching...
No Matches
bse_operator_uks.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2026 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_UKS_H
22#define VOTCA_XTP_BSE_OPERATOR_UKS_H
23
24#include "eigen.h"
25#include "matrixfreeoperator.h"
26#include "threecenter.h"
27
28namespace votca {
29namespace xtp {
30
39
40template <Index cqp, Index cx, Index cd, Index cd2>
42 public:
43 BSE_OPERATOR_UKS(const Eigen::VectorXd& epsilon_0_inv,
44 const TCMatrix_gwbse_spin& Mmn,
45 const Eigen::MatrixXd& Hqp_alpha,
46 const Eigen::MatrixXd& Hqp_beta)
47 : epsilon_0_inv_(epsilon_0_inv),
48 Mmn_(Mmn),
49 Hqp_alpha_(Hqp_alpha),
50 Hqp_beta_(Hqp_beta) {}
51
53
54 Eigen::VectorXd diagonal() const override;
55 Eigen::MatrixXd matmul(const Eigen::MatrixXd& input) const override;
56 Eigen::MatrixXd dense_matrix() const;
57
58 private:
68
70 SpinBlockInfo alpha_;
71 SpinBlockInfo beta_;
72
74
75 const Eigen::VectorXd& epsilon_0_inv_;
77 const Eigen::MatrixXd& Hqp_alpha_;
78 const Eigen::MatrixXd& Hqp_beta_;
79
80 void setup_block(SpinBlockInfo& blk, Index homo, Index offset);
81
82 Eigen::VectorXd Hqp_row(const Eigen::MatrixXd& Hqp, const SpinBlockInfo& blk,
83 Index v1, Index c1) const;
84
85 void add_qp_block(Eigen::MatrixXd& y, const Eigen::MatrixXd& x,
86 const SpinBlockInfo& blk, const Eigen::MatrixXd& Hqp) const;
87
88 void add_exchange_block(Eigen::MatrixXd& y, const Eigen::MatrixXd& x,
89 const SpinBlockInfo& out_blk,
90 const SpinBlockInfo& in_blk,
91 const TCMatrix_gwbse& Mout, const TCMatrix_gwbse& Min,
92 double prefactor) const;
93
94 void add_direct_block(Eigen::MatrixXd& y, const Eigen::MatrixXd& x,
95 const SpinBlockInfo& out_blk,
96 const SpinBlockInfo& in_blk, const TCMatrix_gwbse& Mout,
97 const TCMatrix_gwbse& Min, double prefactor) const;
98
99 void add_direct2_block(Eigen::MatrixXd& y, const Eigen::MatrixXd& x,
100 const SpinBlockInfo& out_blk,
101 const SpinBlockInfo& in_blk,
102 const TCMatrix_gwbse& Mout, const TCMatrix_gwbse& Min,
103 double prefactor) const;
104
105 void add_direct_cross_tda_block(Eigen::MatrixXd& y, const Eigen::MatrixXd& x,
106 const SpinBlockInfo& out_blk,
107 const SpinBlockInfo& in_blk,
108 const TCMatrix_gwbse& Mout,
109 const TCMatrix_gwbse& Min,
110 double prefactor) const;
111};
112
113// TDA A block: Hqp + Hx - Hd
115
116// full BSE B block: Hx - Hd2
118
119// direct-only operators used for perturbative dynamical screening
122
125
126} // namespace xtp
127} // namespace votca
128
129#endif // VOTCA_XTP_BSE_OPERATOR_UKS_H
void add_exchange_block(Eigen::MatrixXd &y, const Eigen::MatrixXd &x, const SpinBlockInfo &out_blk, const SpinBlockInfo &in_blk, const TCMatrix_gwbse &Mout, const TCMatrix_gwbse &Min, double prefactor) const
Eigen::VectorXd diagonal() const override
void configure(BSEOperatorUKS_Options opt)
void add_direct_block(Eigen::MatrixXd &y, const Eigen::MatrixXd &x, const SpinBlockInfo &out_blk, const SpinBlockInfo &in_blk, const TCMatrix_gwbse &Mout, const TCMatrix_gwbse &Min, double prefactor) const
void add_direct_cross_tda_block(Eigen::MatrixXd &y, const Eigen::MatrixXd &x, const SpinBlockInfo &out_blk, const SpinBlockInfo &in_blk, const TCMatrix_gwbse &Mout, const TCMatrix_gwbse &Min, double prefactor) const
void add_direct2_block(Eigen::MatrixXd &y, const Eigen::MatrixXd &x, const SpinBlockInfo &out_blk, const SpinBlockInfo &in_blk, const TCMatrix_gwbse &Mout, const TCMatrix_gwbse &Min, double prefactor) const
BSE_OPERATOR_UKS(const Eigen::VectorXd &epsilon_0_inv, const TCMatrix_gwbse_spin &Mmn, const Eigen::MatrixXd &Hqp_alpha, const Eigen::MatrixXd &Hqp_beta)
Eigen::VectorXd Hqp_row(const Eigen::MatrixXd &Hqp, const SpinBlockInfo &blk, Index v1, Index c1) const
void setup_block(SpinBlockInfo &blk, Index homo, Index offset)
Eigen::MatrixXd matmul(const Eigen::MatrixXd &input) const override
Eigen::MatrixXd dense_matrix() const
void add_qp_block(Eigen::MatrixXd &y, const Eigen::MatrixXd &x, const SpinBlockInfo &blk, const Eigen::MatrixXd &Hqp) const
BSE_OPERATOR_UKS< 0, 1, 0, 1 > ExcitonUKSOperator_BTDA_B
BSE_OPERATOR_UKS< 1, 1, 1, 0 > ExcitonUKSOperator_TDA
BSE_OPERATOR_UKS< 0, 0, 0, 1 > Hd2UKSOperator
BSE_OPERATOR_UKS< 0, 0, 1, 0 > HdUKSOperator
BSE_OPERATOR_UKS< 0, 1, 0, 0 > HxUKSOperator
BSE_OPERATOR_UKS< 1, 0, 0, 0 > QpUKSOperator
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26