[DFTB-Plus-User] consistency check for "makeDensityMatrix"

Sharma SRK Chaitanya Yamijala sharmajncasr at gmail.com
Thu Nov 30 18:47:55 CET 2017


Hi Ben,

Thank you. I have checked. It is exactly 30 (as expected) when I have
considered the lower triangle. Thanks a lot for pointing the mistake.

I am getting the sqrHamSize from the code itself (from the below block).
  ! Nr. of independent spin Hamiltonians
  select case (nSpin)
  case (1)
    nSpinHams = 1
    sqrHamSize = nOrb
  case (2)
    nSpinHams = 2
    sqrHamSize = nOrb
  case (4)
    nSpinHams = 1
    sqrHamSize = 2 * nOrb
  end select

Thanks for pointing me towards spin. Currently, I am doing everything for
the restricted case only. But, I will remember this point when I code the
unrestricted parts.

Thanks again,
Sincerely,
Sharma.





--------------------------------------------------------------------
Sharma
http://www.chem.rochester.edu/groups/huo/people/


On Thu, Nov 30, 2017 at 12:26 PM, Ben Hourahine <
benjamin.hourahine at strath.ac.uk> wrote:

> Hello Sharma,
>
> is most cases only the lower triangle of large symmetric matrices are
> stored (not the upper triangle as your loops run). The blockSymmetrizeHS
> only modifies the atomic on-site part of such matrices. Try something like:
> sumFilling = 0.0d0
> do iSurf = 1, sqrHamSize
>     sumFilling = sumFilling + SSqrReal(iSurf,iSurf) *
> currentOverlap(iSurf,iSurf)
>     do jSurf = iSurf+1, sqrHamSize
>              sumFilling = sumFilling + 2.0_dp * SSqrReal(jSurf,iSurf) *
> currentOverlap(jSurf,iSurf)
>     end do
> end do
>
> Where do you get the value of sqrHamSize?
>
> Also be careful about spin as you've suppressed that index.
> Regards
>
> Ben
>
>
> On 30/11/17 16:46, Sharma SRK Chaitanya Yamijala wrote:
>
> Dear DFTB+ developers,
>
> I am using the output of "makeDensityMatrix" subroutine from a nonSCC
> calculation (on Benzene) for testing our lab code on PLDM.
>
> As
>
> \rho_{\mu \nu} = \sum_{i}^{Occ} n_{i} c_{\mu i}c_{\nu i} \\
>
> and
>
> \langle \Phi_{i}| \Phi_{i}\rangle = \sum_{\mu \nu}^{nAOs} c_{\mu i}c_{\nu
> i} S_{\mu \nu} = 1 ~ \forall ~ i\\
>
> I was expecting that sum over (rho(mu,nu) * S(mu,nu)) will give me the
> summation over the fillings of all occupied orbitals (as given by the below
> equation).
> =
> \sum_{\mu \nu}^{nAOs} \rho_{\mu \nu} S_{\mu \nu} = \sum_{\mu \nu}
> \sum_{i}^{Occ} n_{i} c_{\mu i}c_{\nu i} S_{\mu \nu} = \sum_{i}^{Occ} n_{i}
> \left(\sum_{\mu \nu} c_{\mu i}c_{\nu i} S_{\mu \nu} \right) =
> \sum_{i}^{Occ} n_{i}
>
> If the above equation is right, then for the case of benzene and for
> the nonSCC calculation this should give me 30 as the answer (15 occupied
> levels and with an occupation number of 2 for each). With the minimum
> modifications to the DFTB+ code, I was not able to get this answer.
>
> *The minimal modifications I did were below:*
>
> Before the scc loop, I saved the overlap matrix and I have used this and
> the output of "makeDensityMatrix" routine
>
>     if (tCheck) then
>       call unpackHS(overlapReal, over, neighborList%iNeighbor, nNeighbor,&
>         & iAtomStart, iPair, img2CentCell)
>       call blockSymmetrizeHS(overlapReal, iAtomStart)
>       currentOverlap = overlapReal
>     end if
>
> !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
> !!!!!!!!!!!!!!!
>     !! SCC-loop
> !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
> !!!!!!!!!!!!!!!
>
> ....
> ....
>
> if (tDensON2) then
>               call makeDensityMatrix(SSqrReal, HSqrReal(:,:,iSpin2),
> filling(:,1,iSpin),&
>                   & neighborlist%iNeighbor, nNeighbor, orb, iAtomStart,
> img2CentCell)
>             else
>               call makeDensityMatrix(SSqrReal, HSqrReal(:,:,iSpin2),
> filling(:,1,iSpin))
>
> !********************************************************************
>               sumFilling=0.0d0
>               do iSurf=1,sqrHamSize
>                  sumFilling=sumFilling+SSqrReal(iSurf,iSurf)*
> currentOverlap(iSurf,iSurf)
>               end do
>               do iSurf=1,sqrHamSize-1
>                 do jSurf=iSurf+1,sqrHamSize
>                    sumFilling=sumFilling+2.0d0*SSqrReal(iSurf,jSurf)*
> currentOverlap(iSurf,jSurf)
>                 end do
>               end do
>               write(*,*)"sumFilling",sumFilling
> !********************************************************************
>
>             end if
>
>
> I believe I am doing something wrong (either in the code or in the
> equation) and any advice on your end is highly appreciated.
>
> Thanking you,
> Sincerely,
> Sharma.
>
> P.S.: You can use https://www.codecogs.com/latex/eqneditor.php to
> visualize the above equations (by copying and pasting them)
>
>
> --------------------------------------------------------------------
> Sharma
> http://www.chem.rochester.edu/groups/huo/people/
>
>
>
> _______________________________________________
> DFTB-Plus-User mailing listDFTB-Plus-User at mailman.zfn.uni-bremen.dehttps://mailman.zfn.uni-bremen.de/cgi-bin/mailman/listinfo/dftb-plus-user
>
>
> --
>       Dr. B. Hourahine, SUPA, Department of Physics,
>     University of Strathclyde, John Anderson Building,
>             107 Rottenrow, Glasgow G4 0NG, UK.
>     +44 141 548 2325 <+44%20141%20548%202325>, benjamin.hourahine at strath.ac.uk
>
> 2013/4 THE Awards Entrepreneurial University of the Year
>       2012/13 THE Awards UK University of the Year
>
>    The University of Strathclyde is a charitable body,
>         registered in Scotland, number SC015263
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://mailman.zfn.uni-bremen.de/pipermail/dftb-plus-user/attachments/20171130/2c5ce97c/attachment.html>


More information about the DFTB-Plus-User mailing list