Skip to content

Commit 02bf6a9

Browse files
author
Roberto Di Remigio
committed
Update CHANGELOG.md
1 parent c35b437 commit 02bf6a9

File tree

3 files changed

+21
-26
lines changed

3 files changed

+21
-26
lines changed

CHANGELOG.md

Lines changed: 5 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -10,20 +10,15 @@
1010
- For all solvers, the symmetrization needed to obtain the polarization weights
1111
happens directly on the computed charges, instead of symmetrizing the system
1212
matrices.
13-
- The `IEFSolver` object now stores the unsymmetrized T and R matrices.
14-
The explicit calculation of the inverse of T is thus avoided, the drawback is
15-
that two square matrices of size equal to the cavity size are stored instead
16-
of just one. A [partially pivoted LU decomposition](http://eigen.tuxfamily.org/dox/classEigen_1_1PartialPivLU.html)
17-
is used to solve the linear system.
18-
In the isotropic case no inverses are ever explicitly calculated. In the
19-
anisotropic case, the inverse of the internal S matrix is still needed to
20-
form R.
13+
- The `IEFSolver` object stores the unsymmetrized T^-1R matrices.
14+
A [partially pivoted LU decomposition](http://eigen.tuxfamily.org/dox/classEigen_1_1PartialPivLU.html)
15+
is used to compute T^-1R.
16+
A [robust Cholesky decomposition](http://eigen.tuxfamily.org/dox/classEigen_1_1LDLT.html) is
17+
used to form the R matrix in the anisotropic IEF case.
2118
- The `CPCMSolver` object now stores the scaled, unsymmetrized S matrix. The
2219
explicit calculation and storage of its inverse is thus avoided.
2320
A [robust Cholesky decomposition](http://eigen.tuxfamily.org/dox/classEigen_1_1LDLT.html) is
2421
used to solve the linear equation system.
25-
The symmetrization needed to obtain the polarization weights happens directly
26-
on the computed charges.
2722
- The `hermitivitize` function can now correctly symmetrize vectors.
2823

2924
### Fixed

src/solver/CPCMSolver.cpp

Lines changed: 10 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -31,8 +31,8 @@
3131

3232
#include "Config.hpp"
3333

34-
#include <Eigen/Core>
3534
#include <Eigen/Cholesky>
35+
#include <Eigen/Core>
3636

3737
#include "cavity/Cavity.hpp"
3838
#include "cavity/Element.hpp"
@@ -43,11 +43,11 @@
4343
void CPCMSolver::buildSystemMatrix_impl(const Cavity & cavity, const IGreensFunction & gf_i, const IGreensFunction & gf_o)
4444
{
4545
if (!isotropic_) PCMSOLVER_ERROR("C-PCM is defined only for isotropic environments!", BOOST_CURRENT_FUNCTION);
46-
TIMER_ON("Computing SI");
46+
TIMER_ON("Computing S");
4747
double epsilon = profiles::epsilon(gf_o.permittivity());
48-
SI_ = gf_i.singleLayer(cavity.elements());
49-
SI_ /= (epsilon - 1.0)/(epsilon + correction_);
50-
TIMER_OFF("Computing SI");
48+
S_ = gf_i.singleLayer(cavity.elements());
49+
S_ /= (epsilon - 1.0)/(epsilon + correction_);
50+
TIMER_OFF("Computing S");
5151

5252
// Symmetry-pack
5353
// The number of irreps in the group
@@ -59,10 +59,10 @@ void CPCMSolver::buildSystemMatrix_impl(const Cavity & cavity, const IGreensFunc
5959
// into "block diagonal" when all other manipulations are done.
6060
if (cavity.pointGroup().nrGenerators() != 0) {
6161
TIMER_ON("Symmetry blocking");
62-
symmetryBlocking(SI_, cavity.size(), dimBlock, nrBlocks);
62+
symmetryBlocking(S_, cavity.size(), dimBlock, nrBlocks);
6363
TIMER_OFF("Symmetry blocking");
6464
}
65-
symmetryPacking(blockSI_, SI_, dimBlock, nrBlocks);
65+
symmetryPacking(blockS_, S_, dimBlock, nrBlocks);
6666

6767
built_ = true;
6868
}
@@ -72,12 +72,12 @@ Eigen::VectorXd CPCMSolver::computeCharge_impl(const Eigen::VectorXd & potential
7272
// The potential and charge vector are of dimension equal to the
7373
// full dimension of the cavity. We have to select just the part
7474
// relative to the irrep needed.
75-
int fullDim = SI_.rows();
75+
int fullDim = S_.rows();
7676
Eigen::VectorXd charge = Eigen::VectorXd::Zero(fullDim);
77-
int nrBlocks = blockSI_.size();
77+
int nrBlocks = blockS_.size();
7878
int irrDim = fullDim/nrBlocks;
7979
charge.segment(irrep*irrDim, irrDim) =
80-
- blockSI_[irrep].ldlt().solve(potential.segment(irrep*irrDim, irrDim));
80+
- blockS_[irrep].ldlt().solve(potential.segment(irrep*irrDim, irrDim));
8181
// Symmetrize ASC charge := (charge + charge*)/2
8282
if (hermitivitize_) hermitivitize(charge);
8383

src/solver/CPCMSolver.hpp

Lines changed: 6 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -45,10 +45,10 @@ class IGreensFunction;
4545
* \author Roberto Di Remigio
4646
* \date 2013, 2016
4747
*
48-
* \note We store the unsymmetrized SI matrix and use a robust
48+
* \note We store the unsymmetrized S matrix and use a robust
4949
* Cholesky decomposition to solve for the ASC. The ASC is then
5050
* symmetrized. This avoids computing and storing the inverse explicitly.
51-
* The SI matrix is already scaled by the dielectric factor entering the
51+
* The S matrix is already scaled by the dielectric factor entering the
5252
* definition of the conductor model!
5353
*/
5454

@@ -71,10 +71,10 @@ class CPCMSolver : public PCMSolver
7171
bool hermitivitize_;
7272
/*! Correction for the conductor results */
7373
double correction_;
74-
/*! SI matrix, not symmetry blocked */
75-
Eigen::MatrixXd SI_;
76-
/*! SI matrix, symmetry blocked form */
77-
std::vector<Eigen::MatrixXd> blockSI_;
74+
/*! S matrix, not symmetry blocked */
75+
Eigen::MatrixXd S_;
76+
/*! S matrix, symmetry blocked form */
77+
std::vector<Eigen::MatrixXd> blockS_;
7878

7979
/*! \brief Calculation of the PCM matrix
8080
* \param[in] cavity the cavity to be used

0 commit comments

Comments
 (0)