Merge lp:~nickpapior/siesta/merge-OSSO-class into lp:~albertog/siesta/merge-OSSO

Proposed by Nick Papior
Status: Merged
Merged at revision: 701
Proposed branch: lp:~nickpapior/siesta/merge-OSSO-class
Merge into: lp:~albertog/siesta/merge-OSSO
Diff against target: 704 lines (+173/-177)
8 files modified
Src/compute_energies.F90 (+28/-28)
Src/final_H_f_stress.F (+12/-7)
Src/nlefsm.f (+13/-11)
Src/setup_H0.F (+25/-26)
Src/setup_hamiltonian.F (+76/-85)
Src/sparse_matrices.F (+7/-7)
Src/state_init.F (+11/-12)
version.info (+1/-1)
To merge this branch: bzr merge lp:~nickpapior/siesta/merge-OSSO-class
Reviewer Review Type Date Requested Status
Alberto Garcia Approve
Review via email: mp+343801@code.launchpad.net

Commit message

Transferred H0_offsiteSO to H_so_off_2D as a class-object.

Fixed some possible inconsistencies in the conversion between
complex and real/imaginary part without kind specifications.

Clarified usage in setup_hamiltonian. Now the structure has been
re-assigned as prior to the OSSO implementation.

Description of the change

Transferred H0_offsiteSO to H_so_off_2D as a class-object.

Fixed some possible inconsistencies in the conversion between
complex and real/imaginary part without kind specifications.

Clarified usage in setup_hamiltonian. Now the structure has been
re-assigned as prior to the OSSO implementation.

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

Great. Thanks!

review: Approve

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
=== modified file 'Src/compute_energies.F90'
--- Src/compute_energies.F90 2018-04-23 07:36:55 +0000
+++ Src/compute_energies.F90 2018-04-23 11:10:59 +0000
@@ -48,8 +48,7 @@
48 rmaxo, no_l48 rmaxo, no_l
49 use m_ntm, only: ntm49 use m_ntm, only: ntm
5050
51 use m_spin, only: spin 51 use m_spin, only: spin
52 use sparse_matrices, only: H0_offsiteSO
53 use parallel, only: IONode 52 use parallel, only: IONode
5453
55 use m_dipol, only: dipol54 use m_dipol, only: dipol
@@ -300,10 +299,15 @@
300 use files, only : filesOut_t ! derived type for output file names299 use files, only : filesOut_t ! derived type for output file names
301 use class_dSpData1D, only : val300 use class_dSpData1D, only : val
302 use class_dSpData2D, only : val301 use class_dSpData2D, only : val
303 use sparse_matrices, only: H_kin_1D, H_vkb_1D, H_so_2D302 use class_zSpData2D, only : val
303 use sparse_matrices, only: H_kin_1D, H_vkb_1D
304 use sparse_matrices, only: H_so_on_2D, H_so_off_2D
305
304306
305 type(filesOut_t) :: filesOut ! blank output file names307 type(filesOut_t) :: filesOut ! blank output file names
306 real(dp), pointer :: H_vkb(:), H_kin(:), H_so(:,:)308 real(dp), pointer :: H_vkb(:), H_kin(:), H_so_on(:,:)
309 complex(dp), pointer :: H_so_off(:,:)
310
307311
308 complex(dp) :: Hc, Dc312 complex(dp) :: Hc, Dc
309 real(dp) :: dummy_stress(3,3), dummy_fa(1,1)313 real(dp) :: dummy_stress(3,3), dummy_fa(1,1)
@@ -356,41 +360,37 @@
356 Eso = 0._dp360 Eso = 0._dp
357361
358 if ( spin%SO_offsite ) then362 if ( spin%SO_offsite ) then
363 H_so_off => val(H_so_off_2D)
359364
360 ! The computation of the trace is different here, as H0_offsiteSO has365 ! The computation of the trace is different here, as H_so_off has
361 ! a different structure from H and the DM.366 ! a different structure from H and the DM.
362 do io = 1, maxnh367 do io = 1, maxnh
363368
364!-------- Eso(u,u)369!-------- Eso(u,u)
365 Dc = cmplx(Dscf(io,1),Dscf(io,5),kind=dp)370 Dc = cmplx(Dscf(io,1),Dscf(io,5), dp)
366 Hc = H0_offsiteSO(io,1)371 Eso = Eso + real( H_so_off(io,1)*Dc, dp)
367 Eso = Eso + real( Hc*Dc )
368!-------- Eso(d,d)372!-------- Eso(d,d)
369 Dc = cmplx(Dscf(io,2),Dscf(io,6),kind=dp)373 Dc = cmplx(Dscf(io,2),Dscf(io,6), dp)
370 Hc = H0_offsiteSO(io,2)374 Eso = Eso + real( H_so_off(io,2)*Dc, dp)
371 Eso = Eso + real( Hc*Dc )
372!-------- Eso(u,d)375!-------- Eso(u,d)
373 Dc = cmplx(Dscf(io,3),Dscf(io,4),kind=dp)376 Dc = cmplx(Dscf(io,3),Dscf(io,4), dp)
374 Hc = H0_offsiteSO(io,4)377 Eso = Eso + real( H_so_off(io,4)*Dc, dp)
375 Eso = Eso + real( Hc*Dc )
376!-------- Eso(d,u)378!-------- Eso(d,u)
377 Dc = cmplx(Dscf(io,7),-Dscf(io,8),kind=dp)379 Dc = cmplx(Dscf(io,7),-Dscf(io,8), dp)
378 Hc = H0_offsiteSO(io,3)380 Eso = Eso + real( H_so_off(io,3)*Dc, dp)
379 Eso = Eso + real( Hc*Dc )381
380382 end do
381 enddo383
382384 else if ( spin%SO_onsite ) then
383 endif
384
385 if ( spin%SO_onsite ) then
386 ! Sadly some compilers (g95), does385 ! Sadly some compilers (g95), does
387 ! not allow bounds for pointer assignments :(386 ! not allow bounds for pointer assignments :(
388 H_so => val(H_so_2D)387 H_so_on => val(H_so_on_2D)
389 do io = 1,maxnh388 do io = 1,maxnh
390 Eso = Eso + H_so(io,1)*Dscf(io,7) + H_so(io,2)*Dscf(io,8) &389 Eso = Eso + H_so_on(io,1)*Dscf(io,7) + H_so_on(io,2)*Dscf(io,8) &
391 + H_so(io,5)*Dscf(io,3) + H_so(io,6)*Dscf(io,4) &390 + H_so_on(io,5)*Dscf(io,3) + H_so_on(io,6)*Dscf(io,4) &
392 - H_so(io,3)*Dscf(io,5) - H_so(io,4)*Dscf(io,6)391 - H_so_on(io,3)*Dscf(io,5) - H_so_on(io,4)*Dscf(io,6)
393 end do392 end do
393
394 end if394 end if
395395
396#ifdef MPI396#ifdef MPI
397397
=== modified file 'Src/final_H_f_stress.F'
--- Src/final_H_f_stress.F 2018-04-20 09:41:55 +0000
+++ Src/final_H_f_stress.F 2018-04-23 11:10:59 +0000
@@ -21,8 +21,9 @@
21 use siesta_options, only: recompute_H_after_scf21 use siesta_options, only: recompute_H_after_scf
22 use sparse_matrices, only: numh, listh, listhptr22 use sparse_matrices, only: numh, listh, listhptr
23 use sparse_matrices, only: H, S, Dscf, Escf, maxnh, xijo23 use sparse_matrices, only: H, S, Dscf, Escf, maxnh, xijo
24 use sparse_matrices, only: H_ldau_2D24 use sparse_matrices, only: H_ldau_2D, H_so_off_2D
25 use class_dSpData2D, only: val25 use class_dSpData2D, only: val
26 use class_zSpData2D, only: val
2627
27 use siesta_geom28 use siesta_geom
28 use atomlist, only: no_u, iaorb, iphkb, qtot, indxuo, datm, 29 use atomlist, only: no_u, iaorb, iphkb, qtot, indxuo, datm,
@@ -30,7 +31,7 @@
30 & rmaxo, no_l, iza31 & rmaxo, no_l, iza
31 use metaforce, only: lMetaForce, meta32 use metaforce, only: lMetaForce, meta
32 use molecularmechanics, only : twobody33 use molecularmechanics, only : twobody
33 use m_nlefsm, only: nlefsm, nlefsm_offsiteSO34 use m_nlefsm, only: nlefsm, nlefsm_SO_off
34 use m_overfsm, only: overfsm35 use m_overfsm, only: overfsm
35 use m_kinefsm, only: kinefsm36 use m_kinefsm, only: kinefsm
36 use m_naefs, only: naefs37 use m_naefs, only: naefs
@@ -91,6 +92,8 @@
91 real(dp), pointer :: H_tmp(:,:) => null()92 real(dp), pointer :: H_tmp(:,:) => null()
92 real(dp), pointer :: S_tmp(:) => null() 93 real(dp), pointer :: S_tmp(:) => null()
93 real(dp), pointer :: H_ldau(:,:)94 real(dp), pointer :: H_ldau(:,:)
95 complex(dp), pointer :: H_so_off(:,:)
96
94#ifdef FINAL_CHECK_HS97#ifdef FINAL_CHECK_HS
95 real(dp) :: diff_H, diff_S98 real(dp) :: diff_H, diff_S
96#endif99#endif
@@ -176,13 +179,14 @@
176 call globalize_sum(Enl,buffer1)179 call globalize_sum(Enl,buffer1)
177 Enl = buffer1180 Enl = buffer1
178#endif181#endif
179 else 182 else
180 call nlefsm_offsiteSO(scell, na_u, na_s, isa, xa, indxua,183 H_so_off => val(H_so_off_2D)
184 call nlefsm_SO_off(scell, na_u, na_s, isa, xa, indxua,
181 & maxnh, maxnh, lasto, lastkb, iphorb, iphKB,185 & maxnh, maxnh, lasto, lastkb, iphorb, iphKB,
182 & numh, listhptr, listh, numh, listhptr, listh,186 & numh, listhptr, listh, numh, listhptr, listh,
183 & spin%Grid,187 & spin%Grid,
184 & Enl, Eso, fal,188 & Enl, Eso, fal,
185 & stressl, H_tmp,189 & stressl, H_tmp, H_so_off,
186 & matrix_elements_only=.false.)190 & matrix_elements_only=.false.)
187#ifdef MPI191#ifdef MPI
188! Global reduction of energy terms192! Global reduction of energy terms
@@ -193,7 +197,7 @@
193#endif197#endif
194 endif198 endif
195199
196 if ( spin%SO_onsite) then200 if ( spin%SO_onsite ) then
197 call spinorb(no_u,no_l,iaorb,iphorb,isa,indxuo,201 call spinorb(no_u,no_l,iaorb,iphorb,isa,indxuo,
198 & maxnh,numh,listhptr,listh,Dscf,H_tmp(:,3:),Eso)202 & maxnh,numh,listhptr,listh,Dscf,H_tmp(:,3:),Eso)
199 endif203 endif
@@ -292,7 +296,8 @@
292 & diff_S296 & diff_S
293 end if297 end if
294#endif298#endif
295 299
300 ! TODO I am not sure this works with SO_offsite?
296 if (recompute_H_after_scf) then301 if (recompute_H_after_scf) then
297 if (ionode) then302 if (ionode) then
298 write(6,"(a)") ":!: Updating H after scf cycle" //303 write(6,"(a)") ":!: Updating H after scf cycle" //
299304
=== modified file 'Src/nlefsm.f'
--- Src/nlefsm.f 2018-04-23 09:26:04 +0000
+++ Src/nlefsm.f 2018-04-23 11:10:59 +0000
@@ -8,11 +8,10 @@
8 module m_nlefsm8 module m_nlefsm
99
10 use precision, only : dp10 use precision, only : dp
11 use sparse_matrices, only: H0_offsiteSO
1211
13 implicit none12 implicit none
1413
15 public :: nlefsm, nlefsm_offsiteSO, calc_Vj_offsiteSO14 public :: nlefsm, nlefsm_SO_off
1615
17 private16 private
1817
@@ -445,11 +444,12 @@
445C nlefsm_offsiteSO calculates the KB elements to the total Hamiltonian 444C nlefsm_offsiteSO calculates the KB elements to the total Hamiltonian
446C (including the SO contribution)445C (including the SO contribution)
447C when Off-Site Spin Orbit is included in the calculation 446C when Off-Site Spin Orbit is included in the calculation
448 subroutine nlefsm_offsiteSO( scell, nua, na, isa, xa, indxua,447 subroutine nlefsm_SO_off( scell, nua, na, isa, xa, indxua,
449 . maxnh, maxnd, lasto, lastkb, iphorb, 448 . maxnh, maxnd, lasto, lastkb, iphorb,
450 . iphKB, numd, listdptr, listd, numh, 449 . iphKB, numd, listdptr, listd, numh,
451 . listhptr, listh, nspin, Enl, Enl_offsiteSO, 450 . listhptr, listh, nspin, Enl, Enl_offsiteSO,
452 . fa, stress, H0 , matrix_elements_only)451 . fa, stress, H0 , H0_off,
452 & matrix_elements_only)
453453
454454
455C *********************************************************************455C *********************************************************************
@@ -496,7 +496,8 @@
496C ******************* INPUT and OUTPUT *********************************496C ******************* INPUT and OUTPUT *********************************
497C real*8 fa(3,na) : NL forces (added to input fa)497C real*8 fa(3,na) : NL forces (added to input fa)
498C real*8 stress(3,3) : NL stress (added to input stress)498C real*8 stress(3,3) : NL stress (added to input stress)
499C real*8 H(maxnh,nspin) : NL Hamiltonian (added to input H)499C real*8 H0(maxnh,nspin) : NL Hamiltonian (added to input H)
500C complex*16 H0_off(maxnh,nspin) : NL off-site Hamiltonian (added to input H)
500C **************************** OUTPUT *********************************501C **************************** OUTPUT *********************************
501C real*8 Enl : NL energy502C real*8 Enl : NL energy
502C *********************************************************************503C *********************************************************************
@@ -525,7 +526,8 @@
525526
526 real(dp), intent(in) :: scell(3,3), xa(3,na)527 real(dp), intent(in) :: scell(3,3), xa(3,na)
527 real(dp), intent(inout) :: fa(3,nua), stress(3,3)528 real(dp), intent(inout) :: fa(3,nua), stress(3,3)
528 real(dp), intent(inout) :: H0(maxnh) 529 real(dp), intent(inout) :: H0(maxnh)
530 complex(dp), intent(inout) :: H0_off(maxnh,4)
529 531
530 real(dp), intent(out) :: Enl, Enl_offsiteSO 532 real(dp), intent(out) :: Enl, Enl_offsiteSO
531 logical, intent(in) :: matrix_elements_only533 logical, intent(in) :: matrix_elements_only
@@ -853,10 +855,10 @@
853 ind = listhptr(iio)+j855 ind = listhptr(iio)+j
854 jo = listh(ind)856 jo = listh(ind)
855 H0(ind) = H0(ind) + Vi(jo)857 H0(ind) = H0(ind) + Vi(jo)
856 H0_offsiteSO(ind,1) = H0_offsiteSO(ind,1) + V_so(1,1,jo)858 H0_off(ind,1) = H0_off(ind,1) + V_so(1,1,jo)
857 H0_offsiteSO(ind,2) = H0_offsiteSO(ind,2) + V_so(2,2,jo)859 H0_off(ind,2) = H0_off(ind,2) + V_so(2,2,jo)
858 H0_offsiteSO(ind,3) = H0_offsiteSO(ind,3) + V_so(1,2,jo)860 H0_off(ind,3) = H0_off(ind,3) + V_so(1,2,jo)
859 H0_offsiteSO(ind,4) = H0_offsiteSO(ind,4) + V_so(2,1,jo)861 H0_off(ind,4) = H0_off(ind,4) + V_so(2,1,jo)
860862
861863
862C Careful with this Vi()864C Careful with this Vi()
@@ -894,7 +896,7 @@
894896
895 call timer( 'nlefsm', 2 )897 call timer( 'nlefsm', 2 )
896898
897 end subroutine nlefsm_offsiteSO899 end subroutine nlefsm_SO_off
898900
899c-----------------------------------------------------------------------901c-----------------------------------------------------------------------
900c902c
901903
=== modified file 'Src/setup_H0.F'
--- Src/setup_H0.F 2018-04-19 12:47:26 +0000
+++ Src/setup_H0.F 2018-04-23 11:10:59 +0000
@@ -17,11 +17,10 @@
17 17
18 USE siesta_options, only: g2cut18 USE siesta_options, only: g2cut
19 use sparse_matrices, only: H_kin_1D, H_vkb_1D19 use sparse_matrices, only: H_kin_1D, H_vkb_1D
20 use sparse_matrices, only: H_so_2D20 use sparse_matrices, only: H_so_on_2D, H_so_off_2D
21 use sparse_matrices, only: Dscf21 use sparse_matrices, only: Dscf
2222
23 use sparse_matrices, only: H0_offsiteSO23 use m_nlefsm, only: nlefsm_SO_off
24 use m_nlefsm, only: nlefsm_offsiteSO
25 use m_spin, only: spin24 use m_spin, only: spin
2625
27 use sparse_matrices, only: listh, listhptr, numh, maxnh26 use sparse_matrices, only: listh, listhptr, numh, maxnh
@@ -44,6 +43,7 @@
44 use alloc, only: re_alloc, de_alloc43 use alloc, only: re_alloc, de_alloc
45 use class_dSpData1D, only: val44 use class_dSpData1D, only: val
46 use class_dSpData2D, only: val45 use class_dSpData2D, only: val
46 use class_zSpData2D, only: val
4747
48#ifdef MPI48#ifdef MPI
49 use m_mpi_utils, only: globalize_sum49 use m_mpi_utils, only: globalize_sum
@@ -57,13 +57,15 @@
57 integer :: ia, is57 integer :: ia, is
5858
59 real(dp) :: dummy_Eso59 real(dp) :: dummy_Eso
60 integer :: io, ispin, i, j 60 integer :: ispin, i, j
61 complex(dp) :: Hc, Dc61 complex(dp) :: Dc
62#ifdef MPI62#ifdef MPI
63 real(dp) :: buffer163 real(dp) :: buffer1
64#endif64#endif
6565
66 real(dp), pointer :: H_val(:), H_so(:,:)66 real(dp), pointer :: H_val(:), H_so_on(:,:)
67 complex(dp), pointer :: H_so_off(:,:)
68
6769
68#ifdef DEBUG70#ifdef DEBUG
69 call write_debug( ' PRE setup_H0' )71 call write_debug( ' PRE setup_H0' )
@@ -136,13 +138,14 @@
136 & H_val,138 & H_val,
137 & matrix_elements_only=.true.) 139 & matrix_elements_only=.true.)
138 else140 else
139 H0_offsiteSO=dcmplx(0.0d0,0.0d0)141 H_so_off => val(H_so_off_2D)
140 call nlefsm_offsiteSO(scell, na_u, na_s, isa, xa, indxua,142 H_so_off = dcmplx(0._dp, 0._dp)
143 call nlefsm_SO_off(scell, na_u, na_s, isa, xa, indxua,
141 & maxnh, maxnh, lasto, lastkb, iphorb, iphKB,144 & maxnh, maxnh, lasto, lastkb, iphorb, iphKB,
142 & numh, listhptr, listh, numh, listhptr, listh,145 & numh, listhptr, listh, numh, listhptr, listh,
143 & spin%Grid,146 & spin%Grid,
144 & dummy_E, dummy_Eso, dummy_fa,147 & dummy_E, dummy_Eso, dummy_fa,
145 & dummy_stress, H_val,148 & dummy_stress, H_val, H_so_off,
146 & matrix_elements_only=.true.)149 & matrix_elements_only=.true.)
147150
148151
@@ -153,24 +156,20 @@
153! DM in a such way that the result gives Re{Tr[H_SO*DM]}.156! DM in a such way that the result gives Re{Tr[H_SO*DM]}.
154!157!
155158
156 do io = 1, maxnh159 do i = 1, maxnh
157160
158!-------- Eso(u,u)161!-------- Eso(u,u)
159 Dc = cmplx(Dscf(io,1),Dscf(io,5),kind=dp)162 Dc = cmplx(Dscf(i,1),Dscf(i,5), dp)
160 Hc = H0_offsiteSO(io,1)163 Eso = Eso + real( H_so_off(i,1)*Dc, dp)
161 Eso = Eso + real( Hc*Dc )
162!-------- Eso(d,d)164!-------- Eso(d,d)
163 Dc = cmplx(Dscf(io,2),Dscf(io,6),kind=dp)165 Dc = cmplx(Dscf(i,2),Dscf(i,6),dp)
164 Hc = H0_offsiteSO(io,2)166 Eso = Eso + real( H_so_off(i,2)*Dc, dp)
165 Eso = Eso + real( Hc*Dc )
166!-------- Eso(u,d)167!-------- Eso(u,d)
167 Dc = cmplx(Dscf(io,3),Dscf(io,4),kind=dp)168 Dc = cmplx(Dscf(i,3),Dscf(i,4), dp)
168 Hc = H0_offsiteSO(io,4)169 Eso = Eso + real( H_so_off(i,4)*Dc, dp)
169 Eso = Eso + real( Hc*Dc )
170!-------- Eso(d,u)170!-------- Eso(d,u)
171 Dc = cmplx(Dscf(io,7),-Dscf(io,8),kind=dp)171 Dc = cmplx(Dscf(i,7),-Dscf(i,8), dp)
172 Hc = H0_offsiteSO(io,3)172 Eso = Eso + real( H_so_off(i,3)*Dc, dp)
173 Eso = Eso + real( Hc*Dc )
174173
175 enddo174 enddo
176175
@@ -188,12 +187,12 @@
188! should be enough187! should be enough
189!188!
190 if ( spin%SO_onsite ) then189 if ( spin%SO_onsite ) then
191 H_so => val(H_so_2D)190 H_so_on => val(H_so_on_2D)
192!$OMP parallel workshare default(shared)191!$OMP parallel workshare default(shared)
193 H_so = 0._dp192 H_so_on(:,:) = 0._dp
194!$OMP end parallel workshare193!$OMP end parallel workshare
195 call spinorb(no_u,no_l,iaorb,iphorb,isa,indxuo,194 call spinorb(no_u,no_l,iaorb,iphorb,isa,indxuo,
196 & maxnh,numh,listhptr,listh,Dscf,H_so,Eso)195 & maxnh,numh,listhptr,listh,Dscf,H_so_on,Eso)
197 end if196 end if
198197
199C This will take care of possible changes to the mesh and atomic-related198C This will take care of possible changes to the mesh and atomic-related
200199
=== modified file 'Src/setup_hamiltonian.F'
--- Src/setup_hamiltonian.F 2018-04-20 09:41:55 +0000
+++ Src/setup_hamiltonian.F 2018-04-23 11:10:59 +0000
@@ -14,12 +14,14 @@
1414
15 USE siesta_options15 USE siesta_options
16 use sparse_matrices, only: H_kin_1D, H_vkb_1D16 use sparse_matrices, only: H_kin_1D, H_vkb_1D
17 use sparse_matrices, only: H_ldau_2D, H_so_2D17 use sparse_matrices, only: H_ldau_2D
18 use sparse_matrices, only: H_so_on_2D, H_so_off_2D
18 use sparse_matrices, only: listh, listhptr, numh, maxnh19 use sparse_matrices, only: listh, listhptr, numh, maxnh
19 use sparse_matrices, only: H, S, Hold20 use sparse_matrices, only: H, S, Hold
20 use sparse_matrices, only: Dscf, Escf, xijo21 use sparse_matrices, only: Dscf, Escf, xijo
21 use class_dSpData1D, only: val22 use class_dSpData1D, only: val
22 use class_dSpData2D, only: val23 use class_dSpData2D, only: val
24 use class_zSpData2D, only: val
2325
24 use siesta_geom26 use siesta_geom
25 use atmfuncs, only: uion27 use atmfuncs, only: uion
@@ -42,7 +44,6 @@
42 use m_ntm44 use m_ntm
4345
44 use m_spin, only: spin 46 use m_spin, only: spin
45 use sparse_matrices, only: H0_offsiteSO
4647
47 use m_dipol48 use m_dipol
48 use alloc, only: re_alloc, de_alloc49 use alloc, only: re_alloc, de_alloc
@@ -71,9 +72,11 @@
71 type(filesOut_t) :: filesOut ! blank output file names72 type(filesOut_t) :: filesOut ! blank output file names
72 logical :: use_rhog_in73 logical :: use_rhog_in
7374
74 real(dp), pointer :: H_vkb(:), H_kin(:), H_ldau(:,:), H_so(:,:)75 real(dp), pointer :: H_vkb(:), H_kin(:), H_ldau(:,:)
76 real(dp), pointer :: H_so_on(:,:)
77 complex(dp), pointer:: H_so_off(:,:)
7578
76 complex(dp):: Hc, Dc79 complex(dp):: Dc
77 integer :: ind, i, j80 integer :: ind, i, j
7881
79!------------------------------------------------------------------------- BEGIN82!------------------------------------------------------------------------- BEGIN
@@ -90,53 +93,53 @@
90 do ispin = 1, spin%H93 do ispin = 1, spin%H
91 do io = 1,maxnh94 do io = 1,maxnh
92 Hold(io,ispin) = H(io,ispin)95 Hold(io,ispin) = H(io,ispin)
93 enddo96 end do
94 enddo97 end do
95!$OMP end do98!$OMP end do
9699
97!$OMP single100!$OMP single
98 H_kin => val(H_kin_1D)101 H_kin => val(H_kin_1D)
99 H_vkb => val(H_vkb_1D)102 H_vkb => val(H_vkb_1D)
103
100 if ( spin%SO_onsite ) then104 if ( spin%SO_onsite ) then
101 ! Sadly some compilers (g95), does105 ! Sadly some compilers (g95), does
102 ! not allow bounds for pointer assignments :(106 ! not allow bounds for pointer assignments :(
103 H_so => val(H_so_2D)107 H_so_on => val(H_so_on_2D)
108
109 else if ( spin%SO_offsite ) then
110 H_so_off => val(H_so_off_2D)
111
104 end if112 end if
105!$OMP end single ! keep wait113!$OMP end single ! keep wait
106114
107 if ( .not. spin%SO_offsite ) then115 ! Initialize diagonal Hamiltonian
108 do ispin = 1, spin%spinor116 do ispin = 1, spin%spinor
109 if (ispin .le. 2) then117!$OMP do
110!$OMP do118 do io = 1,maxnh
111 do io = 1,maxnh119 H(io,ispin) = H_kin(io) + H_vkb(io)
112 H(io,ispin) = H_kin(io) + H_vkb(io)120 end do
113 end do121!$OMP end do nowait
114!$OMP end do nowait122 end do
115 else
116!$OMP do
117 do io = 1,maxnh
118 H(io,ispin) = 0.0_dp
119 end do
120!$OMP end do nowait
121 end if
122 end do
123 end if
124123
125 if ( spin%SO_onsite ) then124 if ( spin%SO_onsite ) then
126!$OMP do collapse(2)125!$OMP do collapse(2)
127 do ispin = 3 , spin%H126 do ispin = 3 , spin%H
128 do io = 1,maxnh127 do io = 1,maxnh
129 H(io,ispin) = H_so(io,ispin-2)128 H(io,ispin) = H_so_on(io,ispin-2)
130 end do129 end do
131 end do130 end do
132!$OMP end do nowait131!$OMP end do nowait
133 else if ( spin%NCol ) then132
134!$OMP do133 else
135 do io = 1,maxnh134
136 H(io,3) = 0._dp135!$OMP do collapse(2)
137 H(io,4) = 0._dp136 do ispin = 3 , spin%H
138 end do137 do io = 1,maxnh
139!$OMP end do nowait138 H(io,ispin) = 0._dp
139 end do
140 end do
141!$OMP end do nowait
142
140 end if143 end if
141 144
142! ..................145! ..................
@@ -171,40 +174,36 @@
171! DM in a such way that the result gives Re{Tr[H_SO*DM]}.174! DM in a such way that the result gives Re{Tr[H_SO*DM]}.
172!175!
173176
174 if( spin%SO_offsite ) then177 if ( spin%SO_offsite ) then
175
176 do io = 1, maxnh178 do io = 1, maxnh
177179
178!-------- Eso(u,u)180!-------- Eso(u,u)
179 Dc = cmplx(Dscf(io,1),Dscf(io,5),kind=dp)181 Dc = cmplx(Dscf(io,1),Dscf(io,5), dp)
180 Hc = H0_offsiteSO(io,1)182 Eso = Eso + real( H_so_off(io,1)*Dc, dp)
181 Eso = Eso + real( Hc*Dc )
182!-------- Eso(d,d)183!-------- Eso(d,d)
183 Dc = cmplx(Dscf(io,2),Dscf(io,6),kind=dp)184 Dc = cmplx(Dscf(io,2),Dscf(io,6), dp)
184 Hc = H0_offsiteSO(io,2)185 Eso = Eso + real( H_so_off(io,2)*Dc, dp)
185 Eso = Eso + real( Hc*Dc )
186!-------- Eso(u,d)186!-------- Eso(u,d)
187 Dc = cmplx(Dscf(io,3),Dscf(io,4),kind=dp)187 Dc = cmplx(Dscf(io,3),Dscf(io,4), dp)
188 Hc = H0_offsiteSO(io,4)188 Eso = Eso + real( H_so_off(io,4)*Dc, dp)
189 Eso = Eso + real( Hc*Dc )
190!-------- Eso(d,u)189!-------- Eso(d,u)
191 Dc = cmplx(Dscf(io,7),-Dscf(io,8),kind=dp)190 Dc = cmplx(Dscf(io,7),-Dscf(io,8), dp)
192 Hc = H0_offsiteSO(io,3)191 Eso = Eso + real( H_so_off(io,3)*Dc, dp)
193 Eso = Eso + real( Hc*Dc )192
194193 end do
195 enddo194
196195 else if ( spin%SO_onsite ) then
197 endif196
198
199 if ( spin%SO_onsite ) then
200!$OMP do reduction(+:Eso)197!$OMP do reduction(+:Eso)
201 do io = 1, maxnh198 do io = 1, maxnh
202 Eso = Eso + H_so(io,1)*Dscf(io,7) + H_so(io,2)*Dscf(io,8)199 Eso = Eso + H_so_on(io,1)*Dscf(io,7) +
203 . + H_so(io,5)*Dscf(io,3) + H_so(io,6)*Dscf(io,4)200 & H_so_on(io,2)*Dscf(io,8)+ H_so_on(io,5)*Dscf(io,3) +
204 . - H_so(io,3)*Dscf(io,5) - H_so(io,4)*Dscf(io,6)201 & H_so_on(io,6)*Dscf(io,4)- H_so_on(io,3)*Dscf(io,5) -
202 & H_so_on(io,4)*Dscf(io,6)
205 end do203 end do
206!$OMP end do nowait204!$OMP end do nowait
207 end if205
206 end if
208207
209!$OMP end parallel208!$OMP end parallel
210209
@@ -224,10 +223,6 @@
224! Non-SCF part of total energy223! Non-SCF part of total energy
225 call update_E0()224 call update_E0()
226225
227 if ( spin%SO_offsite ) then
228 H(:,:) = 0.0_dp
229 endif
230
231! Hubbard term for LDA+U: energy, forces, stress and matrix elements ....226! Hubbard term for LDA+U: energy, forces, stress and matrix elements ....
232 if( switch_ldau ) then227 if( switch_ldau ) then
233 if ( spin%NCol ) then228 if ( spin%NCol ) then
@@ -285,24 +280,20 @@
285280
286 if ( spin%SO_offsite ) then281 if ( spin%SO_offsite ) then
287282
288 do is = 1 , spin%spinor283! H(:, [5, 6]) are not updated in dhscf, see vmat for details.
289 do io = 1, maxnh284
290 H(io,is) = H(io,is) + H_kin(io) + H_vkb(io)285!------- H(u,u)
291 enddo286 H(:,1) = H(:,1) + real(H_so_off(:,1), dp)
292 enddo287 H(:,5) = dimag(H_so_off(:,1))
293288!------- H(d,d)
294C------- H(u,u)289 H(:,2) = H(:,2) + real(H_so_off(:,2), dp)
295 H(:,1) = H(:,1) + real(H0_offsiteSO(:,1))290 H(:,6) = dimag(H_so_off(:,2))
296 H(:,5) = imag(H0_offsiteSO(:,1))291!------- H(u,d)
297C------- H(d,d)292 H(:,3) = H(:,3) + real(H_so_off(:,3), dp)
298 H(:,2) = H(:,2) + real(H0_offsiteSO(:,2))293 H(:,4) = H(:,4) +dimag(H_so_off(:,3))
299 H(:,6) = imag(H0_offsiteSO(:,2))294!------- H(d,u)
300C------- H(u,d)295 H(:,7) = H(:,7) + real(H_so_off(:,4), dp)
301 H(:,3) = H(:,3) + real(H0_offsiteSO(:,3))296 H(:,8) = H(:,8) -dimag(H_so_off(:,4))
302 H(:,4) = H(:,4) + imag(H0_offsiteSO(:,3))
303C------- H(d,u)
304 H(:,7) = H(:,7) + real(H0_offsiteSO(:,4))
305 H(:,8) = H(:,8) - imag(H0_offsiteSO(:,4))
306297
307 endif298 endif
308299
309300
=== modified file 'Src/sparse_matrices.F'
--- Src/sparse_matrices.F 2018-04-19 12:47:26 +0000
+++ Src/sparse_matrices.F 2018-04-23 11:10:59 +0000
@@ -9,6 +9,7 @@
9 use precision9 use precision
10 use class_dSpData1D10 use class_dSpData1D
11 use class_dSpData2D11 use class_dSpData2D
12 use class_zSpData2D
12 use class_Sparsity13 use class_Sparsity
13 use class_OrbitalDistribution14 use class_OrbitalDistribution
14 use class_Fstack_Pair_Geometry_dSpData2D15 use class_Fstack_Pair_Geometry_dSpData2D
@@ -37,7 +38,10 @@
37 ! Formerly there was a single array H0 for this38 ! Formerly there was a single array H0 for this
38 type(dSpData1D), public, save :: H_vkb_1D, H_kin_1D39 type(dSpData1D), public, save :: H_vkb_1D, H_kin_1D
39 ! LDA+U and spin-orbit coupling Hamiltonian40 ! LDA+U and spin-orbit coupling Hamiltonian
40 type(dSpData2D), public, save :: H_ldau_2D, H_so_2D41 type(dSpData2D), public, save :: H_ldau_2D, H_so_on_2D
42 ! Spin-orbit off-site
43 type(zSpData2D), public, save :: H_so_off_2D
44
4145
42 ! New interface data46 ! New interface data
43 type(Sparsity), public, save :: sparse_pattern47 type(Sparsity), public, save :: sparse_pattern
@@ -56,7 +60,6 @@
5660
57 subroutine resetSparseMatrices( )61 subroutine resetSparseMatrices( )
58 use alloc, only : de_alloc62 use alloc, only : de_alloc
59 use m_spin, only : spin
6063
61 implicit none64 implicit none
6265
@@ -67,7 +70,8 @@
67 call delete( H_kin_1D )70 call delete( H_kin_1D )
68 call delete( H_vkb_1D )71 call delete( H_vkb_1D )
69 call delete( H_ldau_2D )72 call delete( H_ldau_2D )
70 call delete( H_so_2D )73 call delete( H_so_on_2D )
74 call delete( H_so_off_2D )
7175
72 call delete( DM_2D ) ; nullify(Dscf)76 call delete( DM_2D ) ; nullify(Dscf)
73 call delete( EDM_2D ) ; nullify(Escf)77 call delete( EDM_2D ) ; nullify(Escf)
@@ -76,10 +80,6 @@
76 call delete( xij_2D ) ; nullify(xijo)80 call delete( xij_2D ) ; nullify(xijo)
7781
78 call delete( DM_history )82 call delete( DM_history )
79
80 if( spin%SO_offsite ) then
81 call de_alloc( H0_offsiteSO, 'H0_offsiteSO', 'sparseMat' )
82 endif
83 83
84 ! Using MixH is bad as several utilities84 ! Using MixH is bad as several utilities
85 ! are not relying on FoX, better to leave out85 ! are not relying on FoX, better to leave out
8686
=== modified file 'Src/state_init.F'
--- Src/state_init.F 2018-04-19 12:47:26 +0000
+++ Src/state_init.F 2018-04-23 11:10:59 +0000
@@ -26,10 +26,9 @@
26 use sparse_matrices, only: xijo, xij_2D26 use sparse_matrices, only: xijo, xij_2D
27 use sparse_matrices, only: S , S_1D27 use sparse_matrices, only: S , S_1D
2828
29 use sparse_matrices, only: H0_offsiteSO
30
31 use sparse_matrices, only: H_kin_1D, H_vkb_1D29 use sparse_matrices, only: H_kin_1D, H_vkb_1D
32 use sparse_matrices, only: H_ldau_2D, H_so_2D30 use sparse_matrices, only: H_ldau_2D
31 use sparse_matrices, only: H_so_on_2D, H_so_off_2D
3332
34 use sparse_matrices, only: sparse_pattern33 use sparse_matrices, only: sparse_pattern
35 use sparse_matrices, only: block_dist, single_dist34 use sparse_matrices, only: block_dist, single_dist
@@ -100,6 +99,7 @@
100 use class_Sparsity99 use class_Sparsity
101 use class_dSpData1D100 use class_dSpData1D
102 use class_dSpData2D101 use class_dSpData2D
102 use class_zSpData2D
103 use class_dData2D103 use class_dData2D
104#ifdef TEST_IO104#ifdef TEST_IO
105 use m_test_io105 use m_test_io
@@ -564,15 +564,14 @@
564 end if564 end if
565 end if565 end if
566 566
567 if ( spin%SO_onsite ) then567 if ( spin%SO_onsite ) then
568 write(oname,"(a,i0)") "H_so (onsite) at geom step ", istep568 write(oname,"(a,i0)") "H_so (onsite) at geom step ", istep
569 call newdSpData2D(sparse_pattern,spin%H - 2,569 call newdSpData2D(sparse_pattern,spin%H - 2,
570 & block_dist,H_so_2D,name=oname)570 & block_dist,H_so_on_2D,name=oname)
571 end if571 else if ( spin%SO_offsite ) then
572572 write(oname,"(a,i0)") "H_so (offsite) at geom step ", istep
573 if ( spin%SO_offsite ) then573 call newzSpData2D(sparse_pattern,4,
574 call re_alloc(H0_offsiteSO, 1,maxnh, 1,4, 'H0_offsiteSO',574 & block_dist,H_so_off_2D,name=oname)
575 $ 'state_init')
576 endif575 endif
577576
578 write(oname,"(a,i0)") "S at geom step ", istep577 write(oname,"(a,i0)") "S at geom step ", istep
579578
=== modified file 'version.info'
--- version.info 2018-04-23 09:26:04 +0000
+++ version.info 2018-04-23 11:10:59 +0000
@@ -1,2 +1,2 @@
1trunk-688--merge-OSSO-7001trunk-688--merge-OSSO-700--class-1
22

Subscribers

People subscribed via source and target branches

to all changes: