2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-03-24 23:25:17 +08:00
|
|
|
! Copyright (C) 2002-2003 PWSCF group
|
2003-01-20 05:58:50 +08:00
|
|
|
! This file is distributed under the terms of the
|
|
|
|
! GNU General Public License. See the file `License'
|
|
|
|
! in the root directory of the present distribution,
|
|
|
|
! or http://www.gnu.org/copyleft/gpl.txt .
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-06-12 15:43:14 +08:00
|
|
|
#include "machine.h"
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
!----------------------------------------------------------------------------
|
|
|
|
SUBROUTINE iosys()
|
|
|
|
!-----------------------------------------------------------------------------
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... this subroutine reads input data from standard input ( unit 5 )
|
|
|
|
! ... ---------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! ... access the modules renaming the variables that have the same name
|
|
|
|
! ... as the input parameters, this is required in order to use a code
|
|
|
|
! ... independent input parser
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
!
|
2003-11-04 18:53:05 +08:00
|
|
|
USE constants, ONLY : AU, eV_to_kelvin
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY : stdout
|
2003-10-17 00:30:12 +08:00
|
|
|
USE bp, ONLY : nppstr_ => nppstr, &
|
|
|
|
gdir_ => gdir, &
|
|
|
|
lberry_ => lberry
|
|
|
|
USE brilz, ONLY : at, alat, omega, &
|
|
|
|
celldm_ => celldm, &
|
|
|
|
ibrav_ => ibrav
|
|
|
|
USE basis, ONLY : nat_ => nat, &
|
|
|
|
ntyp_ => ntyp, &
|
|
|
|
ityp, tau, atomic_positions, atm, &
|
|
|
|
startingwfc_ => startingwfc, &
|
|
|
|
startingpot_ => startingpot, &
|
|
|
|
startingconfig
|
|
|
|
USE char, ONLY : title_ => title, &
|
|
|
|
crystal
|
|
|
|
USE cellmd, ONLY : cmass, ttol, omega_old, at_old, ntcheck, &
|
|
|
|
cell_factor_ => cell_factor , &
|
|
|
|
press_ => press, &
|
|
|
|
calc, lmovecell
|
|
|
|
USE constants, ONLY : pi, rytoev, uakbar, amconv, bohr_radius_angs
|
|
|
|
USE dynam, ONLY : dt_ => dt, &
|
|
|
|
temperature, amass, delta_t, nraise
|
|
|
|
USE extfield, ONLY : tefield_ => tefield, &
|
|
|
|
dipfield_ => dipfield, &
|
|
|
|
edir_ => edir, &
|
|
|
|
emaxpos_ => emaxpos, &
|
|
|
|
eopreg_ => eopreg, &
|
|
|
|
eamp_ => eamp, &
|
|
|
|
forcefield
|
2003-11-10 02:30:08 +08:00
|
|
|
USE io_files, ONLY : input_drho, output_drho
|
2003-10-17 00:30:12 +08:00
|
|
|
USE force_mod, ONLY : lforce, lstres, force
|
|
|
|
USE gvect, ONLY : dual, &
|
|
|
|
nr1_ => nr1, &
|
|
|
|
nr2_ => nr2, &
|
|
|
|
nr3_ => nr3, &
|
|
|
|
ecutwfc_ => ecutwfc, &
|
|
|
|
ecfixed_ => ecfixed, &
|
|
|
|
qcutz_ => qcutz, &
|
|
|
|
q2sigma_ => q2sigma
|
|
|
|
USE gsmooth, ONLY : nr1s_ => nr1s, &
|
|
|
|
nr2s_ => nr2s, &
|
|
|
|
nr3s_ => nr3s
|
|
|
|
USE klist, ONLY : xk, wk, nks, ngauss,&
|
|
|
|
xqq_ => xqq, &
|
|
|
|
degauss_ => degauss, &
|
|
|
|
nelec_ => nelec
|
|
|
|
USE ktetra, ONLY : nk1, nk2, nk3, k1, k2, k3, ltetra
|
|
|
|
USE ldaU, ONLY : Hubbard_U_ => hubbard_u, &
|
|
|
|
Hubbard_alpha_ => hubbard_alpha, &
|
2004-02-14 16:39:34 +08:00
|
|
|
niter_with_fixed_ns, starting_ns, &
|
2003-10-17 00:30:12 +08:00
|
|
|
lda_plus_u_ => lda_plus_u
|
|
|
|
USE lsda_mod, ONLY : nspin_ => nspin, &
|
|
|
|
starting_magnetization_ => starting_magnetization, &
|
|
|
|
lsda
|
|
|
|
USE io_files, ONLY : tmp_dir, &
|
|
|
|
prefix_ => prefix, &
|
|
|
|
pseudo_dir_ => pseudo_dir, &
|
|
|
|
psfile
|
|
|
|
USE relax, ONLY : if_pos, epsf, starting_scf_threshold, &
|
|
|
|
restart_bfgs, epse
|
2004-02-09 19:15:33 +08:00
|
|
|
USE control_flags, ONLY : diis_ethr_cg, diis_ndim, diis_wfc_keep, isolve, &
|
2003-10-17 00:30:12 +08:00
|
|
|
max_cg_iter, diis_buff, david, imix, nmix, &
|
2004-01-23 23:08:03 +08:00
|
|
|
iverbosity, tr2, niter, order, iswitch, &
|
2003-10-17 00:30:12 +08:00
|
|
|
upscale_ => upscale, &
|
|
|
|
mixing_beta_ => mixing_beta, &
|
|
|
|
nstep_ => nstep, &
|
|
|
|
iprint_ => iprint, &
|
|
|
|
nosym_ => nosym, &
|
|
|
|
modenum_ => modenum, &
|
2004-01-23 17:50:00 +08:00
|
|
|
reduce_io, ethr, lscf, lbfgs, lmd, lneb, lphonon, &
|
2004-02-26 19:50:36 +08:00
|
|
|
noinv, time_max, restart, loldbfgs, lconstrain, &
|
|
|
|
ldamped
|
2003-12-10 22:57:07 +08:00
|
|
|
USE wvfct, ONLY : ibm_baco2, &
|
|
|
|
nbnd_ => nbnd
|
2003-10-17 00:30:12 +08:00
|
|
|
USE fixed_occ, ONLY : tfixed_occ
|
|
|
|
USE control_flags, ONLY : twfcollect
|
2003-12-10 22:57:07 +08:00
|
|
|
USE neb_variables, ONLY : lsteep_des, lquick_min , ldamped_dyn, lmol_dyn, &
|
|
|
|
num_of_images_ => num_of_images, &
|
|
|
|
CI_scheme_ => CI_scheme, &
|
|
|
|
VEC_scheme_ => VEC_scheme, &
|
|
|
|
optimization_ => optimization, &
|
|
|
|
damp_ => damp, &
|
|
|
|
temp_req_ => temp_req, &
|
|
|
|
ds_ => ds, &
|
|
|
|
k_max_ => k_max, &
|
|
|
|
k_min_ => k_min, &
|
|
|
|
neb_thr_ => neb_thr
|
|
|
|
USE noncollin_module, ONLY : baco_ibm_xlf, &
|
|
|
|
noncolin_ => noncolin, &
|
|
|
|
lambda_ => lambda, &
|
|
|
|
i_cons_ => i_cons, &
|
|
|
|
mcons_ => mcons, &
|
|
|
|
angle1_ => angle1, &
|
|
|
|
angle2_ => angle2, &
|
|
|
|
report_ => report
|
2004-01-24 18:01:55 +08:00
|
|
|
USE bfgs_module, ONLY : bfgs_xlf_bug, &
|
|
|
|
lbfgs_ndim_ => lbfgs_ndim, &
|
2003-12-10 22:57:07 +08:00
|
|
|
trust_radius_max_ => trust_radius_max, &
|
|
|
|
trust_radius_min_ => trust_radius_min, &
|
|
|
|
trust_radius_ini_ => trust_radius_ini, &
|
|
|
|
trust_radius_end_ => trust_radius_end, &
|
|
|
|
w_1_ => w_1, &
|
|
|
|
w_2_ => w_2
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! CONTROL namelist
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
USE input_parameters, ONLY : title, calculation, verbosity, &
|
|
|
|
restart_mode, nstep, iprint, tstress, tprnfor, &
|
|
|
|
dt, outdir, prefix, max_seconds, &
|
|
|
|
etot_conv_thr, forc_conv_thr, pseudo_dir, &
|
|
|
|
disk_io, tefield, dipfield, lberry, gdir, &
|
|
|
|
nppstr, wf_collect
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
! SYSTEM namelist
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
USE input_parameters, ONLY : ibrav, celldm, a, b, c, cosab, cosac, cosbc, &
|
|
|
|
nat, ntyp, nbnd, nelec, ecutwfc, ecutrho, &
|
|
|
|
nr1, nr2, nr3, nr1s, nr2s, nr3s, &
|
|
|
|
nosym, starting_magnetization, &
|
|
|
|
occupations, degauss, smearing, &
|
|
|
|
nspin, ecfixed, qcutz, q2sigma, &
|
|
|
|
lda_plus_U, Hubbard_U, Hubbard_alpha, &
|
2004-02-14 16:39:34 +08:00
|
|
|
edir, emaxpos, eopreg, eamp, starting_ns_eigenvalue, &
|
2003-10-30 21:56:34 +08:00
|
|
|
noncolin, lambda, i_cons, mcons, angle1, &
|
|
|
|
angle2, report
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
! ELECTRONS namelist
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
USE input_parameters, ONLY : electron_maxstep, electron_dynamics, &
|
|
|
|
mixing_mode, mixing_beta, mixing_ndim, &
|
|
|
|
mixing_fixed_ns, diago_cg_maxiter, &
|
|
|
|
diago_david_ndim, diago_diis_buff, &
|
|
|
|
diago_diis_ndim, diago_diis_keep, &
|
|
|
|
diago_diis_ethr, diagonalization, startingwfc, &
|
2004-01-23 17:50:00 +08:00
|
|
|
startingpot, conv_thr, diago_thr_init
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
! IONS namelist
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
USE input_parameters, ONLY : ion_dynamics, ion_positions, ion_temperature, &
|
2003-12-10 22:57:07 +08:00
|
|
|
tempw, tolp, upscale, potential_extrapolation, &
|
|
|
|
CI_scheme, VEC_scheme, minimization_scheme, &
|
|
|
|
num_of_images, optimization, damp, temp_req, &
|
|
|
|
ds, k_max, k_min, neb_thr, &
|
|
|
|
trust_radius_max, trust_radius_min, &
|
|
|
|
trust_radius_ini, trust_radius_end, &
|
|
|
|
w_1, w_2, lbfgs_ndim
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
! CELL namelist
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
USE input_parameters, ONLY : cell_parameters, cell_dynamics, press, &
|
|
|
|
wmass, cell_temperature, cell_dofree, &
|
|
|
|
cell_factor
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
! PHONON namelist
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
USE input_parameters, ONLY : phonon, modenum, xqq
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... NEB specific
|
|
|
|
!
|
|
|
|
USE input_parameters, ONLY : pos
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
USE read_namelists_module, ONLY: read_namelists
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-18 22:28:27 +08:00
|
|
|
#if defined (__PARA)
|
|
|
|
USE para, ONLY : npool
|
|
|
|
#endif
|
2003-10-17 00:30:12 +08:00
|
|
|
IMPLICIT NONE
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... local variables
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-16 23:43:57 +08:00
|
|
|
INTEGER :: unit = 5, &! standard input unit
|
2004-02-14 16:39:34 +08:00
|
|
|
i, iiarg, nargs, ia, ios, ierr, ilen, is, image, nt
|
2003-12-16 23:43:57 +08:00
|
|
|
INTEGER :: iargc
|
|
|
|
EXTERNAL iargc
|
|
|
|
CHARACTER (LEN=80) :: input_file
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
!
|
|
|
|
CALL getenv( 'HOME', pseudo_dir )
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
pseudo_dir = TRIM( pseudo_dir ) // '/pw/pseudo/'
|
2003-04-03 23:35:36 +08:00
|
|
|
!
|
2003-12-16 23:43:57 +08:00
|
|
|
! ... Input from file ?
|
|
|
|
!
|
|
|
|
nargs = iargc()
|
|
|
|
!
|
|
|
|
DO iiarg = 1, ( nargs - 1 )
|
|
|
|
!
|
|
|
|
CALL getarg( iiarg, input_file )
|
|
|
|
IF ( TRIM( input_file ) == '-input' .OR. &
|
|
|
|
TRIM( input_file ) == '-inp' .OR. &
|
|
|
|
TRIM( input_file ) == '-in' ) THEN
|
|
|
|
!
|
|
|
|
CALL getarg( ( iiarg + 1 ) , input_file )
|
|
|
|
OPEN ( UNIT = unit, FILE = input_file, FORM = 'FORMATTED', &
|
|
|
|
STATUS = 'OLD', IOSTAT = ierr )
|
|
|
|
CALL errore( 'iosys', 'input file ' // TRIM( input_file ) // &
|
|
|
|
& ' not found' , ierr )
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
|
|
|
!
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
CALL read_namelists( 'PW' )
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
nraise = 100
|
|
|
|
delta_t = 1.D0
|
|
|
|
!
|
|
|
|
! ... translate from input to internals of PWscf, various checks
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
time_max = max_seconds
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
IF ( tefield .AND. ( .NOT. nosym ) ) THEN
|
|
|
|
nosym = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
WRITE( stdout, &
|
|
|
|
'(5x,"Presently no symmetry can be used with electric field",/)' )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( tefield .AND. tstress ) THEN
|
2003-10-17 00:30:12 +08:00
|
|
|
tstress = .FALSE.
|
2003-12-10 22:57:07 +08:00
|
|
|
WRITE( stdout, &
|
|
|
|
'(5x,"Presently stress not available with electric field",/)' )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
IF ( tefield .AND. ( nspin == 2 ) ) THEN
|
2004-02-18 22:28:27 +08:00
|
|
|
CALL errore( 'iosys', 'LSDA not available with electric field' , 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
!
|
2003-09-02 17:16:09 +08:00
|
|
|
twfcollect = wf_collect
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
! ... Set Values for electron and bands
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
tfixed_occ = .FALSE.
|
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( occupations ) )
|
|
|
|
CASE ( 'fixed' )
|
2003-01-20 05:58:50 +08:00
|
|
|
ngauss = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
ltetra = .FALSE.
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( degauss /= 0.D0 ) THEN
|
|
|
|
CALL errore( ' iosys ', &
|
|
|
|
& ' fixed occupations, gauss. broadening ignored', -1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
degauss = 0.D0
|
2003-01-20 05:58:50 +08:00
|
|
|
END IF
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'smearing' )
|
|
|
|
ltetra = .FALSE.
|
|
|
|
IF ( degauss == 0.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', &
|
|
|
|
& ' smearing requires gaussian broadening', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END IF
|
2003-10-17 00:30:12 +08:00
|
|
|
SELECT CASE ( TRIM( smearing ) )
|
|
|
|
CASE ( 'gaussian' , 'gauss' )
|
2003-01-20 05:58:50 +08:00
|
|
|
ngauss = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'methfessel-paxton' , 'm-p' , 'mp' )
|
2003-01-20 05:58:50 +08:00
|
|
|
ngauss = 1
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'marzari-vanderbilt' , 'cold' , 'm-v' , 'mv' )
|
|
|
|
ngauss = -1
|
|
|
|
CASE ( 'fermi-dirac' , 'f-d' , 'fd' )
|
2003-01-20 05:58:50 +08:00
|
|
|
ngauss = -99
|
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'tetrahedra' )
|
2003-01-20 05:58:50 +08:00
|
|
|
ngauss = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
ltetra = .TRUE.
|
|
|
|
CASE ( 'from_input' )
|
|
|
|
ngauss = 0
|
|
|
|
ltetra = .FALSE.
|
|
|
|
tfixed_occ = .TRUE.
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ',' occupations ' // TRIM( occupations ) // &
|
|
|
|
& 'not implemented', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
IF( nbnd < 1 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' nbnd less than 1 ', nbnd )
|
2003-01-20 05:58:50 +08:00
|
|
|
END IF
|
|
|
|
IF( nelec < 0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' nelec less than 0 ', nelec )
|
2003-01-20 05:58:50 +08:00
|
|
|
END IF
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
lsda = ( nspin == 2 )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... starting_magnetization(ia) = -2.D0 means "not set" -- set it to 0
|
|
|
|
!
|
2003-12-02 22:28:22 +08:00
|
|
|
DO ia = 1, ntyp
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( starting_magnetization(ia) == -2.D0 ) &
|
|
|
|
starting_magnetization(ia) = 0.d0
|
|
|
|
END DO
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
IF ( ecutrho <= 0.D0 ) THEN
|
|
|
|
dual = 4.D0
|
|
|
|
ELSE
|
|
|
|
dual = ecutrho / ecutwfc
|
|
|
|
IF( dual <= 1.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' invalid dual? ', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( restart_mode ) )
|
|
|
|
CASE ( 'from_scratch' )
|
|
|
|
restart = .FALSE.
|
|
|
|
restart_bfgs = .FALSE.
|
2003-01-20 05:58:50 +08:00
|
|
|
startingconfig = 'input'
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'restart' )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... NEB specific
|
|
|
|
!
|
|
|
|
IF ( calculation == 'neb' ) THEN
|
|
|
|
restart = .FALSE.
|
|
|
|
restart_bfgs = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE
|
2003-12-10 22:57:07 +08:00
|
|
|
restart = .TRUE.
|
|
|
|
restart_bfgs = .TRUE.
|
|
|
|
startingpot = 'file'
|
|
|
|
startingwfc = 'file'
|
|
|
|
IF ( TRIM( ion_positions ) == 'from_input' ) THEN
|
|
|
|
startingconfig = 'input'
|
|
|
|
ELSE
|
|
|
|
startingconfig = 'file'
|
|
|
|
END IF
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' unknown restart_mode ' // &
|
|
|
|
& TRIM( restart_mode ), 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( disk_io ) )
|
|
|
|
CASE ( 'high' )
|
|
|
|
reduce_io = .FALSE.
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-10-17 00:30:12 +08:00
|
|
|
reduce_io = .TRUE.
|
|
|
|
restart = .FALSE.
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
Hubbard_U(:) = Hubbard_U(:) / rytoev
|
|
|
|
Hubbard_alpha(:)= Hubbard_alpha(:) / rytoev
|
|
|
|
!
|
2004-01-23 17:50:00 +08:00
|
|
|
ethr = diago_thr_init
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( calculation ) )
|
|
|
|
CASE ( 'scf' )
|
|
|
|
lscf = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
lbfgs = .FALSE.
|
|
|
|
lmd = .FALSE.
|
|
|
|
lneb = .FALSE.
|
|
|
|
iswitch = 0 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
2003-10-17 00:30:12 +08:00
|
|
|
lforce = tprnfor
|
|
|
|
lmovecell = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
nstep = 1
|
|
|
|
CASE ( 'nscf' )
|
|
|
|
lscf = .FALSE.
|
2003-12-10 22:57:07 +08:00
|
|
|
lbfgs = .FALSE.
|
|
|
|
lmd = .FALSE.
|
|
|
|
lneb = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
iswitch = -1
|
|
|
|
lforce = .FALSE.
|
|
|
|
lmovecell = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
nstep = 1
|
|
|
|
CASE ( 'relax' )
|
|
|
|
lscf = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
lbfgs = .TRUE.
|
|
|
|
lmd = .FALSE.
|
|
|
|
lneb = .FALSE.
|
|
|
|
iswitch = 1 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
2003-10-17 00:30:12 +08:00
|
|
|
lforce = .TRUE.
|
|
|
|
lmovecell = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'md' )
|
|
|
|
lscf = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
lbfgs = .FALSE.
|
|
|
|
lmd = .TRUE.
|
|
|
|
lneb = .FALSE.
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
2003-10-17 00:30:12 +08:00
|
|
|
lforce = .TRUE.
|
|
|
|
lmovecell = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'vc-relax' , 'vc-md' )
|
|
|
|
lscf = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
lbfgs = .FALSE.
|
|
|
|
lmd = .TRUE.
|
|
|
|
lneb = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
iswitch = 3
|
|
|
|
lmovecell = .TRUE.
|
|
|
|
lforce = .TRUE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'phonon' )
|
|
|
|
lscf = .FALSE.
|
2003-12-10 22:57:07 +08:00
|
|
|
lbfgs = .FALSE.
|
|
|
|
lmd = .FALSE.
|
|
|
|
lneb = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
iswitch = -2 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
2003-10-17 00:30:12 +08:00
|
|
|
lforce = .FALSE.
|
|
|
|
lmovecell = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .TRUE.
|
2003-10-17 00:30:12 +08:00
|
|
|
nstep = 1
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... NEB specific
|
|
|
|
!
|
|
|
|
CASE ( 'neb' )
|
|
|
|
lscf = .TRUE.
|
|
|
|
lbfgs = .FALSE.
|
|
|
|
lmd = .FALSE.
|
|
|
|
lneb = .TRUE.
|
|
|
|
iswitch = 1 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
lforce = tprnfor
|
|
|
|
lmovecell = .FALSE.
|
2004-01-23 17:50:00 +08:00
|
|
|
lphonon = .FALSE.
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' calculation ' // &
|
|
|
|
& TRIM( calculation ) // ' not implemented', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
IF ( modenum /= 0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = - 4
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
!
|
|
|
|
IF ( startingpot /= 'atomic' .AND. startingpot /= 'file' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys', 'wrong startingpot: use default', -1 )
|
|
|
|
IF ( lscf ) startingpot = 'atomic'
|
|
|
|
IF ( .NOT. lscf ) startingpot = 'file'
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
!
|
|
|
|
IF ( .NOT. lscf .AND. startingpot /= 'file' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys', 'wrong startingpot: use default', -1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
startingpot = 'file'
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
!
|
|
|
|
IF ( startingwfc /= 'atomic' .AND. startingwfc /= 'random' .AND. &
|
|
|
|
startingwfc /= 'file' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys', 'wrong startingwfc: use default', -1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
startingwfc = 'atomic'
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( electron_dynamics ) )
|
|
|
|
CASE ( 'none' )
|
|
|
|
CONTINUE
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ',' unknown electron_dynamics '// &
|
|
|
|
& TRIM( electron_dynamics ), 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( diagonalization ) )
|
|
|
|
CASE ( 'cg' )
|
2003-01-20 05:58:50 +08:00
|
|
|
isolve = 1
|
2003-10-17 00:30:12 +08:00
|
|
|
max_cg_iter = diago_cg_maxiter
|
|
|
|
CASE ( 'diis' )
|
2003-01-20 05:58:50 +08:00
|
|
|
isolve = 2
|
2003-10-17 00:30:12 +08:00
|
|
|
max_cg_iter = diago_cg_maxiter
|
|
|
|
diis_buff = diago_diis_buff
|
|
|
|
diis_ethr_cg = diago_diis_ethr ! ... SF
|
|
|
|
diis_ndim = diago_diis_ndim ! ... SF
|
|
|
|
diis_wfc_keep = diago_diis_keep
|
|
|
|
CASE ( 'david' )
|
2003-01-20 05:58:50 +08:00
|
|
|
isolve = 0
|
|
|
|
david = diago_david_ndim
|
|
|
|
CASE DEFAULT
|
|
|
|
isolve = 0
|
|
|
|
david = diago_david_ndim
|
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
tr2 = conv_thr
|
|
|
|
niter = electron_maxstep
|
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( potential_extrapolation ) )
|
|
|
|
CASE ( 'none' )
|
2003-01-20 05:58:50 +08:00
|
|
|
order = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'atomic' )
|
2003-01-20 05:58:50 +08:00
|
|
|
order = 1
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'wfc' )
|
2003-01-20 05:58:50 +08:00
|
|
|
order = 2
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'wfc2' )
|
2003-01-20 05:58:50 +08:00
|
|
|
order = 3
|
|
|
|
CASE DEFAULT
|
|
|
|
order = 1
|
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
IF ( occupations == 'fixed' .AND. nspin == 2 .AND. lscf ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', &
|
|
|
|
& ' fixed occupations and lsda not implemented ', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END IF
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
! ... initialization
|
|
|
|
!
|
|
|
|
calc = ' '
|
|
|
|
lbfgs = .FALSE.
|
|
|
|
loldbfgs = .FALSE.
|
|
|
|
lmd = .FALSE.
|
|
|
|
ldamped = .FALSE.
|
|
|
|
lconstrain = .FALSE.
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
IF ( TRIM( calculation ) == 'relax' ) THEN
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
SELECT CASE ( TRIM( ion_dynamics ) )
|
|
|
|
CASE ( 'bfgs' )
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = 1 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
lbfgs = .TRUE.
|
2004-02-26 19:50:36 +08:00
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2003-12-18 21:39:53 +08:00
|
|
|
IF ( epse <= 20.D0 * ( tr2 / upscale ) ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' required etot_conv_thr is too small:' // &
|
|
|
|
& ' conv_thr must be reduced', 1 )
|
|
|
|
!
|
2004-02-09 19:15:33 +08:00
|
|
|
CASE ( 'old-bfgs' )
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = 1 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
2004-02-09 19:15:33 +08:00
|
|
|
loldbfgs = .TRUE.
|
2004-02-26 19:50:36 +08:00
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2003-12-18 21:39:53 +08:00
|
|
|
IF ( epse <= 20.D0 * ( tr2 / upscale ) ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' required etot_conv_thr is too small:' // &
|
|
|
|
& ' conv_thr must be reduced', 1 )
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
CASE ( 'constrained-damp' )
|
|
|
|
iswitch = 4 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
lmd = .TRUE.
|
|
|
|
ldamped = .TRUE.
|
|
|
|
lconstrain = .TRUE.
|
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'damp' )
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
lmd = .TRUE.
|
2004-02-26 19:50:36 +08:00
|
|
|
ldamped = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
|
|
|
ntcheck = nstep + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': ion_dynamics=' // TRIM( ion_dynamics ) // &
|
|
|
|
& ' not supported', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
|
|
|
ELSE IF ( TRIM( calculation ) == 'md' ) THEN
|
|
|
|
!
|
|
|
|
lmd = .TRUE.
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
SELECT CASE ( TRIM( ion_dynamics ) )
|
|
|
|
CASE ( 'verlet' )
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
CASE ( 'constrained-verlet' )
|
2004-02-26 19:50:36 +08:00
|
|
|
lconstrain = .TRUE.
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = 4 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
CASE ( 'beeman' )
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
2003-01-20 05:58:50 +08:00
|
|
|
calc = 'md'
|
2003-12-10 22:57:07 +08:00
|
|
|
ntcheck = nstep + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': ion_dynamics=' // TRIM( ion_dynamics ) // &
|
|
|
|
& ' not supported', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
|
|
|
ELSE IF ( TRIM( calculation ) == 'vc-relax' ) THEN
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
SELECT CASE ( TRIM( cell_dynamics ) )
|
|
|
|
CASE ( 'none' )
|
2003-12-10 22:57:07 +08:00
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
calc = 'mm'
|
2003-10-17 00:30:12 +08:00
|
|
|
ntcheck = nstep + 1
|
|
|
|
CASE ( 'damp-pr' )
|
2003-12-10 22:57:07 +08:00
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
calc = 'cm'
|
2003-10-17 00:30:12 +08:00
|
|
|
ntcheck = nstep + 1
|
|
|
|
CASE ( 'damp-w' )
|
2003-12-10 22:57:07 +08:00
|
|
|
epse = etot_conv_thr
|
|
|
|
epsf = forc_conv_thr
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
calc = 'nm'
|
2003-10-17 00:30:12 +08:00
|
|
|
ntcheck = nstep + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': cell_dynamics=' // TRIM( cell_dynamics ) // &
|
|
|
|
& ' not supported', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( TRIM( ion_dynamics ) /= 'damp' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': ion_dynamics=' // TRIM( ion_dynamics ) // &
|
|
|
|
& ' not supported', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
|
|
|
ELSE IF ( TRIM( calculation ) == 'vc-md' ) THEN
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
SELECT CASE ( TRIM( cell_dynamics ) )
|
2003-12-10 22:57:07 +08:00
|
|
|
CASE ( 'none' )
|
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
calc = 'md'
|
2003-10-17 00:30:12 +08:00
|
|
|
ntcheck = nstep + 1
|
|
|
|
CASE ( 'pr' )
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
calc = 'cd'
|
2003-10-17 00:30:12 +08:00
|
|
|
ntcheck = nstep + 1
|
|
|
|
CASE ( 'w' )
|
2003-12-10 22:57:07 +08:00
|
|
|
iswitch = 3 ! ... obsolescent: do not use in new code ( 29/10/2003 C.S.)
|
|
|
|
calc = 'nd'
|
2003-10-17 00:30:12 +08:00
|
|
|
ntcheck = nstep + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': ion_dynamics=' // TRIM( ion_dynamics ) // &
|
|
|
|
& ' not supported', 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( TRIM( ion_dynamics ) /= 'beeman' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': ion_dynamics=' // TRIM( ion_dynamics ) // &
|
|
|
|
& ' not supported', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
|
|
|
ELSE IF ( TRIM( calculation ) == 'neb' ) THEN
|
|
|
|
!
|
|
|
|
! ... NEB specific
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
IF ( num_of_images < 2 ) THEN
|
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': num_of_images must be at least 2', 1 )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
IF ( ( CI_scheme /= "no-CI" ) .AND. &
|
|
|
|
( CI_scheme /= "highest-TS" ) .AND. &
|
|
|
|
( CI_scheme /= "all-SP" ) .AND. &
|
|
|
|
( CI_scheme /= "manual" ) ) THEN
|
|
|
|
!
|
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': unknown CI_scheme', 1 )
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
IF ( ( VEC_scheme /= "energy-weighted" ) .AND. &
|
|
|
|
( VEC_scheme /= "gradient-weighted" ) ) THEN
|
|
|
|
!
|
|
|
|
CALL errore( ' iosys ', 'calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': unknown VEC_scheme', 1 )
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
! ... initialization
|
|
|
|
!
|
|
|
|
lsteep_des = .FALSE.
|
|
|
|
lquick_min = .FALSE.
|
|
|
|
ldamped_dyn = .FALSE.
|
|
|
|
lmol_dyn = .FALSE.
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
SELECT CASE ( minimization_scheme )
|
2003-12-10 22:57:07 +08:00
|
|
|
CASE ( "sd" )
|
|
|
|
lsteep_des = .TRUE.
|
|
|
|
CASE ( "quick-min" )
|
|
|
|
lquick_min = .TRUE.
|
|
|
|
CASE ( "damped-dyn" )
|
|
|
|
ldamped_dyn = .TRUE.
|
|
|
|
CASE ( "mol-dyn" )
|
|
|
|
lmol_dyn = .TRUE.
|
|
|
|
IF ( temp_req == 0 ) &
|
|
|
|
WRITE( stdout,'(/,T2,"WARNING: tepm_req has not been set" )')
|
|
|
|
temp_req = temp_req / ( eV_to_kelvin * AU )
|
|
|
|
CASE default
|
|
|
|
CALL errore( ' iosys ','calculation=' // TRIM( calculation ) // &
|
|
|
|
& ': unknown minimization_scheme', 1 )
|
|
|
|
END SELECT
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
SELECT CASE ( TRIM( ion_temperature ) )
|
|
|
|
CASE ( 'not_controlled' )
|
|
|
|
CONTINUE
|
|
|
|
CASE ( 'rescaling' )
|
2003-01-20 05:58:50 +08:00
|
|
|
temperature = tempw
|
2003-10-17 00:30:12 +08:00
|
|
|
ttol = tolp
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' unknown ion_temperature ' // &
|
|
|
|
& TRIM( ion_temperature ), 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( cell_temperature ) )
|
|
|
|
CASE ( 'not_controlled' )
|
|
|
|
CONTINUE
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' unknown cell_temperature ' // &
|
|
|
|
& TRIM( cell_temperature ), 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( cell_dofree ) )
|
|
|
|
CASE ( 'all' )
|
|
|
|
CONTINUE
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', &
|
|
|
|
& ' unknown cell_dofree ' // TRIM( cell_dofree ), 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( mixing_mode ) )
|
|
|
|
CASE ( 'plain' )
|
2003-01-20 05:58:50 +08:00
|
|
|
imix = 0
|
|
|
|
starting_scf_threshold = tr2
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'TF' )
|
2003-01-20 05:58:50 +08:00
|
|
|
imix = 1
|
|
|
|
starting_scf_threshold = tr2
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'local-TF' )
|
2003-01-20 05:58:50 +08:00
|
|
|
imix = 2
|
|
|
|
starting_scf_threshold = tr2
|
2003-10-17 00:30:12 +08:00
|
|
|
CASE ( 'potential' )
|
2003-01-20 05:58:50 +08:00
|
|
|
imix = -1
|
2003-12-10 22:57:07 +08:00
|
|
|
starting_scf_threshold = SQRT( tr2 )
|
2003-01-20 05:58:50 +08:00
|
|
|
CASE DEFAULT
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' unknown mixing ' // TRIM( mixing_mode ), 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
END SELECT
|
2004-02-16 19:47:26 +08:00
|
|
|
|
|
|
|
IF ( dipfield .AND. imix.ne.-1 ) THEN
|
2004-02-18 22:28:27 +08:00
|
|
|
CALL errore( 'iosys', 'use mixing_mod=potential with dipfield' , 1 )
|
2004-02-16 19:47:26 +08:00
|
|
|
END IF
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
nmix = mixing_ndim
|
|
|
|
niter_with_fixed_ns = mixing_fixed_ns
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
SELECT CASE ( TRIM( verbosity ) )
|
|
|
|
CASE ( 'high' )
|
2003-01-20 05:58:50 +08:00
|
|
|
iverbosity = 1
|
|
|
|
CASE DEFAULT
|
|
|
|
iverbosity = 0
|
|
|
|
END SELECT
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
tmp_dir = TRIM( outdir )
|
|
|
|
lstres = ( tstress .AND. lscf )
|
|
|
|
!
|
|
|
|
! ... Copy values from input module to PW internals
|
|
|
|
!
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
nppstr_ = nppstr
|
|
|
|
gdir_ = gdir
|
|
|
|
lberry_ = lberry
|
2004-02-18 22:28:27 +08:00
|
|
|
#ifdef __PARA
|
|
|
|
if ( lberry_ .AND. npool > 1 ) &
|
|
|
|
CALL errore (' iosys ', ' Berry Phase not implemented with pools ', 1)
|
|
|
|
#endif
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
title_ = title
|
|
|
|
dt_ = dt
|
|
|
|
tefield_ = tefield
|
|
|
|
dipfield_ = dipfield
|
2003-10-17 00:30:12 +08:00
|
|
|
prefix_ = TRIM( prefix )
|
|
|
|
pseudo_dir_ = TRIM( pseudo_dir )
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
nstep_ = nstep
|
|
|
|
iprint_ = iprint
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
celldm_ = celldm
|
|
|
|
ibrav_ = ibrav
|
|
|
|
nat_ = nat
|
|
|
|
ntyp_ = ntyp
|
|
|
|
edir_ = edir
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
emaxpos_ = emaxpos
|
2003-10-17 00:30:12 +08:00
|
|
|
eopreg_ = eopreg
|
|
|
|
eamp_ = eamp
|
|
|
|
nr1_ = nr1
|
|
|
|
nr2_ = nr2
|
|
|
|
nr3_ = nr3
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
ecutwfc_ = ecutwfc
|
|
|
|
ecfixed_ = ecfixed
|
2003-10-17 00:30:12 +08:00
|
|
|
qcutz_ = qcutz
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
q2sigma_ = q2sigma
|
2003-10-17 00:30:12 +08:00
|
|
|
nr1s_ = nr1s
|
|
|
|
nr2s_ = nr2s
|
|
|
|
nr3s_ = nr3s
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
degauss_ = degauss
|
2003-10-17 00:30:12 +08:00
|
|
|
nelec_ = nelec
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2003-10-30 21:56:34 +08:00
|
|
|
noncolin_ = noncolin
|
2003-12-10 22:57:07 +08:00
|
|
|
angle1_ = angle1
|
|
|
|
angle2_ = angle2
|
|
|
|
report_ = report
|
|
|
|
i_cons_ = i_cons
|
|
|
|
mcons_ = mcons
|
|
|
|
lambda_ = lambda
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
Hubbard_U_( 1 : ntyp ) = hubbard_u( 1 : ntyp )
|
|
|
|
Hubbard_alpha_( 1 : ntyp ) = hubbard_alpha( 1 : ntyp )
|
|
|
|
lda_plus_u_ = lda_plus_u
|
|
|
|
nspin_ = nspin
|
|
|
|
starting_magnetization_ = starting_magnetization
|
2004-02-14 16:39:34 +08:00
|
|
|
starting_ns = starting_ns_eigenvalue
|
2003-10-17 00:30:12 +08:00
|
|
|
nosym_ = nosym
|
|
|
|
nbnd_ = nbnd
|
|
|
|
!
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
startingwfc_ = startingwfc
|
|
|
|
startingpot_ = startingpot
|
|
|
|
mixing_beta_ = mixing_beta
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
upscale_ = upscale
|
|
|
|
press_ = press
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
cell_factor_ = cell_factor
|
2003-10-17 00:30:12 +08:00
|
|
|
modenum_ = modenum
|
|
|
|
xqq_ = xqq
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... NEB specific
|
|
|
|
!
|
|
|
|
num_of_images_ = num_of_images
|
|
|
|
CI_scheme_ = CI_scheme
|
|
|
|
VEC_scheme_ = VEC_scheme
|
|
|
|
optimization_ = optimization
|
|
|
|
damp_ = damp
|
|
|
|
temp_req_ = temp_req
|
|
|
|
ds_ = ds
|
|
|
|
k_max_ = k_max
|
|
|
|
k_min_ = k_min
|
|
|
|
neb_thr_ = neb_thr
|
|
|
|
!
|
|
|
|
! ... new BFGS specific
|
|
|
|
!
|
|
|
|
lbfgs_ndim_ = lbfgs_ndim
|
|
|
|
trust_radius_max_ = trust_radius_max
|
|
|
|
trust_radius_min_ = trust_radius_min
|
|
|
|
trust_radius_ini_ = trust_radius_ini
|
|
|
|
trust_radius_end_ = trust_radius_end
|
|
|
|
w_1_ = w_1
|
|
|
|
w_2_ = w_2
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... read following cards
|
|
|
|
!
|
|
|
|
ALLOCATE( tau( 3, nat_ ) )
|
|
|
|
ALLOCATE( ityp( nat_ ) )
|
|
|
|
ALLOCATE( force( 3, nat_ ) ) ! ... compatibility with old readin
|
|
|
|
ALLOCATE( if_pos( 3, nat_ ) )
|
|
|
|
!
|
|
|
|
IF ( tefield ) ALLOCATE( forcefield( 3, nat_ ) )
|
|
|
|
!
|
|
|
|
CALL read_cards( psfile, atomic_positions )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... set up atomic positions and crystal lattice
|
2003-07-11 23:16:00 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( celldm_(1) == 0.D0 .AND. a /= 0.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( ibrav_ == 0 ) ibrav = 14
|
|
|
|
celldm_(1) = a / bohr_radius_angs
|
|
|
|
celldm_(2) = b / a
|
|
|
|
celldm_(3) = c / a
|
|
|
|
celldm_(4) = cosab
|
|
|
|
celldm_(5) = cosac
|
|
|
|
celldm_(6) = cosbc
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE IF ( celldm_(1) /= 0.D0 .AND. a /= 0.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( 'input', ' do not specify both celldm and a,b,c!', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( ibrav_ == 0 .AND. celldm_(1) /= 0.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... input at are in units of alat
|
|
|
|
!
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
alat = celldm_(1)
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE IF ( ibrav_ == 0 .AND. celldm_(1) == 0.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... input at are in atomic units: define alat
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
celldm_(1) = SQRT( at(1,1)**2 + at(1,2)**2 + at(1,3)**2 )
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
alat = celldm_(1)
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... bring at to alat units
|
|
|
|
!
|
2003-07-28 23:03:32 +08:00
|
|
|
at(:,:) = at(:,:) / alat
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... generate at (atomic units)
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
CALL latgen( ibrav, celldm_, at(1,1), at(1,2), at(1,3), omega )
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
alat = celldm_(1)
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... bring at to alat units
|
|
|
|
!
|
2003-07-28 23:03:32 +08:00
|
|
|
at(:,:) = at(:,:) / alat
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
CALL volume( alat, at(1,1), at(1,2), at(1,3), omega )
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( calculation == 'neb' ) THEN
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... NEB specific
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
DO image = 1, num_of_images_
|
|
|
|
!
|
|
|
|
tau = RESHAPE( SOURCE = pos(1:3*nat_,image), SHAPE = (/ 3 , nat_ /) )
|
|
|
|
!
|
|
|
|
SELECT CASE ( atomic_positions )
|
|
|
|
!
|
|
|
|
! ... convert input atomic positions to internally used format:
|
|
|
|
! ... tau in a0 units
|
|
|
|
!
|
|
|
|
CASE ( 'alat' )
|
|
|
|
!
|
|
|
|
! ... input atomic positions are divided by a0: do nothing
|
|
|
|
!
|
|
|
|
CASE ( 'bohr' )
|
|
|
|
!
|
|
|
|
! ... input atomic positions are in a.u.: divide by alat
|
|
|
|
!
|
|
|
|
tau = tau / alat
|
|
|
|
CASE ( 'crystal' )
|
|
|
|
!
|
|
|
|
! ... input atomic positions are in crystal axis
|
|
|
|
!
|
|
|
|
CALL cryst_to_cart( nat_, tau, at, 1 )
|
|
|
|
CASE ( 'angstrom' )
|
|
|
|
!
|
|
|
|
! ... atomic positions in A: convert to a.u. and divide by alat
|
|
|
|
!
|
|
|
|
tau = tau / bohr_radius_angs / alat
|
|
|
|
CASE DEFAULT
|
|
|
|
CALL errore( ' iosys ',' atomic_positions=' // &
|
|
|
|
& TRIM( atomic_positions ) // ' not implemented ', 1 )
|
|
|
|
END SELECT
|
|
|
|
!
|
|
|
|
pos(1:3*nat_,image) = RESHAPE( SOURCE = tau, SHAPE = (/ 3 * nat_ /) )
|
|
|
|
!
|
|
|
|
END DO
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
ELSE
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
SELECT CASE ( atomic_positions )
|
|
|
|
!
|
|
|
|
! ... convert input atomic positions to internally used format:
|
|
|
|
! ... tau in a0 units
|
|
|
|
!
|
|
|
|
CASE ( 'alat' )
|
|
|
|
!
|
|
|
|
! ... input atomic positions are divided by a0: do nothing
|
|
|
|
!
|
|
|
|
CASE ( 'bohr' )
|
|
|
|
!
|
|
|
|
! ... input atomic positions are in a.u.: divide by alat
|
|
|
|
!
|
|
|
|
tau = tau / alat
|
|
|
|
CASE ( 'crystal' )
|
|
|
|
!
|
|
|
|
! ... input atomic positions are in crystal axis
|
|
|
|
!
|
|
|
|
CALL cryst_to_cart( nat_, tau, at, 1 )
|
|
|
|
CASE ( 'angstrom' )
|
|
|
|
!
|
|
|
|
! ... atomic positions in A: convert to a.u. and divide by alat
|
|
|
|
!
|
|
|
|
tau = tau / bohr_radius_angs / alat
|
|
|
|
CASE DEFAULT
|
|
|
|
CALL errore( ' iosys ',' atomic_positions=' // &
|
|
|
|
& TRIM( atomic_positions ) // ' not implemented ', 1 )
|
|
|
|
END SELECT
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... set default value of wmass
|
2003-01-29 22:29:11 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( wmass == 0.D0 ) THEN
|
|
|
|
IF ( calc == 'nd' .OR. calc == 'nm' ) THEN
|
|
|
|
DO ia = 1, nat
|
2003-01-29 22:29:11 +08:00
|
|
|
wmass = wmass + amass(ityp(ia))
|
2003-10-17 00:30:12 +08:00
|
|
|
END DO
|
2003-12-10 22:57:07 +08:00
|
|
|
wmass = 0.75D0 * wmass / pi / pi / omega**( 2.D0 / 3.D0 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
IF ( calc == 'cd' .OR. calc == 'cm' ) THEN
|
|
|
|
DO ia = 1, nat
|
2003-01-29 22:29:11 +08:00
|
|
|
wmass = wmass + amass(ityp(ia))
|
2003-10-17 00:30:12 +08:00
|
|
|
END DO
|
2003-12-10 22:57:07 +08:00
|
|
|
wmass = 0.75D0 * wmass / pi / pi
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
END IF
|
2003-01-29 22:29:11 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... unit conversion for cell mass and pressure
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
cmass = wmass * amconv
|
|
|
|
press_ = press_ / uakbar
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... read pseudopotentials
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL readpp()
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... Renata's dynamics uses masses in atomic units
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( calc /= ' ' ) THEN
|
2003-01-20 05:58:50 +08:00
|
|
|
amass = amass * amconv
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... In the case of variable cell dynamics save old cell variables
|
|
|
|
! ... and initialize a few other variables
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( lmovecell ) THEN
|
|
|
|
at_old = at
|
2003-01-20 05:58:50 +08:00
|
|
|
omega_old = omega
|
2003-10-17 00:30:12 +08:00
|
|
|
lstres = .TRUE.
|
|
|
|
IF ( cell_factor_ <= 0.D0 ) cell_factor_ = 1.2D0
|
|
|
|
IF ( cmass <= 0.D0 ) &
|
2004-01-22 00:40:09 +08:00
|
|
|
CALL errore( 'readin', &
|
2003-12-10 22:57:07 +08:00
|
|
|
& 'vcsmd: a positive value for cell mass is required', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE
|
|
|
|
cell_factor_ = 1.D0
|
|
|
|
END IF
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL verify_tmpdir()
|
|
|
|
!
|
|
|
|
CALL restart_from_file()
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( startingconfig == 'file' ) CALL read_config_from_file()
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL write_config_to_file()
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
! ... Files
|
|
|
|
!
|
|
|
|
input_drho = ' '
|
|
|
|
output_drho = ' '
|
|
|
|
crystal = ' '
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE iosys
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
2003-12-10 22:57:07 +08:00
|
|
|
SUBROUTINE read_cards( psfile, atomic_positions_ )
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
USE wvfct, ONLY : gamma_only
|
2004-02-26 19:50:36 +08:00
|
|
|
USE brilz, ONLY : at, ibrav, symm_type, celldm
|
2003-10-17 00:30:12 +08:00
|
|
|
USE basis, ONLY : nat, ntyp, ityp, tau, atm
|
|
|
|
USE klist, ONLY : nks
|
|
|
|
USE ktetra, ONLY : nk1_ => nk1, &
|
|
|
|
nk2_ => nk2, &
|
|
|
|
nk3_ => nk3, &
|
|
|
|
k1_ => k1, &
|
|
|
|
k2_ => k2, &
|
|
|
|
k3_ => k3
|
|
|
|
USE klist, ONLY : lxkcry, &
|
|
|
|
xk_ => xk, &
|
|
|
|
wk_ => wk
|
|
|
|
USE fixed_occ, ONLY : tfixed_occ, &
|
|
|
|
f_inp_ => f_inp
|
|
|
|
USE relax, ONLY : fixatom, &
|
|
|
|
if_pos_ => if_pos
|
|
|
|
USE dynam, ONLY : amass
|
2004-02-26 19:50:36 +08:00
|
|
|
USE control_flags, ONLY : lconstrain
|
|
|
|
USE constrains_module, ONLY : nconstr, constr_tol, constr, target
|
2003-12-10 22:57:07 +08:00
|
|
|
USE input_parameters, ONLY : atom_label, atom_pfile, atom_mass, &
|
|
|
|
atom_ptyp, taspc, tapos, rd_pos, &
|
|
|
|
atomic_positions, if_pos, sp_pos, &
|
|
|
|
k_points, xk, wk, nk1, nk2, nk3, &
|
|
|
|
k1, k2, k3, nkstot, cell_symmetry, rd_ht, &
|
2004-02-26 19:50:36 +08:00
|
|
|
trd_ht, f_inp, calculation,&
|
|
|
|
nconstr_inp, constr_tol_inp, constr_inp
|
2003-10-17 00:30:12 +08:00
|
|
|
USE read_cards_module, ONLY : read_cards_base => read_cards
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
USE parser
|
|
|
|
USE basic_algebra_routines, ONLY : norm
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER (LEN=80) :: psfile(ntyp)
|
|
|
|
CHARACTER (LEN=30) :: atomic_positions_
|
|
|
|
!
|
|
|
|
LOGICAL :: tcell = .FALSE.
|
|
|
|
INTEGER :: i, is, ns, ia, ik
|
|
|
|
!
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
amass = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
CALL read_cards_base( 'PW' )
|
|
|
|
!
|
|
|
|
IF ( .NOT. taspc ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' cards ', ' atomic species info missing', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( .NOT. tapos ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' cards ', ' atomic position info missing', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
DO is = 1, ntyp
|
|
|
|
amass(is) = atom_mass(is)
|
|
|
|
psfile(is) = atom_pfile(is)
|
|
|
|
atm(is) = atom_label(is)
|
|
|
|
IF( amass(is) <= 0.D0 ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' iosys ', ' invalid mass ', is )
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
END IF
|
2003-10-17 00:30:12 +08:00
|
|
|
END DO
|
|
|
|
!
|
|
|
|
DO ia = 1, nat
|
2003-12-10 22:57:07 +08:00
|
|
|
tau(:,ia) = rd_pos(:,ia)
|
|
|
|
ityp(ia) = sp_pos(ia)
|
|
|
|
END DO
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... TEMP: calculate fixatom (to be removed)
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
!
|
|
|
|
fixatom = 0
|
|
|
|
fix1: DO ia = nat, 1, -1
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( if_pos(1,ia) /= 0 .OR. &
|
|
|
|
if_pos(2,ia) /= 0 .OR. &
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
if_pos(3,ia) /= 0 ) EXIT fix1
|
|
|
|
fixatom = fixatom + 1
|
|
|
|
END DO fix1
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
! ... The constrain on fixed coordinates is implemented using the array
|
|
|
|
! ... if_pos whose value is 0 when the coordinate is to be kept fixed, 1
|
2004-01-07 00:53:55 +08:00
|
|
|
! ... otherwise. fixatom is maintained for compatibility. ( C.S. 15/10/2003 )
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
if_pos_ = 1
|
|
|
|
if_pos_(:,:) = if_pos(:,1:nat)
|
|
|
|
!
|
|
|
|
atomic_positions_ = TRIM( atomic_positions )
|
|
|
|
!
|
|
|
|
IF ( k_points == 'automatic' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... automatic generation of k-points
|
|
|
|
!
|
|
|
|
gamma_only = .FALSE.
|
|
|
|
lxkcry = .FALSE.
|
|
|
|
nks = 0
|
2004-01-07 00:53:55 +08:00
|
|
|
! nk1,nk2,nk3 and k1,k2,k3 are initialized even when not used
|
2003-12-10 22:57:07 +08:00
|
|
|
nk1_ = nk1
|
|
|
|
nk2_ = nk2
|
|
|
|
nk3_ = nk3
|
|
|
|
k1_ = k1
|
|
|
|
k2_ = k2
|
|
|
|
k3_ = k3
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE IF ( k_points == 'tpiba' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... input k-points are in 2pi/a units
|
|
|
|
!
|
|
|
|
gamma_only = .FALSE.
|
|
|
|
lxkcry = .FALSE.
|
|
|
|
nks = nkstot
|
2003-10-17 00:30:12 +08:00
|
|
|
xk_(:,1:nks) = xk(:,1:nks)
|
|
|
|
wk_(1:nks) = wk(1:nks)
|
2004-01-07 00:53:55 +08:00
|
|
|
nk1_ = 0
|
|
|
|
nk2_ = 0
|
|
|
|
nk3_ = 0
|
|
|
|
k1_ = 0
|
|
|
|
k2_ = 0
|
|
|
|
k3_ = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE IF ( k_points == 'crystal' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... input k-points are in crystal (reciprocal lattice) axis
|
|
|
|
!
|
|
|
|
gamma_only = .FALSE.
|
|
|
|
lxkcry = .TRUE.
|
|
|
|
nks = nkstot
|
2003-10-17 00:30:12 +08:00
|
|
|
xk_(:,1:nks) = xk(:,1:nks)
|
|
|
|
wk_(1:nks) = wk(1:nks)
|
2004-01-07 00:53:55 +08:00
|
|
|
nk1_ = 0
|
|
|
|
nk2_ = 0
|
|
|
|
nk3_ = 0
|
|
|
|
k1_ = 0
|
|
|
|
k2_ = 0
|
|
|
|
k3_ = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE IF ( k_points == 'gamma' ) THEN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... Only Gamma (k=0) is used
|
|
|
|
! ... specialized routines are used
|
|
|
|
!
|
|
|
|
gamma_only = .TRUE.
|
|
|
|
lxkcry = .FALSE.
|
|
|
|
nks = 1
|
|
|
|
xk_(:,1) = 0.0
|
|
|
|
wk_(1) = 1.0
|
2004-01-07 00:53:55 +08:00
|
|
|
nk1_ = 0
|
|
|
|
nk2_ = 0
|
|
|
|
nk3_ = 0
|
|
|
|
k1_ = 0
|
|
|
|
k2_ = 0
|
|
|
|
k3_ = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... default: input k-points are in 2pi/a units
|
|
|
|
!
|
|
|
|
gamma_only = .FALSE.
|
|
|
|
lxkcry = .FALSE.
|
|
|
|
nks = nkstot
|
2003-10-17 00:30:12 +08:00
|
|
|
xk_(:,1:nks) = xk(:,1:nks)
|
|
|
|
wk_(1:nks) = wk(1:nks)
|
2004-01-07 00:53:55 +08:00
|
|
|
nk1_ = 0
|
|
|
|
nk2_ = 0
|
|
|
|
nk3_ = 0
|
|
|
|
k1_ = 0
|
|
|
|
k2_ = 0
|
|
|
|
k3_ = 0
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
!
|
|
|
|
IF ( tfixed_occ ) THEN
|
|
|
|
IF ( nks > 1 .OR. ( nk1 * nk2 * nk3 ) > 1 ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( 'read_cards', &
|
|
|
|
& 'only one k point with fixed occupations', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
f_inp_ = f_inp
|
|
|
|
ENDIF
|
|
|
|
!
|
|
|
|
IF ( trd_ht ) THEN
|
All namelists and cards moved to Modules/input_parameters.f90 .
From now on, all new input variables should be added
to this module, and then copied to the code internal
variables in the input.f90 subroutine
The namelists and cards parsers are in :
Modules/read_namelists.f90 and Modules/read_cards.f90
files input_parameters.f90 read_namelists.f90 read_cards.f90
are shared by all codes, while each code has its own version
of input.f90 ( used to copy input values into internals variables ).
EXAMPLE:
suppose you need to add a new input variable called "pippo"
to the namelist control, then:
1) add pippo to the input_parameters.f90 file containing the
namelist control
INTEGER :: pippo = 0
NAMELIST / control / ....., pippo
remember: always set an initialization value!
2) add pippo to the control_default subroutine
( cantained in module read_namelists.f90 )
subroutine control_default( prog )
...
IF( prog == 'PW' ) pippo = 10
...
end subroutine
this routine set the default value for pippo,
that could vary with the code
3) add pippo to the control_bcast subroutine
( cantained in module read_namelists.f90 )
subroutine control_bcast( )
...
call mp_bcast( pippo )
...
end subroutine
4) add pippo to the control_checkin subroutine
( cantained in module read_namelists.f90 )
subroutine control_checking( prog )
...
IF( pippo < 0 ) &
CALL error(' control_checkin ',' variable pippo less than 0 ', 1 )
...
end subroutine
5) Copy the value of pippo in the code internal variables
( file input.f90 )
subroutine iosys()
use input_parameters, only: ...., pippo
use pwcom, only: ....., myvar
...
call read_namelists( 'PW' )
...
myvar = pippo
...
end subroutine
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@282 c92efa57-630b-4861-b058-cf58834340f0
2003-07-31 21:24:20 +08:00
|
|
|
symm_type = cell_symmetry
|
2003-10-17 00:30:12 +08:00
|
|
|
at = TRANSPOSE( rd_ht )
|
|
|
|
tcell = .TRUE.
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
IF ( ibrav == 0 .AND. .NOT. tcell ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' cards ', ' ibrav=0: must read cell parameters', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
IF ( ibrav /= 0 .AND. tcell ) &
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( ' cards ', ' redundant data for cell parameters', 2 )
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
IF ( lconstrain ) THEN
|
|
|
|
!
|
|
|
|
nconstr = nconstr_inp
|
|
|
|
constr_tol = constr_tol_inp
|
|
|
|
!
|
|
|
|
ALLOCATE( constr(2,nconstr) )
|
|
|
|
ALLOCATE( target(nconstr) )
|
|
|
|
!
|
|
|
|
constr(:,:) = constr_inp(:,1:nconstr)
|
|
|
|
!
|
|
|
|
! ... target value of the constrain ( in bohr )
|
|
|
|
!
|
|
|
|
DO ia = 1, nconstr
|
|
|
|
!
|
|
|
|
target(ia) = norm( tau(:,constr(1,ia)) - &
|
|
|
|
tau(:,constr(2,ia)) ) * celldm(1)
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
|
|
|
! ... "if_pos" constrains are not compatible with other constrains
|
|
|
|
!
|
|
|
|
if_pos_ = 1
|
|
|
|
!
|
|
|
|
IF ( fixatom > 0 ) &
|
|
|
|
CALL errore( ' cards ', &
|
|
|
|
& ' constrains are not compatible with fixed atoms', 1 )
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
RETURN
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
END SUBROUTINE read_cards
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
2003-12-10 22:57:07 +08:00
|
|
|
SUBROUTINE verify_tmpdir()
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
USE input_parameters, ONLY : restart_mode
|
2004-01-23 23:08:03 +08:00
|
|
|
USE control_flags, ONLY : lneb
|
2003-12-10 22:57:07 +08:00
|
|
|
USE io_files, ONLY : prefix, tmp_dir, nd_nmbr
|
|
|
|
USE neb_variables, ONLY : num_of_images
|
|
|
|
#if defined (__PARA)
|
|
|
|
USE para, ONLY : me, mypool
|
|
|
|
USE mp, ONLY : mp_barrier
|
|
|
|
#endif
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
USE parser, ONLY : int_to_char, delete_if_present
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
INTEGER :: l, ios, image
|
|
|
|
CHARACTER (LEN=80) :: tmp_dir_saved
|
|
|
|
INTEGER :: c_mkdir
|
|
|
|
EXTERNAL c_mkdir
|
2003-10-17 00:30:12 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
! ... verify if tmp_dir ends with /, add one if needed
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
l = LEN_TRIM( tmp_dir )
|
|
|
|
IF ( tmp_dir(l:l) /= '/' ) THEN
|
|
|
|
IF ( l > 0 .AND. l < LEN( tmp_dir ) ) THEN
|
2003-01-20 05:58:50 +08:00
|
|
|
tmp_dir(l+1:l+1)='/'
|
2003-10-17 00:30:12 +08:00
|
|
|
ELSE
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL errore( 'outdir: ', tmp_dir // ' truncated or empty', 1 )
|
2003-10-17 00:30:12 +08:00
|
|
|
END IF
|
|
|
|
END IF
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
ios = 0
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
OPEN( UNIT = 4, FILE = TRIM( tmp_dir ) // 'pwscf' // nd_nmbr, &
|
2003-10-17 00:30:12 +08:00
|
|
|
STATUS = 'UNKNOWN', FORM = 'UNFORMATTED', IOSTAT = ios )
|
|
|
|
CLOSE( UNIT = 4, STATUS = 'DELETE' )
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( ios /= 0 ) CALL errore( 'outdir: ', TRIM( tmp_dir ) // &
|
|
|
|
& ' non existent or non writable', 1 )
|
|
|
|
!
|
|
|
|
! ... if starting from scratch all temporary files are removed
|
|
|
|
! ... from tmp_dir ( only by the master node )
|
|
|
|
!
|
|
|
|
IF ( restart_mode == 'from_scratch' ) THEN
|
|
|
|
!
|
|
|
|
#if defined (__PARA)
|
|
|
|
IF ( me == 1 .AND. mypool == 1 ) THEN
|
2004-02-26 19:50:36 +08:00
|
|
|
#endif
|
|
|
|
!
|
|
|
|
! ... wfc-extrapolation file is removed
|
|
|
|
!
|
|
|
|
CALL delete_if_present( TRIM( tmp_dir ) // TRIM( prefix ) // '.update' )
|
|
|
|
!
|
|
|
|
! ... MD restart file is removed
|
|
|
|
!
|
|
|
|
CALL delete_if_present( TRIM( tmp_dir ) // TRIM( prefix ) // '.md' )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... BFGS rstart file is removed
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
CALL delete_if_present( TRIM( tmp_dir ) // TRIM( prefix ) // '.bfgs' )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
#if defined (__PARA)
|
|
|
|
END IF
|
2004-01-06 02:11:01 +08:00
|
|
|
#endif
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
! ... NEB specific
|
|
|
|
! ... in the scratch directory the tree of subdirectories needed by neb are
|
|
|
|
! ... created
|
|
|
|
!
|
|
|
|
IF ( lneb ) THEN
|
|
|
|
!
|
|
|
|
tmp_dir_saved = tmp_dir
|
|
|
|
!
|
|
|
|
DO image = 1, num_of_images
|
|
|
|
!
|
|
|
|
ios = 0
|
|
|
|
tmp_dir = TRIM( tmp_dir_saved ) // TRIM( prefix ) //"_" // &
|
|
|
|
TRIM( int_to_char( image ) ) // '/'
|
|
|
|
!
|
|
|
|
#if defined (__PARA)
|
|
|
|
IF ( me == 1 .AND. mypool == 1 ) THEN
|
|
|
|
#endif
|
|
|
|
!
|
|
|
|
! ... a scratch directory for this image of the elastic band is
|
|
|
|
! ... created ( only by the master node )
|
|
|
|
!
|
|
|
|
ios = c_mkdir( TRIM( tmp_dir ), LEN_TRIM( tmp_dir ) )
|
|
|
|
!
|
|
|
|
#if defined (__PARA)
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
! ... all jobs are syncronized
|
|
|
|
!
|
|
|
|
CALL mp_barrier()
|
|
|
|
#endif
|
|
|
|
!
|
|
|
|
! ... each job checks whether the scratch directory is accessible
|
|
|
|
! ... or not
|
|
|
|
!
|
|
|
|
OPEN( UNIT = 4, FILE = TRIM( tmp_dir ) // 'pwscf' // nd_nmbr, &
|
|
|
|
STATUS = 'UNKNOWN', FORM = 'UNFORMATTED', IOSTAT = ios )
|
|
|
|
CLOSE( UNIT = 4, STATUS = 'DELETE' )
|
|
|
|
!
|
|
|
|
IF ( ios /= 0 ) &
|
|
|
|
CALL errore( 'outdir: ', TRIM( tmp_dir ) // &
|
|
|
|
& ' non existent or non writable', 1 )
|
|
|
|
!
|
|
|
|
! ... if starting from scratch all temporary files are removed
|
|
|
|
! ... from tmp_dir ( only by the master node )
|
|
|
|
!
|
|
|
|
IF ( restart_mode == 'from_scratch' ) THEN
|
|
|
|
!
|
|
|
|
#if defined (__PARA)
|
|
|
|
IF ( me == 1 .AND. mypool == 1 ) THEN
|
|
|
|
#endif
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
! ... standard output of the self-consistency is removed
|
2004-01-06 02:11:01 +08:00
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
CALL delete_if_present( TRIM( tmp_dir ) // 'PW.out' )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
#if defined (__PARA)
|
|
|
|
END IF
|
2004-01-06 02:11:01 +08:00
|
|
|
#endif
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
|
|
|
tmp_dir = tmp_dir_saved
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
2003-10-17 00:30:12 +08:00
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE verify_tmpdir
|