~jose-soler/siesta/unfolding

« back to all changes in this revision

Viewing changes to Src/final_H_f_stress.F

  • Committer: Jose M Soler
  • Date: 2018-10-29 18:09:36 UTC
  • mfrom: (718.1.26 trunk)
  • Revision ID: jose.soler@uam.es-20181029180936-7p8z2sbde0mkvxoz
Merged with trunk-744

modified:
  Docs/Contributors.txt
  Docs/siesta.bib
  Docs/siesta.tex
  Docs/tbtrans.tex
  Docs/tex/fdf.tex
  Docs/tex/setup.tex
  RELEASE_NOTES
  Src/final_H_f_stress.F
  Src/grdsam.F
  Src/m_ts_options.F90
  Src/m_ts_tri_init.F90
  Src/m_ts_weight.F90
  Src/save_density_matrix.F90
  Src/state_init.F
  Util/TS/TBtrans/Makefile
  Util/TS/TBtrans/m_tbt_tri_init.F90

Show diffs side-by-side

added added

removed removed

Lines of Context:
17
17
#ifdef NCDF_4
18
18
      use siesta_options, only: write_cdf
19
19
#endif
 
20
      use siesta_options, only: fixspin
20
21
      use siesta_options, only: g2cut, savehs, temp, nmove
21
22
      use siesta_options, only: recompute_H_after_scf
22
23
      use sparse_matrices, only: numh, listh, listhptr
23
24
      use sparse_matrices, only: H, S, Dscf, Escf, maxnh, xijo
24
25
      use sparse_matrices, only: H_ldau_2D, H_so_off_2D
 
26
      use class_dSpData1D, only: val
25
27
      use class_dSpData2D, only: val
26
28
      use class_zSpData2D, only: val
27
29
 
86
88
      integer               :: istr    ! Calculate stress? 0=>no, 1=>yes
87
89
      real(dp)              :: g2max
88
90
      ! To avoid overwriting the current Hamiltonian and S
89
 
      real(dp), pointer     :: H_tmp(:,:) => null()
90
 
      real(dp), pointer     :: H_ldau(:,:)
 
91
      real(dp), pointer :: H_tmp(:,:) => null()
 
92
      real(dp), pointer :: H_ldau(:,:), S_tmp(:)
91
93
      complex(dp), pointer  :: H_so_off(:,:)
92
94
 
93
95
#ifdef FINAL_CHECK_HS
339
341
      end if
340
342
#endif
341
343
#endif
 
344
 
 
345
      if ( fixspin .and. (write_tshs_history .or. TS_HS_save) ) then
 
346
        ! For fixed spin calculations we shift the Hamiltonian according
 
347
        ! to the first spin such that we have a "common" E_F, after writing
 
348
        ! we shift back again.
 
349
        H_tmp => val(H_2D)
 
350
        S_tmp => val(S_1D)
 
351
        call daxpy(size(S_tmp),Efs(1)-Efs(2),S_tmp(1),1,H_tmp(1,2),1)
 
352
        Ef = Efs(1)
 
353
      end if
 
354
      
342
355
      if ( write_tshs_history ) then
343
356
         ! This is "pure" MD and we only write consecutive numbers
344
357
         ! Together with this you cannot also save FC
382
395
         end if
383
396
      endif
384
397
 
 
398
      if ( fixspin .and. (write_tshs_history .or. TS_HS_save) ) then
 
399
        ! Shift back
 
400
        call daxpy(size(S_tmp),Efs(2)-Efs(1),S_tmp(1),1,H_tmp(1,2),1)
 
401
      end if
 
402
 
385
403
!----------------------------------------------------------------------- END
386
404
      END subroutine final_H_f_stress
387
405
      END module m_final_H_f_stress