Skip to content

Commit 671e758

Browse files
author
Roberto Di Remigio
committed
Solve two linear systems of equations to get ASC when IEFPCM is used
1 parent fb141c2 commit 671e758

File tree

5 files changed

+57
-9
lines changed

5 files changed

+57
-9
lines changed

CHANGELOG.md

Lines changed: 13 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,6 +6,19 @@
66

77
- The `CPCMSolver` object now stores the scaled, Hermitian, symmetry-adapted S matrix.
88
Polarization weights are then directly computed from the incoming MEP.
9+
- The `IEFSolver` object now stores the non-Hermitian, symmetry-adapted T and R matrices.
10+
The explicit calculation of the inverse of T is thus avoided.
11+
However, two square matrices of size equal to the cavity size are stored instead
12+
of just one. To obtain the polarization weights _two_ linear systems of equations are solved.
13+
A [partially pivoted LU decomposition](http://eigen.tuxfamily.org/dox/classEigen_1_1PartialPivLU.html)
14+
is used to solve the linear system(s).
15+
The strategy used in v1.1.3 suffered from a reduced numerical accuracy, due to the fact that
16+
the polarization weights were not correctly defined.
17+
18+
### Removed
19+
20+
- The `hermitivitize` function will only work correctly on matrices. This
21+
reverts modifications in the previous release.
922

1023
## [v1.1.3] (2016-07-03)
1124

doc/pcmsolver.bib

Lines changed: 17 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -299,3 +299,20 @@ @ARTICLE{DiRemigio2016-nn
299299
doi = "10.1063/1.4943782",
300300
addendum = "Times cited: 0"
301301
}
302+
303+
@ARTICLE{Barone2004-ae,
304+
title = "{Achieving linear-scaling computational cost for the polarizable
305+
continuum model of solvation}",
306+
author = "Barone, Vincenzo and Frisch, Michael J and Scalmani, Giovanni and
307+
Kudin, Konstantin N and Scuseria, Gustavo E and Pomelli, Christian
308+
S",
309+
journal = "Theor. Chem. Acc.",
310+
volume = 111,
311+
number = "2-6",
312+
pages = "90--100",
313+
month = "1~" # mar,
314+
year = 2004,
315+
issn = "1432-881X, 1432-2234",
316+
doi = "10.1007/s00214-003-0527-2"
317+
}
318+

src/solver/IEFSolver.cpp

Lines changed: 11 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -93,8 +93,17 @@ Eigen::VectorXd IEFSolver::computeCharge_impl(const Eigen::VectorXd & potential,
9393
int irrDim = fullDim/nrBlocks;
9494
charge.segment(irrep*irrDim, irrDim) =
9595
- blockTepsilon_[irrep].partialPivLu().solve(blockRinfinity_[irrep] * potential.segment(irrep*irrDim, irrDim));
96-
// Symmetrize ASC charge := (charge + charge*)/2
97-
if (hermitivitize_) hermitivitize(charge);
96+
97+
// Obtain polarization weights
98+
if (hermitivitize_) {
99+
Eigen::VectorXd adj_asc = Eigen::VectorXd::Zero(fullDim);
100+
// Form T^\dagger * v = c
101+
adj_asc.segment(irrep*irrDim, irrDim) = blockTepsilon_[irrep].adjoint().partialPivLu().solve(potential.segment(irrep*irrDim, irrDim));
102+
// Form R^\dagger * c = q^* ("transposed" polarization charges)
103+
adj_asc.segment(irrep*irrDim, irrDim) = - blockRinfinity_[irrep].adjoint() * (adj_asc.segment(irrep*irrDim, irrDim).eval());
104+
// Get polarization weights
105+
charge = 0.5 * (adj_asc + charge.eval());
106+
}
98107

99108
return charge;
100109
}

src/solver/IEFSolver.hpp

Lines changed: 14 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -48,7 +48,20 @@ class IGreensFunction;
4848
* \note We store the non-Hermitian, symmetry-blocked T(epsilon) and Rinfinity matrices.
4949
* The ASC is obtained by multiplying the MEP by Rinfinity and then using a partially
5050
* pivoted LU decomposition of T(epsilon) on the resulting vector.
51-
* This avoids computing and storing the inverse explicitly.
51+
* In case the polarization weights are requested, we use the approach suggested in
52+
* \cite Barone2004-ae.
53+
* First, the adjoint problem is solved:
54+
* \f[
55+
* \mathbf{T}_\varepsilon^\dagger \tilde{v} = v
56+
* \f]
57+
* Also in this case a partially pivoted LU decomposition is used.
58+
* The "transposed" ASC is obtained by the matrix-vector multiplication:
59+
* \f[
60+
* q^* = \mathbf{R}_\infty^\dagger \tilde{v}
61+
* \f]
62+
* Eventually, the two sets of charges are summed and divided by 2
63+
* This avoids computing and storing the inverse explicitly, at the expense of storing
64+
* both T(epsilon) and Rinfinity.
5265
*/
5366

5467
class IEFSolver : public PCMSolver

src/utils/MathUtils.hpp

Lines changed: 2 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -211,17 +211,13 @@ inline void symmetryPacking(std::vector<Eigen::MatrixXd> & blockedMatrix,
211211
* \note We check if a matrix or vector was given, since in the latter
212212
* case we only want the complex conjugation operation to happen.
213213
*/
214-
template <typename Derived>
214+
template <typename Derived>
215215
inline void hermitivitize(Eigen::MatrixBase<Derived> & obj_)
216216
{
217217
// We need to use adjoint().eval() to avoid aliasing issues, see:
218218
// http://eigen.tuxfamily.org/dox/group__TopicAliasing.html
219219
// The adjoint is evaluated explicitly into an intermediate.
220-
if ((obj_.rows() != 1) && (obj_.cols() != 1)) {
221-
obj_ = 0.5 * (obj_ + obj_.adjoint().eval());
222-
} else {
223-
obj_ = 0.5 * (obj_ + obj_.conjugate().eval());
224-
}
220+
obj_ = 0.5 * (obj_ + obj_.adjoint().eval());
225221
}
226222

227223
/*! \fn inline void eulerRotation(Eigen::Matrix3d & R_, const Eigen::Vector3d & eulerAngles_)

0 commit comments

Comments
 (0)