forked from abacusmodeling/abacus-develop
-
Notifications
You must be signed in to change notification settings - Fork 139
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Solving linear equations to evolve the wave function in RT-TDDFT. (#5925
) * Add files via upload * Add files via upload * Add files via upload * Update input-main.md * Update solve_propagation.cpp
- Loading branch information
Showing
7 changed files
with
187 additions
and
12 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
111 changes: 111 additions & 0 deletions
111
source/module_hamilt_lcao/module_tddft/solve_propagation.cpp
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,111 @@ | ||
#include "solve_propagation.h" | ||
|
||
#include <iostream> | ||
|
||
#include "module_base/lapack_connector.h" | ||
#include "module_base/scalapack_connector.h" | ||
|
||
namespace module_tddft | ||
{ | ||
#ifdef __MPI | ||
void solve_propagation(const Parallel_Orbitals* pv, | ||
const int nband, | ||
const int nlocal, | ||
const double dt, | ||
const std::complex<double>* Stmp, | ||
const std::complex<double>* Htmp, | ||
const std::complex<double>* psi_k_laststep, | ||
std::complex<double>* psi_k) | ||
{ | ||
// (1) init A,B and copy Htmp to A & B | ||
std::complex<double>* operator_A = new std::complex<double>[pv->nloc]; | ||
ModuleBase::GlobalFunc::ZEROS(operator_A, pv->nloc); | ||
BlasConnector::copy(pv->nloc, Htmp, 1, operator_A, 1); | ||
|
||
std::complex<double>* operator_B = new std::complex<double>[pv->nloc]; | ||
ModuleBase::GlobalFunc::ZEROS(operator_B, pv->nloc); | ||
BlasConnector::copy(pv->nloc, Htmp, 1, operator_B, 1); | ||
|
||
// ->>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> | ||
// (2) compute operator_A & operator_B by GEADD | ||
// operator_A = Stmp + i*para * Htmp; beta2 = para = 0.25 * dt | ||
// operator_B = Stmp - i*para * Htmp; beta1 = - para = -0.25 * dt | ||
std::complex<double> alpha = {1.0, 0.0}; | ||
std::complex<double> beta1 = {0.0, -0.25 * dt}; | ||
std::complex<double> beta2 = {0.0, 0.25 * dt}; | ||
|
||
ScalapackConnector::geadd('N', | ||
nlocal, | ||
nlocal, | ||
alpha, | ||
Stmp, | ||
1, | ||
1, | ||
pv->desc, | ||
beta2, | ||
operator_A, | ||
1, | ||
1, | ||
pv->desc); | ||
ScalapackConnector::geadd('N', | ||
nlocal, | ||
nlocal, | ||
alpha, | ||
Stmp, | ||
1, | ||
1, | ||
pv->desc, | ||
beta1, | ||
operator_B, | ||
1, | ||
1, | ||
pv->desc); | ||
// ->>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> | ||
// (3) b = operator_B @ psi_k_laststep | ||
std::complex<double>* tmp_b = new std::complex<double>[pv->nloc_wfc]; | ||
ScalapackConnector::gemm('N', | ||
'N', | ||
nlocal, | ||
nband, | ||
nlocal, | ||
1.0, | ||
operator_B, | ||
1, | ||
1, | ||
pv->desc, | ||
psi_k_laststep, | ||
1, | ||
1, | ||
pv->desc_wfc, | ||
0.0, | ||
tmp_b, | ||
1, | ||
1, | ||
pv->desc_wfc); | ||
//get ipiv | ||
int* ipiv = new int[pv->nloc]; | ||
int info = 0; | ||
// (4) solve Ac=b | ||
ScalapackConnector::gesv(nlocal, | ||
nband, | ||
operator_A, | ||
1, | ||
1, | ||
pv->desc, | ||
ipiv, | ||
tmp_b, | ||
1, | ||
1, | ||
pv->desc_wfc, | ||
&info); | ||
|
||
//copy solution to psi_k | ||
BlasConnector::copy(pv->nloc_wfc, tmp_b, 1, psi_k, 1); | ||
|
||
delete []tmp_b; | ||
delete []ipiv; | ||
delete []operator_A; | ||
delete []operator_B; | ||
} | ||
#endif // __MPI | ||
} // namespace module_tddft |
34 changes: 34 additions & 0 deletions
34
source/module_hamilt_lcao/module_tddft/solve_propagation.h
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,34 @@ | ||
#ifndef TD_SOLVE_PROPAGATION_H | ||
#define TD_SOLVE_PROPAGATION_H | ||
|
||
#include "module_basis/module_ao/parallel_orbitals.h" | ||
#include <complex> | ||
|
||
namespace module_tddft | ||
{ | ||
#ifdef __MPI | ||
/** | ||
* @brief solve propagation equation A@c(t+dt) = B@c(t) | ||
* | ||
* @param[in] pv information of parallel | ||
* @param[in] nband number of bands | ||
* @param[in] nlocal number of orbitals | ||
* @param[in] dt time interval | ||
* @param[in] Stmp overlap matrix S(t+dt/2) | ||
* @param[in] Htmp H(t+dt/2) | ||
* @param[in] psi_k_laststep psi of last step | ||
* @param[out] psi_k psi of this step | ||
*/ | ||
void solve_propagation(const Parallel_Orbitals* pv, | ||
const int nband, | ||
const int nlocal, | ||
const double dt, | ||
const std::complex<double>* Stmp, | ||
const std::complex<double>* Htmp, | ||
const std::complex<double>* psi_k_laststep, | ||
std::complex<double>* psi_k); | ||
|
||
#endif | ||
} // namespace module_tddft | ||
|
||
#endif // TD_SOLVE_H |