Merge lp:~nickpapior/siesta/4.1-init-DM into lp:siesta/4.1

Proposed by Nick Papior
Status: Merged
Merge reported by: Nick Papior
Merged at revision: not available
Proposed branch: lp:~nickpapior/siesta/4.1-init-DM
Merge into: lp:siesta/4.1
Diff against target: 278 lines (+189/-19) (has conflicts)
3 files modified
Docs/siesta.tex (+41/-0)
Src/m_new_dm.F90 (+144/-19)
version.info (+4/-0)
Text conflict in version.info
To merge this branch: bzr merge lp:~nickpapior/siesta/4.1-init-DM
Reviewer Review Type Date Requested Status
Alberto Garcia Pending
Review via email: mp+372254@code.launchpad.net

Commit message

Added a new initialization method for DM

The current DM initialization *only* specifies the charge
on the truly diagonal DM entries. I.e. io == col(ind).
However, if one considers that atomic charges are actually
atomic *states* i.e. a state-vector with element in the
diagonal one will also get DM entries in off-diagonal
but io - 1 == col(ind) % no_u which is probably more correct.

I have added two new methods for DM initialization:

1. The DM.Init flag determines one of these methods:
   a) atomic
      the default on-site only atomic density
   b) expand to a state with on-site density
      This enables small auxiliary cell calculations
      to have the density spread amongst overlapping
      orbitals.
      In some systems this seem to behave better
      (graphite/graphene) whereas others it doesn't
      (si64).

2. Added DM.Init.RandomStates <integer> to remove
   N electrons from the initial density matrix
   and add N electrons using a randomized state.
   This state is ensured to be in the non-orthogonal
   basis but is not necessarily orthogonal to any
   other state generated (nor will it be orthogonal
   to the atomic states).
   Currently this method scales orbital coefficients with
   the atomic charges. This will match the atomic charges to
   some degree.
   This is done to limit randoized states with charge in
   high-lying shells which seems dubious in many cases.
   Generally this does not work too good, but may be good
   in some cases I don't know off.

Currently the defaults are as they have always been, but
perhaps the atomic-state would be good in the future?
This needs more heuristics.

Description of the change

I have tried atomic vs. atomic-state on the following tests:

The numbers shows the number of SCF steps needed to converge and the first step dDmax, dHmax

graphite_c6:
atomic: 9, 1.879391, 1.748376
atomic-state: 7, 1.568393, 0.207440

batio3:
atomic: 19, 1.961004, 22.268878
atomic-state: 21, 2.144922, 20.791337

sic-slab:
COMMENT: fun-fact, the SCF cycle indicates that atomic-states is better, but then suddenly atomic converges faster at iSCF 46! History has lots to say here since there is an electric field in this example.
atomic: 69, 1.746005, 13.382736
atomic-state: 144, 2.007246, 13.363516

nanotube-c-5-0:
atomic: 11, 0.741924, 4.242417
atomic-state: 11, 2.037193, 4.270371

To post a comment you must log in.
Revision history for this message
Alberto Garcia (albertog) wrote :

I am trying to understand this a bit better.

The initial DM is used to compute the first H, but *only* the charge density computed from that DM is relevant (H0 does not use the DM, and dhscf only uses rho(r)). The initial charge density then
is the key for initialization. Currently the heuristic is that a superposition of atomic charges is an appropriate method. To improve on that we would need some extra heuristics about how bonding effects in that particular system might influence the charge density. But this cannot be guessed simply from manipulation of 'atomic' concepts.

Does your alternate method give as starting rho(r) the same atomic-charges superposition? I guess not, since the convergence results are different. So you are changing the charge heuristic, but not in a physical way, so I think that the results will be erratic and depend very much on whether your trick leads to a "closer" rho(r) for a particular system.

Besides this, the initial DM will indeed have an influence on the delta_DM found after the first scf step. For example, the 'diagonal' method disregards the nzeta-1 orbitals beyond the first in DM_in, and these are likely to be occupied in DM_out. But delta_DM(scf=1) should be discarded anyway if we are not mixing in the first scf step, as is the case when using atomic initialization.

There are other possibilites for a more targeted initialization: one could extract from a set of calculations in the same or similar system some kind of "deformation charge" estimate and add it to the initial charge density. Note that this is already implemented, but not really well explored.

Or a "physically motivated" initial H could be read. Maybe an extended-Hückel form.

Revision history for this message
Nick Papior (nickpapior) wrote :

1. Agreed on your first paragraph.

2. Regarding the 2nd paragraph.
The *diagonal* Currently the initial DM is made from *only* diagonal elements in the DM. I.e. for supercell systems where orbital io connects with mirror io only io, io is initialized.
For my *state* option the DM is populated as though the identity matrix is the set of eigenvectors (rescaled for the non-orthogonal basis). And thus in this case one will find DM values in both io, io, AND io, mirror io.
My main motivation was that an eigenstate with {1,0,0,0,...} would expand to also mirror io couplings as for calculations without auxiliary supercells.
Do you think this is *very* unphysical? If so then non-auxiliary supercells should always be turned off for diagonal indices since this leads to non atomic initial DM charges. This is substantiated by the sic-slab calculation performed in 3 different settings (in all cases the SCF has only the Gamma-point):
diagonal-forced-auxiliary: 83, 1.820250, 13.346312
diagonal-no-auxiliary: 96, 1.767353, 13.326340
state-forced-auxiliary: 139, 1.767353, 13.326336

Here it can be seen that the initial DM change for the diagonal-no-auxiliary vs. state-forced-auxiliary are the same (as expected via the above).

3. This was my *first* and most basic attempt to try and investigate implications of initial DM. As I have told you previously this could be moved to sisl (I have tried this with success) to do more advanced things quicker, but for now I thought that these would be nice to have for certain cases.
Note that my implementation for *randomized* states was inspired by PW codes which does this (apparently with great success?) and hence I thought it could be used as well here. However, they are much worse, but sometimes using 1 randomized state helps...

4. Agreed that some better heuristic and more thorough investigations could be needed here.

Revision history for this message
Alberto Garcia (albertog) wrote :

I was trying to frame the issue in a more general setting to understand it better.
And my "unphysical" label refers just to not including system-specific 'bonding' effects.

From what you say I might have missed an important more concrete point: is the current
algorithm not giving atomic-charge superpositions in all cases?
I thought we had fixed this (in a zeroth order approximation) in rev 576 of 4.0, which deals with the (arguably pathological) case of severe folding (diagonal images) and gamma-point.

Are there any other cases in which rho(r)[DM_initial] is not a perfect superposition of atomic charges?

Revision history for this message
Nick Papior (nickpapior) wrote :

1. Ok :)

2. No, the current algorithm will *not* give atomic-only-charges IFF 1) the auxiliary is not used AND 2) if the cell is so small that orbitals have interactions with its own mirror orbitals (diagonal images).

If we always *only* want atomic-only charges (i.e. no supercell image overlap charges), then my state solution should not go in.

3. Let me however note that rarely would people do NO k-point sampling for cells needing an auxiliary supercell.

This is probably also a notice for us to consider why we get different results (in SCF cycles) for small cells with/without auxiliary supercells.

lp:~nickpapior/siesta/4.1-init-DM updated
1114. By Nick Papior

Removed the atomic-state method which should never be used

Retained the randomized initialization for experimental
use.

Revision history for this message
Nick Papior (nickpapior) wrote :

Merged on GitLab.

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
=== modified file 'Docs/siesta.tex'
--- Docs/siesta.tex 2019-08-22 07:58:14 +0000
+++ Docs/siesta.tex 2019-09-19 13:04:58 +0000
@@ -5149,6 +5149,25 @@
5149 5149
5150\end{fdflogicalF}5150\end{fdflogicalF}
51515151
5152\begin{fdfentry}{DM!Init}<atomic>
5153 \index{spin!initialization}%
5154
5155 Specify the initial density matrix composition. All these methods
5156 also obey \fdf{DM!InitSpin!AF}.
5157
5158 \begin{fdfoptions}
5159
5160 \option[atomic]%
5161 \fdfindex*{DM!Init!atomic}
5162
5163 Only initialize the diagonal (on-site) density matrices according
5164 to the initial charges on the atomic orbital.
5165
5166 \end{fdfoptions}
5167
5168\end{fdfentry}
5169
5170
5152\begin{fdflogicalF}{DM!InitSpin!AF}5171\begin{fdflogicalF}{DM!InitSpin!AF}
5153 \index{spin!initialization}%5172 \index{spin!initialization}%
5154 \index{ferromagnetic initial DM}%5173 \index{ferromagnetic initial DM}%
@@ -5205,6 +5224,28 @@
52055224
5206\end{fdfentry}5225\end{fdfentry}
52075226
5227\begin{fdfentry}{DM!Init.RandomStates}[integer]<0>
5228
5229 Randomize $N$ number of states in the initial density-matrix. This
5230 may be useful to give the system an initial kick. The randomized
5231 states only has coefficients where the initial atomic charges are
5232 non-zero. This prohibits initial charges in high-lying shells.
5233
5234 This only has an effect if the density matrix is started from an
5235 atomic density and/or when using \fdf{DM!InitSpin}.
5236
5237 In case it is used together with \fdf{DM!InitSpin} it also
5238 randomizes the spin-configuration which may be undesirable.
5239
5240 \note this is currently experimental since the randomized states are
5241 not ensured orthogonal. This flag may be removed in later
5242 revisions with some super seeding functions. If testing this, start
5243 with a value of $1$ to see if it has an effect, any higher numbers
5244 will probably be worse.
5245
5246\end{fdfentry}
5247
5248
5208\begin{fdflogicalT}{DM!AllowReuse}5249\begin{fdflogicalT}{DM!AllowReuse}
52095250
5210 Controls whether density matrix information from previous geometry5251 Controls whether density matrix information from previous geometry
52115252
=== modified file 'Src/m_new_dm.F90'
--- Src/m_new_dm.F90 2019-08-22 07:58:14 +0000
+++ Src/m_new_dm.F90 2019-09-19 13:04:58 +0000
@@ -935,6 +935,8 @@
935935
936 ! *** Local variables936 ! *** Local variables
937 type(block_fdf) :: bfdf937 type(block_fdf) :: bfdf
938 character(len=32) :: method
939 integer :: i_method ! 1 == atomic
938 logical :: init_spin_block940 logical :: init_spin_block
939 type(dData2D) :: arr_2D941 type(dData2D) :: arr_2D
940 ! The pointers for the sparse pattern942 ! The pointers for the sparse pattern
@@ -942,14 +944,21 @@
942 integer :: io, gio, i, ind, jo944 integer :: io, gio, i, ind, jo
943 integer, pointer :: ncol(:), ptr(:), col(:)945 integer, pointer :: ncol(:), ptr(:), col(:)
944 real(dp), pointer :: DM(:,:)946 real(dp), pointer :: DM(:,:)
947
948 method = fdf_get('DM.Init', 'atomic')
949 if ( leqi(method, 'atomic') ) then
950 i_method = 1
951 else
952 call die('m_new_dm: erroneous option: DM.Init [atomic]')
953 end if
945 954
946 ! First we need to figure out the options955 ! First we need to figure out the options
947 if ( spin%DM > 1 ) then956 if ( spin%DM > 1 ) then
948 init_spin_block = fdf_block('DM.InitSpin', bfdf)957 init_spin_block = fdf_block('DM.InitSpin', bfdf)
949 else958 else
950 ! Do not read the initspin block in case of959 ! Do not read the initspin block in case of
951 ! non-polarized calculation960 ! non-polarized calculation
952 init_spin_block = .false.961 init_spin_block = .false.
953 end if962 end if
954963
955 ! Retrieve pointers to data964 ! Retrieve pointers to data
@@ -963,10 +972,10 @@
963 call delete(arr_2D)972 call delete(arr_2D)
964 DM => val(DM_2D)973 DM => val(DM_2D)
965974
966 if ( .not. init_spin_block ) then975 if ( init_spin_block ) then
967 call init_atomic()976 call init_user()
968 else977 else if ( i_method == 1 ) then
969 call init_user()978 call init_atomic()
970 end if979 end if
971980
972 ! We have initialized with atomic information. Correct in case we981 ! We have initialized with atomic information. Correct in case we
@@ -975,18 +984,27 @@
975 ! Gamma-point only. Otherwise S(diagonal) is always 1.0 and the984 ! Gamma-point only. Otherwise S(diagonal) is always 1.0 and the
976 ! simple atomic-orbital superposition works as intended.985 ! simple atomic-orbital superposition works as intended.
977986
987 if ( init_spin_block .or. i_method == 1 ) then
988
978 do io = 1, no_l989 do io = 1, no_l
979 ! Retrieve global orbital index990 ! Retrieve global orbital index
980 gio = index_local_to_global(dit, io)991 gio = index_local_to_global(dit, io)
981 do i = 1 , ncol(io)992 do i = 1 , ncol(io)
982 ind = ptr(io) + i993 ind = ptr(io) + i
983 jo = col(ind)994 jo = col(ind)
984 ! Check for diagonal element995 ! Check for diagonal element
985 if ( gio == jo ) then996 if ( gio == jo ) then
986 DM(ind,:) = DM(ind,:) / S(ind)997 DM(ind,:) = DM(ind,:) / S(ind)
987 endif998 endif
988 end do999 end do
989 end do1000 end do
1001
1002 end if
1003
1004 ! Check if we should randomize some states
1005 ! The default is *no* states are randomized
1006 ! Also, this method does not obey anti-ferro!
1007 call init_randomize()
9901008
991 contains1009 contains
9921010
@@ -1346,9 +1364,116 @@
13461364
1347 end subroutine init_user1365 end subroutine init_user
13481366
1367 subroutine init_randomize()
1368
1369 use parallel, only: Node
1370 use geom_helper, only: ucorb
1371#ifdef MPI
1372 use mpi_siesta
1373#endif
1374
1375 integer :: io, gio, i, ind, jo, i1, i2
1376 integer :: n_random, ir
1377 real(dp), allocatable :: rs(:)
1378 real(dp) :: qtot, qinit, scale
1379#ifdef MPI
1380 integer :: ierr
1381#endif
1382
1383 n_random = fdf_get('DM.Init.RandomStates', 0)
1384
1385 ! If the user does not request any randomization, simply return
1386 if ( n_random == 0 ) return
1387
1388 ! Calculate the number of charges in the current DM (this will be our target!)
1389 qtot = 0._dp
1390 do i2 = 1 , spin%DM
1391 do i1 = 1 , nnz
1392 qtot = qtot + DM(i1,i2) * S(i1)
1393 end do
1394 end do
1395#ifdef MPI
1396 call MPI_AllReduce(qtot, qinit, 1, mpi_double_precision, MPI_SUM, &
1397 MPI_COMM_WORLD, ierr)
1398 qtot = qinit
1399#endif
1400
1401 ! Truncate number of random states (NOT nint)
1402 if ( qtot < n_random ) n_random = int(qtot)
1403
1404 ! Calculate the target *initial DM*
1405 qinit = qtot - n_random
1406 scale = qinit / qtot
1407
1408 ! Rescale DM to remove n_random electrons
1409 do i2 = 1 , spin%spinor
1410 do i1 = 1 , nnz
1411 DM(i1,i2) = DM(i1,i2) * scale
1412 end do
1413 end do
1414
1415 ! Allocate for random states
1416 ! TODO allocate all states and make them orthogonal
1417 allocate(rs(no_u))
1418
1419 do ir = 1, n_random
1420
1421 if ( Node == 0 ) then
1422 call random_number(rs)
1423 ! Scale with atomic weights to not have *random* states
1424 ! on, say, f-shells ;)
1425 do io = 1 , no_u
1426 rs(io) = (rs(io) - 0.5_dp) * DM_atom(io)
1427 end do
1428 end if
1429#ifdef MPI
1430 call MPI_Bcast(rs,no_u,MPI_Double_Precision,0,Mpi_Comm_World,ierr)
1431#endif
1432
1433 ! Calculate inner product to rescale for 1 electron for this state
1434 scale = 0._dp
1435 do io = 1, no_l
1436 gio = index_local_to_global(dit, io)
1437 do ind = ptr(io) + 1 , ptr(io) + ncol(io)
1438 jo = ucorb(col(ind), no_u)
1439 scale = scale + rs(gio) * S(ind) * rs(jo)
1440 end do
1441 end do
1442#ifdef MPI
1443 call MPI_AllReduce(scale, qinit, 1, mpi_double_precision, MPI_SUM, &
1444 MPI_COMM_WORLD, ierr)
1445 scale = qinit
1446#endif
1447 ! Inner product value for getting a charge of 1
1448 scale = sqrt(scale)
1449 if ( spin%spinor > 1 ) then
1450 ! Same state for spin-up/down
1451 rs(:) = rs(:) / (scale * 2)
1452 else
1453 rs(:) = rs(:) / scale
1454 end if
1455
1456 do io = 1, no_l
1457 gio = index_local_to_global(dit, io)
1458 do ind = ptr(io) + 1 , ptr(io) + ncol(io)
1459 jo = ucorb(col(ind), no_u)
1460 DM(ind,1:spin%spinor) = DM(ind,1:spin%spinor) + rs(gio) * rs(jo)
1461 end do
1462 end do
1463 end do
1464
1465 ! Clean-up
1466 deallocate(rs)
1467
1468 if ( IONode ) then
1469 write(*,'(a,i0,a)') 'DM randomized ', n_random, ' states'
1470 end if
1471
1472 end subroutine init_randomize
1473
1349 end subroutine init_DM_atomic1474 end subroutine init_DM_atomic
13501475
1351 1476
1352 subroutine extrapolate_dm_with_coords(DM_history,na_u,xa,sparse_pattern,DMnew)1477 subroutine extrapolate_dm_with_coords(DM_history,na_u,xa,sparse_pattern,DMnew)
1353 use class_Sparsity1478 use class_Sparsity
1354 use class_dData2D1479 use class_dData2D
13551480
=== modified file 'version.info'
--- version.info 2019-09-19 12:58:46 +0000
+++ version.info 2019-09-19 13:04:58 +0000
@@ -1,1 +1,5 @@
1<<<<<<< TREE
1siesta-4.1-11202siesta-4.1-1120
3=======
4siesta-4.1-1112--init-DM-2
5>>>>>>> MERGE-SOURCE

Subscribers

People subscribed via source and target branches