[gpaw-users] GPAW :plane wave mode and optimization use out of memory

Marcin Dulak Marcin.Dulak at fysik.dtu.dk
Fri Dec 14 19:57:05 CET 2012


Hi,

let's continue on the mailing list.

On 12/12/12 13:32, qian li wrote:
> Dear Marcin,
>
> Thank you very much for your help!
> There two questions:
>
> 1.
>
> It seems WFS don't converge in the optimization of benzene.So I use 
> plane wave mode as you suggest me.
>
> I use mode= 'pw' in the input file, but it errors. Do you know how to 
> set plane wave mode in the GPAW?
>
> This is my input file about pw
>
> fromase import*
>
> fromgpaw import*
>
> from gpaw.transport.tools import sort_atoms
>
> from gpaw.mixer import Mixer
>
> from gpaw import GPAW, PW
>
> from ase.io <http://ase.io> import read,write
>
> from ase.visualize import view
>
> from ase.lattice.surface import *
>
> from ase.optimize import QuasiNewton
>
> from ase.data.molecules import molecule
>
> import numpy as np
>
> atoms = read('optimize_gas.traj')
>
> atoms.set_cell([20,20,20])
>
> atoms.center()
>
> vdw = 'vdW-DF'
>
> calc = GPAW(xc=vdw,# XC functional
>
>   mode='pw',
>
this should be rather:
mode=PW(500)

>   nbands=32,
>
>   mixer=Mixer(0.05, 2),
>
>   maxiter=250, #eigensolver='cg',
>
>   convergence={'energy': 0.00001},
>
>             txt="optimize_gas_DF.txt") # Output text file
>
> atoms.set_calculator(calc) # Attach calculator to atoms
>
> relax = QuasiNewton(atoms,
>
>           trajectory='optimize_gas_DF.traj',
>
>           logfile='qn.log')
>
> relax.run(fmax=0.01)
>
> calc.write("optimize_gas_DF.gpw")
>
> This is the error:
>
> Traceback (most recent call last):
>   File "benze.py", line 28, in ?
>     relax.run(fmax=0.01)
>   File "/kemi/strange/ase/ase/optimize/optimize.py", line 114, in run
>     f = self.atoms.get_forces()
>   File "/kemi/strange/ase/ase/atoms.py", line 650, in get_forces
>     forces = self._calc.get_forces(self)
>   File "/users/kemi/strange/gpaw/gpaw/aseinterface.py", line 69, in 
> get_forces
> force_call_to_set_positions=force_call_to_set_positions)
>   File "/users/kemi/strange/gpaw/gpaw/paw.py", line 230, in calculate
>     self.set_positions(atoms)
>   File "/users/kemi/strange/gpaw/gpaw/paw.py", line 308, in set_positions
>     self.wfs.initialize(self.density, self.hamiltonian, spos_ac)
>   File "/kemi/strange/gpaw/gpaw/wavefunctions/fdpw.py", line 67, in 
> initialize
>     self.initialize_wave_functions_from_basis_functions(
>   File "/kemi/strange/gpaw/gpaw/wavefunctions/fdpw.py", line 75, in 
> initialize_wave_functions_from_basis_functions
>     lcaoksl, lcaobd = self.initksl, self.initksl.bd 
> <http://self.initksl.bd>
> AttributeError: 'NoneType' object has no attribute 'bd'
> GPAW CLEANUP (node 4): exceptions.AttributeError occurred.  Calling 
> MPI_Abort!
use basis='dzp' in order to have enough initial bands (see 
https://listserv.fysik.dtu.dk/pipermail/gpaw-users/2012-December/001898.html) 

>
> 2.
>
> In addition, when I changed h from 0.20 to 0.12, they needed a lot 
> memory and the machines were down.
> This is the input file:
in general, use the --dry-run=N, when N is number of cores to estimate 
the "Calculator" memory required per core.
See 
https://wiki.fysik.dtu.dk/gpaw/documentation/parallel_runs/parallel_runs.html
Note however that the vdw implementation does not have memory estimate 
implemented,
resulting in underestimation by factor of 10 in the case of this 
particular calculation.
There is some information about reducing memory usage in vdw calculations:
https://listserv.fysik.dtu.dk/pipermail/gpaw-users/2011-April/000798.html

I attach an example that converged SCF for me, but it's true that vdw 
run for this case shows instabilities, e.g.:
iter:  96  17:45:11  -3.8   -3.5     -103.563981  0      2
iter:  97  17:45:47  -2.7   -3.5     -103.563867  0      6
iter:  98  17:46:23  +1.5   -4.1     -104.711229  0      2
iter:  99  17:46:59  +0.7   -0.9     -105.401991  0      2
iter: 100  17:47:50  -1.9   -0.9     -103.549780  0      71

Best regards,

Marcin
>
> from ase import *
>
> from gpaw import *
>
> from gpaw.transport.tools import sort_atoms
>
> from gpaw.mixer import Mixer
>
> from ase.io <http://ase.io/> import read,write
>
> from ase.visualize import view
>
> from ase.lattice.surface import *
>
> from ase.optimize import QuasiNewton
>
> from ase.data.molecules import molecule
>
> import numpy as np
>
> atoms = read('benzene.xyz')
>
> atoms.set_cell([20,20,20])
>
> atoms.center()
>
> vdw = 'vdW-DF'
>
> calc = GPAW(xc=vdw,# XC functional
>
>   h=0.12,
>
>     basis='dzp',
>
>   mixer=Mixer(0.05, 2),
>
>   maxiter=250, eigensolver='cg',
>
>             txt="optimize_gas_DF.txt") # Output text file
>
> atoms.set_calculator(calc) # Attach calculator to atoms
>
> relax = QuasiNewton(atoms,
>
>           trajectory='optimize_gas_DF.traj',
>
>             logfile='qn.log')
>
> relax.run(fmax=0.01)
>
> calc.write("optimize_gas_DF.gpw")
>
>
> This is the .out for memory:
>
>
> Mem:         96680      64009      32671        0         53       1532
>
> Mem:         96680      15645      81034        0         53       1532
>
> Mem:         96680      39696      56984          0         53       1532
>
> Mem:         96680      50398      46281        0         53       1532
>
> Mem:         96680      56323      40356        0         53       1532
>
> Mem:         96680      62551      34129        0         53       1533
>
> Mem:         96680      73105      23575        0         53       1533
>
> Mem:         96680      67812      28867          0         53       1533
>
> Mem:         96680      70900      25780        0         54       1533
>
> Mem:         96680      79990      16690        0         54       1533
>
> Mem:         96680      76051      20629          0         54       1533
>
> Mem:         96680      79141      17538        0         54       1533
>
> Mem:         96680      72374      24306        0         54       1533
>
> Mem:         96680      88612       8068          0         40       1150
>
> Mem:         96680      96321        359        0         40       1141
>
>
> Do you know why this happens?
> Could you give me some suggestions?
> Thank you very much!
>
> Best regards


-- 
***********************************
  
Marcin Dulak
Technical University of Denmark
Department of Physics
Building 307, Room 229
DK-2800 Kongens Lyngby
Denmark
Tel.: (+45) 4525 3157
Fax.: (+45) 4593 2399
email:Marcin.Dulak at fysik.dtu.dk

***********************************

-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://listserv.fysik.dtu.dk/pipermail/gpaw-users/attachments/20121214/34081809/attachment-0001.html 
-------------- next part --------------
A non-text attachment was scrubbed...
Name: test_default.py
Type: text/x-python
Size: 989 bytes
Desc: not available
Url : http://listserv.fysik.dtu.dk/pipermail/gpaw-users/attachments/20121214/34081809/attachment-0001.py 


More information about the gpaw-users mailing list