[gpaw-users] PDOS with LCAO
Ask Hjorth Larsen
asklarsen at gmail.com
Mon Jan 14 13:11:57 CET 2013
Hi
2012/11/19 Ask Hjorth Larsen <asklarsen at gmail.com>:
> Hi
>
> I will finish it [pdos based on atomic orbitals in LCAO mode] and check it in sooner or later.
>
> Regards
> Ask
The above has been done and I claim that it works, but have only done
very little testing, so there may be cases where it does not.
Basically:
from gpaw.utilities.dos import LCAODOS, RestartLCAODOS, fold
if you're restarting:
calc = GPAW(filename, ...)
dos = RestartLCAODOS(calc)
else:
calc = GPAW(...)
<perform calculation>
dos = LCAODOS(calc)
energies, weights = dos.get_subspace_pdos([7, 8, 9, 12])
# That gets the subspace PDOS on the space of orbitals number 7-9 and 12
e, w = fold(eps_n * Hartree, weights, 800, 0.1)
e and w are continuous energy/DOS arrays. For example the integral of
e dw should be the number of states. I would appreciate if someone
can test any of this. It should in some ways be far superior to the
PDOS calculated from PAW projectors.
Also:
dos.get_atomic_subspace_pdos(a) gets the dos for all basis functions on atom a
dos.get_orbital_pdos(M) gets the dos on only one orbital
Finaly there's the keyword ravel=True. If set to False you get the
energies and weights as multidimensional lists with kpt/spin/band
indices.
Best regards
Ask
>
> 2012/11/18 Silvio a Beccara <s.abeccara at gmail.com>:
>> Dear Ask,
>>
>> I would be interested in having the PDOS with LCAO, too. I would use it to
>> calculate overlap matrices.
>> Best
>>
>> Silvio
>>
>>
>>>
>>>
>>> That's not implemented I'm afraid.
>>>
>>> The best way to calculate PDOS with LCAO is to use the orbitals
>>> directly. The advantage is that you project onto something complete
>>> *and* take the non-orthogonality into account. This makes it the only
>>> quantitiative PDOS. Its disadvantage is that it hasn't been
>>> implemented either, but I claim to have an almost finished version
>>> lying around somewhere, if some day someone expresses sufficient
>>> interest.
>>>
>>> Regards
>>> Ask
>>>
>>>
>>> ------------------------------
>>>
>>> Message: 2
>>> Date: Thu, 15 Nov 2012 14:31:51 +0100
>>> From: Torsten Hahn <torstenhahn at fastmail.fm>
>>> Subject: Re: [gpaw-users] get_all_electron_ldos in lcao
>>> To: Ask Hjorth Larsen <asklarsen at gmail.com>
>>> Cc: gpaw-users <gpaw-users at listserv.fysik.dtu.dk>
>>> Message-ID: <EA094016-2088-424C-944F-0242047446AF at fastmail.fm>
>>> Content-Type: text/plain; charset=windows-1252
>>>
>>> Am 15.11.2012 um 12:31 schrieb Ask Hjorth Larsen <asklarsen at gmail.com>:
>>>
>>> > Hi
>>> >
>>> > 2012/11/15 Torsten Hahn <torstenhahn at fastmail.fm>:
>>> >> Hi,
>>> >>
>>> >> is there a way for LCAO-mode calculations to use the
>>> >> get_all_electron_ldos function?
>>> >>
>>> >> If you give the function-argument wf_k other than 'None', the functions
>>> >> fails (see. pDOS example on the GPAW-Website) due to the fact that
>>> >> kpt.psit_nG is not set in LCAO mode?
>>> >>
>>> >> Is there a workaround ?
>>> >>
>>> >> Best regards,
>>> >> Torsten.
>>> >> _______________________________________________
>>> >> gpaw-users mailing list
>>> >> gpaw-users at listserv.fysik.dtu.dk
>>> >> https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users
>>> >
>>> > That's not implemented I'm afraid.
>>> >
>>> > The best way to calculate PDOS with LCAO is to use the orbitals
>>> > directly. The advantage is that you project onto something complete
>>> > *and* take the non-orthogonality into account. This makes it the only
>>> > quantitiative PDOS. Its disadvantage is that it hasn't been
>>> > implemented either, but I claim to have an almost finished version
>>> > lying around somewhere, if some day someone expresses sufficient
>>> > interest.
>>> >
>>> > Regards
>>> > Ask
>>>
>>> I express my interest now and i could also spend some time to test /
>>> debugging ?
>>>
>>> Regards,
>>> Torsten.
>>>
>>>
>>>
>>>
>>>
>>> ------------------------------
>>>
>>> Message: 3
>>> Date: Fri, 16 Nov 2012 11:52:01 +0100
>>> From: Torsten Hahn <torstenhahn at fastmail.fm>
>>> Subject: Re: [gpaw-users] get parallel hamiltonian
>>> To: Ask Hjorth Larsen <asklarsen at gmail.com>
>>> Cc: gpaw-users <gpaw-users at listserv.fysik.dtu.dk>
>>> Message-ID: <943B1B81-9B5B-4E5F-852B-D24E205EE90F at fastmail.fm>
>>> Content-Type: text/plain; charset=us-ascii
>>>
>>>
>>> Am 26.10.2012 um 23:38 schrieb Ask Hjorth Larsen <asklarsen at gmail.com>:
>>>
>>> > 2012/10/26 Ask Hjorth Larsen <asklarsen at gmail.com>:
>>> > (...)
>>> >>
>>> >> For each s and k, do:
>>> >>
>>> >> Hfull_MM = calc.wfs.ksl.mmdescriptor.collect_on_master(H_MM)
>>> >>
>>> >> ksl is the "Kohn-Sham layout", representing the parallelization of KS
>>> >> states.
>>> >>
>>> >> mmdescriptor is a matrix descriptor (from gpaw/blacs.py) describing
>>> >> the blocked ("mm" rather than "MM") layout of the Hamiltonian.
>>> >>
>>> >> Regards
>>> >> Ask
>>> >
>>> > (I write this partially from memory so if something goes wrong, some
>>> > tinkering may be required)
>>>
>>>
>>> thx, it worked just fine.
>>>
>>>
>>>
>>>
>>> ------------------------------
>>>
>>> _______________________________________________
>>> gpaw-users mailing list
>>> gpaw-users at listserv.fysik.dtu.dk
>>> https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users
>>>
>>> End of gpaw-users Digest, Vol 34, Issue 12
>>> ******************************************
>>
>>
>>
>> _______________________________________________
>> gpaw-users mailing list
>> gpaw-users at listserv.fysik.dtu.dk
>> https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users
More information about the gpaw-users
mailing list