[CP2K-user] [CP2K:18553] Re: Evaluation of basis functions over a grid and computation of overlap integrals

Krack Matthias matthias.krack at psi.ch
Fri Mar 17 16:39:48 UTC 2023


Hi

You want to normalize a contracted Gaussian function (cgf) basis set, which is performed in two steps: (1) normalization of the primitive Cartesian Gaussian functions and thereafter (2) normalization of the contracted Gaussian functions as indicated in line 1075<https://github.com/cp2k/cp2k/blob/f26eaef31a9d3f80ca30d8d2f11790a2a072e370/src/aobasis/basis_set_types.F#L1075>. This is all in Cartesian representation. For the spherical orbital representation, you will need a further Cartesian to spherical orbital transformation step. CP2K calculates integrals internally in the Cartesian representation, but the default printout, e.g. for the overlap integral matrix, is in spherical orbitals.
Firstly, I would try to reproduce the normalization for an uncontracted basis set (one primitive function per l only). In the next step, you can try two primitive functions. Normalization is required to keep for instance the electron count (trace of PS).

HTH

Matthias

From: cp2k at googlegroups.com <cp2k at googlegroups.com> on behalf of Aleksandros Sobczyk <sobczykalek at gmail.com>
Date: Friday, 17 March 2023 at 16:59
To: cp2k <cp2k at googlegroups.com>
Subject: Re: [CP2K:18552] Re: Evaluation of basis functions over a grid and computation of overlap integrals
Below is the precise example we are using.
These un-normalized orbitals:
Sr SZV-MOLOPT-SR-GTH SZV-MOLOPT-SR-GTH-q10
1
2 0 1 6 2 1
     7.290111894735  0.069364270475 -0.016182746349  0.035659445929
     2.536776771327 -0.571246927373  0.158928639982 -0.195822349727
     1.283099546928  0.167836311459 -0.041157757852  0.260320252229
     0.532449841650  0.904733330629 -0.431882417196  0.555386362294
     0.211628059408  0.250907816117 -0.073959284415  0.267635587013
     0.050841303698  0.007135721199  0.825600560103  0.001214128781
after the normalization (cp2k output) become as follows:
Sr SZV-MOLOPT-SR-GTH SZV-MOLOPT-SR-GTH-q10
1
2 0 1 6 2 1
 7.290111894735       0.221325      -0.071449       0.684627
 2.536776771327      -0.825808       0.317914      -1.004793
 1.283099546928       0.145521      -0.049379       0.569764
 0.532449841650       0.405577      -0.267899       0.404860
 0.211628059408       0.056304      -0.022965       0.061570
 0.050841303698       0.000549       0.087969       0.000047

None of the two aforementioned normalization procedures can reproduce this.

On Friday, March 17, 2023 at 4:54:28 PM UTC+1 Aleksandros Sobczyk wrote:
Update on this:
There appear to be two (maybe more) normalization procedures for orbital coefficients:
1) https://github.com/cp2k/cp2k/blob/f26eaef31a9d3f80ca30d8d2f11790a2a072e370/src/aobasis/basis_set_types.F#L1099
2) https://github.com/cp2k/cp2k/blob/f26eaef31a9d3f80ca30d8d2f11790a2a072e370/src/atom_types.F#L2374
The results of 1) seem to be closer to the output of cp2k (for s-orbitals). For p-orbitals it starts to diverge.

It would be great if someone can confirm which normalization method is used internally for contracted Gaussians.
(what is the goal of the normalization? are the primitive Gaussians normalized to integrate to 1? is it something more advanced?)

Best regards,
Aleksandros Sobczyk
On Monday, March 13, 2023 at 11:00:04 AM UTC+1 Aleksandros Sobczyk wrote:
Dear prof. Hutter and CP2K developers,

The example helped to figure out how to reproduce the overlap matrix for our systems. We had to apply the following two steps:

1) For the basis sets, we had to use the normalized coefficients that are found inside the output file of CP2K. The un-normalized coefficients from the original BASIS_SET file did not work (as you suggested). For now, we can directly use the CP2K output to derive these normalized coefficients for our experiments, but it would be also helpful to know how to reproduce them (i.e., what is the normalization procedure).

2) We also had to scale all the exponents of all the basis sets by a factor of (1/0.529)^2 (to convert Angstrom to atomic units). This was an arbitrary guess, but without it the overlaps do not match.

Could you confirm that this is the correct way to compute the orbital overlaps? (and if not, propose the correct way).
It would also be helpful if these details can be documented, e.g. in this page: https://www.cp2k.org/basis_sets

Best regards and thank you again for the feedback,
Aleksandros Sobczyk
On Friday, March 10, 2023 at 6:35:28 PM UTC+1 Aleksandros Sobczyk wrote:
Dear Prof. Hutter,

Thank you very much for your reply and for the example.
Let us investigate it with my colleagues and see if we can resolve our problem.

Best regards,
Aleksandros
On Friday, March 10, 2023 at 3:18:28 PM UTC+1 Jürg Hutter wrote:
Hi

are you assuming normalized or un-normalized Gaussians?
The basis set input in CP2K uses (like all QC codes) normalized Gaussians.
Internally, CP2K works with un-normalized Cartesian Gaussians, i.e. the
coefficients are adapted at the beginning of the calculation.

I have attached a simple example where you can play with the basis set
in the input and the overlap matrix is printed.

regards
JH

________________________________________
From: cp... at googlegroups.com <cp... at googlegroups.com> on behalf of Aleksandros Sobczyk <sobcz... at gmail.com>
Sent: Friday, March 10, 2023 2:26 PM
To: cp2k
Subject: [CP2K:18530] Re: Evaluation of basis functions over a grid and computation of overlap integrals

Update: I have also calculated several overlap integrals analytically (for s-orbitals which are simpler), they still don't match the values of the S matrix.
Any feedback would be greatly appreciated.

On Wednesday, March 8, 2023 at 12:57:02 PM UTC+1 Aleksandros Sobczyk wrote:
Hello,

I have a set of atoms in real-space and the corresponding SZV basis sets.
I want to evaluate each basis function over a grid of points in the cell.
E.g., I have a grid of 3d points [r1, r2, ..., rk] and I want to evaluate each
Φj(r1), Φj(r2), ... Φj(rk)
As a test, I tried to numerically integrate Φj * conj(Φj) over the grid that it was evaluated, and compare the result with the corresponding entry S[j, j] of the overlap matrix that
is returned by CP2K.
Unfortunately my integral differs substantially from the element S[j, j], so I am doing something wrong.

Can we find somewhere more detailed documentation on the precise mathematical formulation of the basis sets, and also on the specific algorithms that are used by CP2K to compute the overlap integrals?
(So far I have followed as precisely as possible the following page: https://www.cp2k.org/basis_sets
but it is still missing information, e.g. are the coefficients normalized? do we assume that the spherical harmonics include the phase factor? etc.)

Thanks a lot in advance!
Aleksandros

--
You received this message because you are subscribed to the Google Groups "cp2k" group.
To unsubscribe from this group and stop receiving emails from it, send an email to cp2k+uns... at googlegroups.com<mailto:cp2k+uns... at googlegroups.com>.
To view this discussion on the web visit https://groups.google.com/d/msgid/cp2k/fa61a2bc-2c9c-4771-8606-445d46af9148n%40googlegroups.com<https://groups.google.com/d/msgid/cp2k/fa61a2bc-2c9c-4771-8606-445d46af9148n%40googlegroups.com?utm_medium=email&utm_source=footer>.
--
You received this message because you are subscribed to the Google Groups "cp2k" group.
To unsubscribe from this group and stop receiving emails from it, send an email to cp2k+unsubscribe at googlegroups.com<mailto:cp2k+unsubscribe at googlegroups.com>.
To view this discussion on the web visit https://groups.google.com/d/msgid/cp2k/073f7bb7-8358-47b8-ba7c-06fe5075ff6cn%40googlegroups.com<https://groups.google.com/d/msgid/cp2k/073f7bb7-8358-47b8-ba7c-06fe5075ff6cn%40googlegroups.com?utm_medium=email&utm_source=footer>.

-- 
You received this message because you are subscribed to the Google Groups "cp2k" group.
To unsubscribe from this group and stop receiving emails from it, send an email to cp2k+unsubscribe at googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/cp2k/ZRAP278MB08271DA6C4EF55654EE56F8AF4BD9%40ZRAP278MB0827.CHEP278.PROD.OUTLOOK.COM.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://lists.cp2k.org/archives/cp2k-user/attachments/20230317/d7314567/attachment-0001.htm>


More information about the CP2K-user mailing list