2003-07-07 05:47:17 +08:00
|
|
|
!
|
2011-12-24 01:23:56 +08:00
|
|
|
! Copyright (C) 2002-2011 Quantum ESPRESSO group
|
2003-07-07 05:47:17 +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 .
|
|
|
|
!
|
2004-11-03 17:53:12 +08:00
|
|
|
!----------------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
MODULE read_namelists_module
|
2004-11-03 17:53:12 +08:00
|
|
|
!----------------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... this module handles the reading of input namelists
|
|
|
|
! ... written by: Carlo Cavazzoni
|
|
|
|
! --------------------------------------------------
|
|
|
|
!
|
2005-09-17 10:14:39 +08:00
|
|
|
USE kinds, ONLY : DP
|
2003-11-17 19:52:06 +08:00
|
|
|
USE input_parameters
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
2005-09-20 23:17:18 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
SAVE
|
|
|
|
!
|
|
|
|
PRIVATE
|
2005-09-20 23:17:18 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
REAL(DP), PARAMETER :: sm_not_set = -20.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2005-09-17 10:14:39 +08:00
|
|
|
PUBLIC :: read_namelists, sm_not_set
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2010-09-01 22:24:41 +08:00
|
|
|
! ... modules needed by read_xml.f90
|
|
|
|
!
|
|
|
|
PUBLIC :: control_defaults, system_defaults, ee_defaults, &
|
|
|
|
electrons_defaults, wannier_ac_defaults, ions_defaults, &
|
|
|
|
cell_defaults, press_ai_defaults, wannier_defaults, control_bcast, &
|
|
|
|
system_bcast, ee_bcast, electrons_bcast, ions_bcast, cell_bcast, &
|
|
|
|
press_ai_bcast, wannier_bcast, wannier_ac_bcast, control_checkin, &
|
|
|
|
system_checkin, electrons_checkin, ions_checkin, cell_checkin, &
|
2010-12-21 06:14:44 +08:00
|
|
|
wannier_checkin, wannier_ac_checkin, fixval
|
2010-09-01 22:24:41 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! ... end of module-scope declarations
|
|
|
|
!
|
|
|
|
! ----------------------------------------------
|
|
|
|
!
|
|
|
|
CONTAINS
|
|
|
|
!
|
|
|
|
!=-----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist CONTROL
|
|
|
|
!
|
|
|
|
!=-----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE control_defaults( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
!
|
|
|
|
IF ( prog == 'PW' ) THEN
|
2007-03-12 21:50:45 +08:00
|
|
|
title = ' '
|
|
|
|
calculation = 'scf'
|
2003-11-17 19:52:06 +08:00
|
|
|
ELSE
|
2007-03-12 21:50:45 +08:00
|
|
|
title = 'MD Simulation'
|
|
|
|
calculation = 'cp'
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
2010-11-02 05:08:30 +08:00
|
|
|
|
2003-11-17 19:52:06 +08:00
|
|
|
verbosity = 'default'
|
2007-03-12 21:50:45 +08:00
|
|
|
IF( prog == 'PW' ) restart_mode = 'from_scratch'
|
|
|
|
IF( prog == 'CP' ) restart_mode = 'restart'
|
2003-11-17 19:52:06 +08:00
|
|
|
nstep = 50
|
|
|
|
IF( prog == 'PW' ) iprint = 100000
|
|
|
|
IF( prog == 'CP' ) iprint = 10
|
2004-07-20 01:19:16 +08:00
|
|
|
IF( prog == 'PW' ) isave = 0
|
|
|
|
IF( prog == 'CP' ) isave = 100
|
2004-12-09 15:31:40 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
tstress = .FALSE.
|
|
|
|
tprnfor = .FALSE.
|
2006-11-21 04:11:43 +08:00
|
|
|
tabps = .FALSE.
|
2004-12-09 15:31:40 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( prog == 'PW' ) dt = 20.0_DP
|
|
|
|
IF( prog == 'CP' ) dt = 1.0_DP
|
2004-12-09 15:31:40 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
ndr = 50
|
|
|
|
ndw = 50
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
! ... use the path specified as outdir and the filename prefix
|
2006-08-09 02:05:16 +08:00
|
|
|
! ... to store output data
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2006-10-11 22:42:42 +08:00
|
|
|
CALL get_env( 'ESPRESSO_TMPDIR', outdir )
|
2006-08-09 02:05:16 +08:00
|
|
|
IF ( TRIM( outdir ) == ' ' ) outdir = './'
|
2007-03-12 21:50:45 +08:00
|
|
|
IF( prog == 'PW' ) prefix = 'pwscf'
|
|
|
|
IF( prog == 'CP' ) prefix = 'cp'
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... directory containing the pseudopotentials
|
|
|
|
!
|
2006-12-11 18:19:53 +08:00
|
|
|
CALL get_env( 'ESPRESSO_PSEUDO', pseudo_dir )
|
2006-08-09 02:05:16 +08:00
|
|
|
IF ( TRIM( pseudo_dir ) == ' ') THEN
|
2006-12-11 18:19:53 +08:00
|
|
|
CALL get_env( 'HOME', pseudo_dir )
|
2006-08-09 02:05:16 +08:00
|
|
|
pseudo_dir = TRIM( pseudo_dir ) // '/espresso/pseudo/'
|
|
|
|
END IF
|
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
refg = 0.05_DP
|
|
|
|
max_seconds = 1.E+7_DP
|
|
|
|
ekin_conv_thr = 1.E-6_DP
|
|
|
|
etot_conv_thr = 1.E-4_DP
|
|
|
|
forc_conv_thr = 1.E-3_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
disk_io = 'default'
|
|
|
|
dipfield = .FALSE.
|
|
|
|
lberry = .FALSE.
|
|
|
|
gdir = 0
|
|
|
|
nppstr = 0
|
|
|
|
wf_collect = .FALSE.
|
2010-06-05 04:37:38 +08:00
|
|
|
IF( prog == 'CP' ) wf_collect = .TRUE. ! default for CP is true
|
2004-12-21 23:48:19 +08:00
|
|
|
printwfc = -1
|
2005-08-16 20:56:49 +08:00
|
|
|
lelfield = .FALSE.
|
2005-10-26 07:11:53 +08:00
|
|
|
nberrycyc = 1
|
2007-02-19 17:21:42 +08:00
|
|
|
lkpoint_dir = .TRUE.
|
2010-08-20 21:10:13 +08:00
|
|
|
lecrpa = .FALSE.
|
2006-06-07 10:01:57 +08:00
|
|
|
!
|
2006-08-11 16:29:52 +08:00
|
|
|
saverho = .TRUE.
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist SYSTEM
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE system_defaults( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
ibrav = -1
|
2007-06-12 01:13:15 +08:00
|
|
|
celldm = (/ 0.0_DP, 0.0_DP, 0.0_DP, 0.0_DP, 0.0_DP, 0.0_DP /)
|
|
|
|
a = 0.0_DP
|
|
|
|
b = 0.0_DP
|
|
|
|
c = 0.0_DP
|
|
|
|
cosab = 0.0_DP
|
|
|
|
cosac = 0.0_DP
|
|
|
|
cosbc = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
nat = 0
|
|
|
|
ntyp = 0
|
|
|
|
nbnd = 0
|
2007-06-12 01:13:15 +08:00
|
|
|
tot_charge = 0.0_DP
|
2006-01-02 20:07:18 +08:00
|
|
|
tot_magnetization = -1
|
2007-06-12 01:13:15 +08:00
|
|
|
ecutwfc = 0.0_DP
|
|
|
|
ecutrho = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
nr1 = 0
|
|
|
|
nr2 = 0
|
|
|
|
nr3 = 0
|
|
|
|
nr1s = 0
|
|
|
|
nr2s = 0
|
|
|
|
nr3s = 0
|
2007-07-18 18:23:06 +08:00
|
|
|
nr1b = 0
|
|
|
|
nr2b = 0
|
|
|
|
nr3b = 0
|
2003-11-17 19:52:06 +08:00
|
|
|
occupations = 'fixed'
|
|
|
|
smearing = 'gaussian'
|
2007-06-12 01:13:15 +08:00
|
|
|
degauss = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
nspin = 1
|
|
|
|
nosym = .FALSE.
|
2008-09-16 21:53:01 +08:00
|
|
|
nosym_evc = .FALSE.
|
2009-06-19 21:40:43 +08:00
|
|
|
force_symmorphic = .FALSE.
|
2011-08-15 22:30:00 +08:00
|
|
|
use_all_frac = .FALSE.
|
2008-03-11 23:38:23 +08:00
|
|
|
noinv = .FALSE.
|
2007-06-12 01:13:15 +08:00
|
|
|
ecfixed = 0.0_DP
|
|
|
|
qcutz = 0.0_DP
|
|
|
|
q2sigma = 0.01_DP
|
2005-11-02 23:42:06 +08:00
|
|
|
input_dft = 'none'
|
2012-04-10 20:38:48 +08:00
|
|
|
ecutfock = -1.0_DP
|
2005-03-24 01:20:26 +08:00
|
|
|
!
|
|
|
|
! ... set starting_magnetization to an invalid value:
|
|
|
|
! ... in PW starting_magnetization MUST be set for at least one atomic type
|
2006-10-26 00:27:50 +08:00
|
|
|
! ... (unless the magnetization is set in other ways)
|
2007-03-12 21:50:45 +08:00
|
|
|
! ... in CP starting_magnetization MUST REMAIN UNSET
|
2005-03-24 01:20:26 +08:00
|
|
|
!
|
2005-09-17 10:14:39 +08:00
|
|
|
starting_magnetization = sm_not_set
|
2005-03-24 01:20:26 +08:00
|
|
|
|
2004-02-14 16:39:34 +08:00
|
|
|
IF ( prog == 'PW' ) THEN
|
2004-07-05 14:50:22 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
starting_ns_eigenvalue = -1.0_DP
|
2004-07-05 14:50:22 +08:00
|
|
|
U_projection_type = 'atomic'
|
|
|
|
!
|
|
|
|
END IF
|
2003-11-17 19:52:06 +08:00
|
|
|
lda_plus_U = .FALSE.
|
2012-02-20 19:01:51 +08:00
|
|
|
lda_plus_u_kind = 0
|
2007-06-12 01:13:15 +08:00
|
|
|
Hubbard_U = 0.0_DP
|
2012-02-20 19:01:51 +08:00
|
|
|
Hubbard_J = 0.0_DP
|
2007-06-12 01:13:15 +08:00
|
|
|
Hubbard_alpha = 0.0_DP
|
2012-02-20 19:01:51 +08:00
|
|
|
step_pen=.false.
|
|
|
|
A_pen=0.0_DP
|
|
|
|
sigma_pen=0.01_DP
|
|
|
|
alpha_pen=0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
edir = 1
|
2007-06-12 01:13:15 +08:00
|
|
|
emaxpos = 0.5_DP
|
|
|
|
eopreg = 0.1_DP
|
2009-06-28 04:45:51 +08:00
|
|
|
eamp = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2005-12-28 18:42:31 +08:00
|
|
|
! ... postprocessing of DOS & phonons & el-ph
|
|
|
|
la2F = .FALSE.
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! ... non collinear program variables
|
|
|
|
!
|
2004-05-03 22:23:10 +08:00
|
|
|
lspinorb = .FALSE.
|
2011-03-16 00:56:08 +08:00
|
|
|
starting_spin_angle=.FALSE.
|
2003-11-17 19:52:06 +08:00
|
|
|
noncolin = .FALSE.
|
2007-06-12 01:13:15 +08:00
|
|
|
lambda = 1.0_DP
|
2005-03-24 01:20:26 +08:00
|
|
|
constrained_magnetization= 'none'
|
2007-06-12 01:13:15 +08:00
|
|
|
fixed_magnetization = 0.0_DP
|
|
|
|
B_field = 0.0_DP
|
|
|
|
angle1 = 0.0_DP
|
|
|
|
angle2 = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
report = 1
|
2006-12-11 22:51:54 +08:00
|
|
|
!
|
2010-08-19 17:33:14 +08:00
|
|
|
no_t_rev = .FALSE.
|
|
|
|
!
|
2009-09-16 04:29:07 +08:00
|
|
|
assume_isolated = 'none'
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2009-07-30 23:47:47 +08:00
|
|
|
one_atom_occupations=.FALSE.
|
|
|
|
!
|
2007-01-25 19:42:52 +08:00
|
|
|
spline_ps = .false.
|
2008-06-11 18:47:40 +08:00
|
|
|
!
|
2009-04-03 00:05:09 +08:00
|
|
|
real_space = .false.
|
2009-06-15 03:05:55 +08:00
|
|
|
!
|
|
|
|
! ... DFT-D
|
|
|
|
!
|
|
|
|
london = .false.
|
|
|
|
london_s6 = 0.75_DP
|
|
|
|
london_rcut = 200.00_DP
|
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
#ifdef __ENVIRON
|
|
|
|
! ... Environ
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
do_environ = .false.
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
|
|
|
#endif
|
2011-04-16 03:17:09 +08:00
|
|
|
! ... ESM
|
|
|
|
!
|
|
|
|
esm_bc='pbc'
|
|
|
|
esm_efield=0.0_DP
|
|
|
|
esm_w=0.0_DP
|
|
|
|
esm_nfit=4
|
|
|
|
esm_debug=.FALSE.
|
|
|
|
esm_debug_gpmax=0
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2011-04-22 00:12:36 +08:00
|
|
|
|
2012-02-15 00:18:50 +08:00
|
|
|
#ifdef __ENVIRON
|
2011-04-22 00:12:36 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
! Variables initialization for Namelist ENVIRON
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
2012-02-15 00:18:50 +08:00
|
|
|
SUBROUTINE environ_defaults( prog )
|
2011-04-22 00:12:36 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
!
|
2012-02-21 23:36:54 +08:00
|
|
|
verbose = 0
|
|
|
|
environ_thr = 1.D-1
|
|
|
|
environ_type = 'input'
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
|
|
|
stype = 1
|
2012-01-20 00:10:44 +08:00
|
|
|
rhomax = 0.005
|
2011-04-22 00:12:36 +08:00
|
|
|
rhomin = 0.0001
|
|
|
|
tbeta = 4.8
|
|
|
|
!
|
2012-02-21 23:36:54 +08:00
|
|
|
env_static_permittivity = 1.D0
|
2011-04-22 00:12:36 +08:00
|
|
|
eps_mode = 'electronic'
|
|
|
|
solvationrad(:) = 3.D0
|
|
|
|
atomicspread(:) = 0.5D0
|
2012-05-01 02:21:54 +08:00
|
|
|
add_jellium = .false.
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2011-09-27 23:46:58 +08:00
|
|
|
ifdtype = 1
|
2012-01-20 00:10:44 +08:00
|
|
|
nfdpoint = 2
|
2011-09-27 23:46:58 +08:00
|
|
|
!
|
2012-04-16 21:24:39 +08:00
|
|
|
mixtype = 'linear'
|
|
|
|
ndiis = 1
|
2011-04-22 00:12:36 +08:00
|
|
|
mixrhopol = 0.5
|
|
|
|
tolrhopol = 1.D-10
|
|
|
|
!
|
2012-02-21 23:36:54 +08:00
|
|
|
env_surface_tension = 0.D0
|
2011-04-22 00:12:36 +08:00
|
|
|
delta = 0.00001D0
|
|
|
|
!
|
2012-02-21 23:36:54 +08:00
|
|
|
env_pressure = 0.D0
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2012-04-16 21:24:39 +08:00
|
|
|
cion = 0.0D0
|
|
|
|
zion = 1.0D0
|
|
|
|
rhopb = 0.0001D0
|
|
|
|
solvent_temperature = 300.0D0
|
|
|
|
!
|
2011-04-22 00:12:36 +08:00
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
!
|
|
|
|
#endif
|
2008-06-11 18:47:40 +08:00
|
|
|
! DCC
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist EE
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
SUBROUTINE ee_defaults( prog )
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
!
|
|
|
|
ncompx = 1
|
|
|
|
ncompy = 1
|
|
|
|
ncompz = 1
|
|
|
|
mr1 = 0
|
|
|
|
mr2 = 0
|
|
|
|
mr3 = 0
|
|
|
|
ecutcoarse = 100.D0
|
|
|
|
errtol = 1.d-22
|
|
|
|
nlev = 2
|
|
|
|
itmax = 1000
|
|
|
|
whichbc = 0
|
|
|
|
! centercompx = 0.D0
|
|
|
|
! centercompy = 0.D0
|
|
|
|
! centercompz = 0.D0
|
|
|
|
! spreadcomp = -9999.D0
|
|
|
|
mixing_charge_compensation = 1.0D0
|
|
|
|
n_charge_compensation = 5
|
|
|
|
comp_thr = 1.D-4
|
|
|
|
! multipole = 'dipole'
|
|
|
|
! poisson_maxiter = 5000
|
|
|
|
! poisson_thr = 1.D-6
|
|
|
|
! comp_thr = 1.D-2
|
|
|
|
! ebc_thr = 1.D-2
|
|
|
|
! rhoionmax = 1.D0
|
|
|
|
! smoothspr = 0.25D0
|
|
|
|
! deltapot = 5.D-1
|
|
|
|
nlev = 2
|
|
|
|
! which_smoothing = 'sphere'
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist ELECTRONS
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE electrons_defaults( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
emass = 400.0_DP
|
|
|
|
emass_cutoff = 2.5_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
orthogonalization = 'ortho'
|
2007-06-12 01:13:15 +08:00
|
|
|
ortho_eps = 1.E-8_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
ortho_max = 20
|
|
|
|
electron_maxstep = 100
|
|
|
|
!
|
|
|
|
! ... ( 'sd' | 'cg' | 'damp' | 'verlet' | 'none' | 'diis' )
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
electron_dynamics = 'none'
|
2007-06-12 01:13:15 +08:00
|
|
|
electron_damping = 0.1_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... ( 'zero' | 'default' )
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
electron_velocities = 'default'
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... ( 'nose' | 'not_controlled' | 'rescaling')
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
electron_temperature = 'not_controlled'
|
2007-06-12 01:13:15 +08:00
|
|
|
ekincw = 0.001_DP
|
|
|
|
fnosee = 1.0_DP
|
|
|
|
ampre = 0.0_DP
|
|
|
|
grease = 1.0_DP
|
|
|
|
conv_thr = 1.E-6_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
diis_size = 4
|
|
|
|
diis_nreset = 3
|
2007-06-12 01:13:15 +08:00
|
|
|
diis_hcut = 1.0_DP
|
|
|
|
diis_wthr = 0.0_DP
|
|
|
|
diis_delt = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
diis_maxstep = 100
|
|
|
|
diis_rot = .FALSE.
|
2007-06-12 01:13:15 +08:00
|
|
|
diis_fthr = 0.0_DP
|
|
|
|
diis_temp = 0.0_DP
|
|
|
|
diis_achmix = 0.0_DP
|
|
|
|
diis_g0chmix = 0.0_DP
|
|
|
|
diis_g1chmix = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
diis_nchmix = 3
|
|
|
|
diis_nrot = 3
|
2007-06-12 01:13:15 +08:00
|
|
|
diis_rothr = 0.0_DP
|
|
|
|
diis_ethr = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
diis_chguess = .FALSE.
|
|
|
|
mixing_mode = 'plain'
|
|
|
|
mixing_fixed_ns = 0
|
2007-06-12 01:13:15 +08:00
|
|
|
mixing_beta = 0.7_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
mixing_ndim = 8
|
2008-11-22 01:36:45 +08:00
|
|
|
diagonalization = 'david'
|
2007-06-12 01:13:15 +08:00
|
|
|
diago_thr_init = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
diago_cg_maxiter = 20
|
|
|
|
diago_david_ndim = 4
|
2006-04-25 01:32:08 +08:00
|
|
|
diago_full_acc = .FALSE.
|
2004-11-03 17:53:12 +08:00
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
sic = 'none'
|
2007-06-12 01:13:15 +08:00
|
|
|
sic_epsilon = 0.0_DP
|
|
|
|
sic_alpha = 0.0_DP
|
2004-04-02 00:30:59 +08:00
|
|
|
force_pairing = .false.
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
fermi_energy = 0.0_DP
|
2005-10-22 06:57:21 +08:00
|
|
|
n_inner = 2
|
2007-03-22 07:45:24 +08:00
|
|
|
niter_cold_restart=1
|
2007-06-12 01:13:15 +08:00
|
|
|
lambda_cold=0.03_DP
|
2004-12-21 23:48:19 +08:00
|
|
|
rotation_dynamics = "line-minimization"
|
|
|
|
occupation_dynamics = "line-minimization"
|
2007-06-12 01:13:15 +08:00
|
|
|
rotmass = 0.0_DP
|
|
|
|
occmass = 0.0_DP
|
|
|
|
rotation_damping = 0.0_DP
|
|
|
|
occupation_damping = 0.0_DP
|
2005-05-27 06:42:05 +08:00
|
|
|
!
|
|
|
|
tcg = .FALSE.
|
2008-07-29 09:46:21 +08:00
|
|
|
maxiter = 100
|
2007-06-12 01:13:15 +08:00
|
|
|
passop = 0.3_DP
|
2007-03-22 07:45:24 +08:00
|
|
|
niter_cg_restart = 20
|
2007-06-12 01:13:15 +08:00
|
|
|
etresh = 1.E-6_DP
|
2005-05-27 06:42:05 +08:00
|
|
|
!
|
|
|
|
epol = 3
|
2007-06-12 01:13:15 +08:00
|
|
|
efield = 0.0_DP
|
2006-02-02 01:57:32 +08:00
|
|
|
epol2 = 3
|
2007-06-12 01:13:15 +08:00
|
|
|
efield2 = 0.0_DP
|
2008-06-10 00:43:53 +08:00
|
|
|
efield_cart(1)=0.d0
|
|
|
|
efield_cart(2)=0.d0
|
|
|
|
efield_cart(3)=0.d0
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2006-12-16 01:58:58 +08:00
|
|
|
occupation_constraints = .false.
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2011-03-11 21:45:04 +08:00
|
|
|
adaptive_thr = .false.
|
|
|
|
conv_thr_init = 0.1E-2_DP
|
|
|
|
conv_thr_multi = 0.1_DP
|
2012-02-02 05:07:04 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2009-02-12 16:07:11 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist WANNIER_AC
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
SUBROUTINE wannier_ac_defaults( prog )
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
!
|
|
|
|
plot_wannier = .FALSE.
|
|
|
|
use_energy_int = .FALSE.
|
|
|
|
print_wannier_coeff = .FALSE.
|
|
|
|
nwan = 0
|
|
|
|
constrain_pot = 0.d0
|
|
|
|
plot_wan_num = 0
|
|
|
|
plot_wan_spin = 1
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
|
2003-11-17 19:52:06 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist IONS
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE ions_defaults( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
! ... ( 'full' | 'coarse-grained' )
|
|
|
|
!
|
|
|
|
phase_space = 'full'
|
|
|
|
!
|
2006-01-29 06:35:48 +08:00
|
|
|
! ... ( 'sd' | 'cg' | 'damp' | 'verlet' | 'none' | 'bfgs' | 'beeman' )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2004-03-20 23:43:47 +08:00
|
|
|
ion_dynamics = 'none'
|
2007-06-12 01:13:15 +08:00
|
|
|
ion_radius = 0.5_DP
|
|
|
|
ion_damping = 0.1_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... ( 'default' | 'from_input' )
|
|
|
|
!
|
|
|
|
ion_positions = 'default'
|
|
|
|
!
|
|
|
|
! ... ( 'zero' | 'default' | 'from_input' )
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
ion_velocities = 'default'
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2007-03-01 02:26:11 +08:00
|
|
|
! ... ( 'nose' | 'not_controlled' | 'rescaling' | 'berendsen' |
|
|
|
|
! 'andersen' | 'langevin' )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2006-01-29 06:35:48 +08:00
|
|
|
ion_temperature = 'not_controlled'
|
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
tempw = 300.0_DP
|
|
|
|
fnosep = -1.0_DP
|
|
|
|
fnosep(1) = 1.0_DP
|
2006-01-29 06:35:48 +08:00
|
|
|
nhpcl = 0
|
|
|
|
nhptyp = 0
|
|
|
|
ndega = 0
|
|
|
|
tranp = .FALSE.
|
2007-06-12 01:13:15 +08:00
|
|
|
amprp = 0.0_DP
|
|
|
|
greasp = 1.0_DP
|
|
|
|
tolp = 100.0_DP
|
2006-01-29 06:35:48 +08:00
|
|
|
ion_nstepe = 1
|
|
|
|
ion_maxstep = 100
|
2007-06-12 01:13:15 +08:00
|
|
|
delta_t = 1.0_DP
|
2007-03-05 18:38:15 +08:00
|
|
|
nraise = 1
|
2006-01-29 06:35:48 +08:00
|
|
|
!
|
2006-06-15 22:27:14 +08:00
|
|
|
refold_pos = .FALSE.
|
|
|
|
remove_rigid_rot = .FALSE.
|
|
|
|
!
|
2011-03-16 17:30:59 +08:00
|
|
|
upscale = 100.0_DP
|
2005-01-10 14:56:14 +08:00
|
|
|
pot_extrapolation = 'atomic'
|
|
|
|
wfc_extrapolation = 'none'
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2004-05-03 16:48:48 +08:00
|
|
|
! ... BFGS defaults
|
|
|
|
!
|
2005-04-21 07:43:15 +08:00
|
|
|
bfgs_ndim = 1
|
2007-06-12 01:13:15 +08:00
|
|
|
trust_radius_max = 0.8_DP ! bohr
|
|
|
|
trust_radius_min = 1.E-4_DP ! bohr
|
|
|
|
trust_radius_ini = 0.5_DP ! bohr
|
|
|
|
w_1 = 0.01_DP
|
|
|
|
w_2 = 0.50_DP
|
2004-05-03 16:48:48 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
sic_rloc = 0.0_DP
|
2004-04-02 00:30:59 +08:00
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
! ... meta-dynamics defaults
|
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
fe_step = 0.4_DP
|
2005-10-16 07:27:47 +08:00
|
|
|
fe_nstep = 100
|
2006-09-13 02:28:57 +08:00
|
|
|
sw_nstep = 10
|
2006-12-11 22:51:54 +08:00
|
|
|
eq_nstep = 0
|
2007-06-12 01:13:15 +08:00
|
|
|
g_amplitude = 0.005_DP
|
2005-09-20 23:17:18 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! Variables initialization for Namelist CELL
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE cell_defaults( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
cell_parameters = 'default'
|
|
|
|
!
|
2007-12-19 22:41:24 +08:00
|
|
|
! ... ( 'sd' | 'pr' | 'none' | 'w' | 'damp-pr' | 'damp-w' | 'bfgs' )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
cell_dynamics = 'none'
|
|
|
|
!
|
|
|
|
! ... ( 'zero' | 'default' )
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
cell_velocities = 'default'
|
2007-06-12 01:13:15 +08:00
|
|
|
press = 0.0_DP
|
|
|
|
wmass = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... ( 'nose' | 'not_controlled' | 'rescaling' )
|
|
|
|
!
|
|
|
|
cell_temperature = 'not_controlled'
|
2007-06-12 01:13:15 +08:00
|
|
|
temph = 0.0_DP
|
|
|
|
fnoseh = 1.0_DP
|
|
|
|
greash = 1.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... ('all'* | 'volume' | 'x' | 'y' | 'z' | 'xy' | 'xz' | 'yz' | 'xyz' )
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
cell_dofree = 'all'
|
2007-06-12 01:13:15 +08:00
|
|
|
cell_factor = 0.0_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
cell_nstepe = 1
|
2007-06-12 01:13:15 +08:00
|
|
|
cell_damping = 0.0_DP
|
|
|
|
press_conv_thr = 0.5_DP
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2006-11-21 04:11:43 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Variables initialization for Namelist PRESS_AI
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
SUBROUTINE press_ai_defaults( prog )
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
abivol = .false.
|
|
|
|
abisur = .false.
|
|
|
|
pvar = .false.
|
|
|
|
fill_vac = .false.
|
|
|
|
cntr = .false.
|
|
|
|
scale_at = .false.
|
|
|
|
t_gauss = .false.
|
|
|
|
jellium = .false.
|
|
|
|
|
2007-06-12 01:13:15 +08:00
|
|
|
P_ext = 0.0_DP
|
|
|
|
P_in = 0.0_DP
|
|
|
|
P_fin = 0.0_DP
|
|
|
|
Surf_t = 0.0_DP
|
|
|
|
rho_thr = 0.0_DP
|
|
|
|
dthr = 0.0_DP
|
|
|
|
step_rad = 0.0_DP
|
|
|
|
delta_eps = 0.0_DP
|
|
|
|
delta_sigma = 0.0_DP
|
|
|
|
R_j = 0.0_DP
|
|
|
|
h_j = 0.0_DP
|
2006-11-21 04:11:43 +08:00
|
|
|
|
|
|
|
n_cntr = 0
|
|
|
|
axis = 3
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
! Variables initialization for Namelist WANNIER
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2005-07-01 22:26:10 +08:00
|
|
|
SUBROUTINE wannier_defaults( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
!
|
|
|
|
!
|
|
|
|
wf_efield = .FALSE.
|
|
|
|
wf_switch = .FALSE.
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
sw_len = 1
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
efx0 = 0.0_DP
|
|
|
|
efy0 = 0.0_DP
|
|
|
|
efz0 = 0.0_DP
|
|
|
|
efx1 = 0.0_DP
|
|
|
|
efy1 = 0.0_DP
|
|
|
|
efz1 = 0.0_DP
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2009-03-04 18:43:48 +08:00
|
|
|
wfsd = 1
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
wfdt = 5.0_DP
|
|
|
|
maxwfdt = 0.30_DP
|
|
|
|
wf_q = 1500.0_DP
|
|
|
|
wf_friction = 0.3_DP
|
2011-08-01 17:57:39 +08:00
|
|
|
!=======================================================================
|
|
|
|
!Lingzhu Kong
|
|
|
|
neigh = 48
|
|
|
|
vnbsp = 0
|
|
|
|
poisson_eps = 1.D-6
|
|
|
|
dis_cutoff = 7.0_DP
|
|
|
|
exx_ps_rcut = 5.0
|
|
|
|
exx_me_rcut = 10.0
|
|
|
|
!=======================================================================
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
nit = 10
|
|
|
|
nsd = 10
|
|
|
|
nsteps = 20
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
tolw = 1.E-8_DP
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
2006-02-14 22:15:28 +08:00
|
|
|
adapt = .TRUE.
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
calwf = 3
|
2005-07-02 04:39:45 +08:00
|
|
|
nwf = 0
|
2005-07-01 22:26:10 +08:00
|
|
|
wffort = 40
|
|
|
|
!
|
|
|
|
writev = .FALSE.
|
2005-07-02 04:39:45 +08:00
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
END SUBROUTINE
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! Broadcast variables values for Namelist CONTROL
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE control_bcast()
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY : ionode_id
|
|
|
|
USE mp, ONLY : mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( title, ionode_id )
|
|
|
|
CALL mp_bcast( calculation, ionode_id )
|
|
|
|
CALL mp_bcast( verbosity, ionode_id )
|
|
|
|
CALL mp_bcast( restart_mode, ionode_id )
|
|
|
|
CALL mp_bcast( nstep, ionode_id )
|
|
|
|
CALL mp_bcast( iprint, ionode_id )
|
|
|
|
CALL mp_bcast( isave, ionode_id )
|
|
|
|
CALL mp_bcast( tstress, ionode_id )
|
|
|
|
CALL mp_bcast( tprnfor, ionode_id )
|
2006-11-21 04:11:43 +08:00
|
|
|
CALL mp_bcast( tabps, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( dt, ionode_id )
|
|
|
|
CALL mp_bcast( ndr, ionode_id )
|
|
|
|
CALL mp_bcast( ndw, ionode_id )
|
|
|
|
CALL mp_bcast( outdir, ionode_id )
|
2006-05-03 22:19:57 +08:00
|
|
|
CALL mp_bcast( wfcdir, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( prefix, ionode_id )
|
|
|
|
CALL mp_bcast( max_seconds, ionode_id )
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL mp_bcast( ekin_conv_thr, ionode_id )
|
|
|
|
CALL mp_bcast( etot_conv_thr, ionode_id )
|
|
|
|
CALL mp_bcast( forc_conv_thr, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( pseudo_dir, ionode_id )
|
|
|
|
CALL mp_bcast( refg, ionode_id )
|
|
|
|
CALL mp_bcast( disk_io, ionode_id )
|
|
|
|
CALL mp_bcast( tefield, ionode_id )
|
2006-05-03 22:19:57 +08:00
|
|
|
CALL mp_bcast( tefield2, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( dipfield, ionode_id )
|
|
|
|
CALL mp_bcast( lberry, ionode_id )
|
|
|
|
CALL mp_bcast( gdir, ionode_id )
|
|
|
|
CALL mp_bcast( nppstr, ionode_id )
|
2007-02-19 17:21:42 +08:00
|
|
|
CALL mp_bcast( lkpoint_dir, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( wf_collect, ionode_id )
|
|
|
|
CALL mp_bcast( printwfc, ionode_id )
|
|
|
|
CALL mp_bcast( lelfield, ionode_id )
|
2005-10-26 07:11:53 +08:00
|
|
|
CALL mp_bcast( nberrycyc, ionode_id )
|
2006-08-11 16:29:52 +08:00
|
|
|
CALL mp_bcast( saverho, ionode_id )
|
2010-08-19 21:24:00 +08:00
|
|
|
CALL mp_bcast( lecrpa, ionode_id )
|
2012-05-20 04:19:27 +08:00
|
|
|
CALL mp_bcast( vdw_table_name,ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Broadcast variables values for Namelist SYSTEM
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE system_bcast()
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY : ionode_id
|
|
|
|
USE mp, ONLY : mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2006-12-11 22:51:54 +08:00
|
|
|
CALL mp_bcast( ibrav, ionode_id )
|
|
|
|
CALL mp_bcast( celldm, ionode_id )
|
|
|
|
CALL mp_bcast( a, ionode_id )
|
|
|
|
CALL mp_bcast( b, ionode_id )
|
|
|
|
CALL mp_bcast( c, ionode_id )
|
|
|
|
CALL mp_bcast( cosab, ionode_id )
|
|
|
|
CALL mp_bcast( cosac, ionode_id )
|
|
|
|
CALL mp_bcast( cosbc, ionode_id )
|
|
|
|
CALL mp_bcast( nat, ionode_id )
|
|
|
|
CALL mp_bcast( ntyp, ionode_id )
|
|
|
|
CALL mp_bcast( nbnd, ionode_id )
|
|
|
|
CALL mp_bcast( tot_charge, ionode_id )
|
|
|
|
CALL mp_bcast( tot_magnetization, ionode_id )
|
|
|
|
CALL mp_bcast( ecutwfc, ionode_id )
|
|
|
|
CALL mp_bcast( ecutrho, ionode_id )
|
|
|
|
CALL mp_bcast( nr1, ionode_id )
|
|
|
|
CALL mp_bcast( nr2, ionode_id )
|
|
|
|
CALL mp_bcast( nr3, ionode_id )
|
|
|
|
CALL mp_bcast( nr1s, ionode_id )
|
|
|
|
CALL mp_bcast( nr2s, ionode_id )
|
|
|
|
CALL mp_bcast( nr3s, ionode_id )
|
|
|
|
CALL mp_bcast( nr1b, ionode_id )
|
|
|
|
CALL mp_bcast( nr2b, ionode_id )
|
|
|
|
CALL mp_bcast( nr3b, ionode_id )
|
|
|
|
CALL mp_bcast( occupations, ionode_id )
|
|
|
|
CALL mp_bcast( smearing, ionode_id )
|
|
|
|
CALL mp_bcast( degauss, ionode_id )
|
|
|
|
CALL mp_bcast( nspin, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
CALL mp_bcast( nosym, ionode_id )
|
2008-09-16 21:53:01 +08:00
|
|
|
CALL mp_bcast( nosym_evc, ionode_id )
|
2008-03-11 23:38:23 +08:00
|
|
|
CALL mp_bcast( noinv, ionode_id )
|
2009-06-19 21:40:43 +08:00
|
|
|
CALL mp_bcast( force_symmorphic, ionode_id )
|
2011-08-15 22:30:00 +08:00
|
|
|
CALL mp_bcast( use_all_frac, ionode_id )
|
2006-12-11 22:51:54 +08:00
|
|
|
CALL mp_bcast( ecfixed, ionode_id )
|
|
|
|
CALL mp_bcast( qcutz, ionode_id )
|
|
|
|
CALL mp_bcast( q2sigma, ionode_id )
|
|
|
|
CALL mp_bcast( input_dft, ionode_id )
|
2009-09-14 00:48:24 +08:00
|
|
|
CALL mp_bcast( nqx1, ionode_id )
|
|
|
|
CALL mp_bcast( nqx2, ionode_id )
|
|
|
|
CALL mp_bcast( nqx3, ionode_id )
|
2010-01-14 21:40:16 +08:00
|
|
|
CALL mp_bcast( exx_fraction, ionode_id )
|
|
|
|
CALL mp_bcast( screening_parameter, ionode_id )
|
2009-09-14 00:48:24 +08:00
|
|
|
CALL mp_bcast( exxdiv_treatment, ionode_id )
|
|
|
|
CALL mp_bcast( x_gamma_extrapolation, ionode_id )
|
|
|
|
CALL mp_bcast( yukawa, ionode_id )
|
|
|
|
CALL mp_bcast( ecutvcut, ionode_id )
|
2011-12-14 02:52:33 +08:00
|
|
|
CALL mp_bcast( ecutfock, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
CALL mp_bcast( starting_magnetization, ionode_id )
|
|
|
|
CALL mp_bcast( starting_ns_eigenvalue, ionode_id )
|
|
|
|
CALL mp_bcast( U_projection_type, ionode_id )
|
|
|
|
CALL mp_bcast( lda_plus_U, ionode_id )
|
2012-02-20 19:01:51 +08:00
|
|
|
CALL mp_bcast( lda_plus_u_kind, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
CALL mp_bcast( Hubbard_U, ionode_id )
|
2012-02-20 19:01:51 +08:00
|
|
|
CALL mp_bcast( Hubbard_J, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
CALL mp_bcast( Hubbard_alpha, ionode_id )
|
2010-05-08 17:58:35 +08:00
|
|
|
CALL mp_bcast( step_pen, ionode_id )
|
|
|
|
CALL mp_bcast( A_pen, ionode_id )
|
|
|
|
CALL mp_bcast( sigma_pen, ionode_id )
|
|
|
|
CALL mp_bcast( alpha_pen, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( edir, ionode_id )
|
|
|
|
CALL mp_bcast( emaxpos, ionode_id )
|
|
|
|
CALL mp_bcast( eopreg, ionode_id )
|
|
|
|
CALL mp_bcast( eamp, ionode_id )
|
2005-12-28 18:42:31 +08:00
|
|
|
CALL mp_bcast( la2F, ionode_id )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! ... non collinear broadcast
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( lspinorb, ionode_id )
|
2011-03-16 00:56:08 +08:00
|
|
|
CALL mp_bcast( starting_spin_angle, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( noncolin, ionode_id )
|
|
|
|
CALL mp_bcast( angle1, ionode_id )
|
|
|
|
CALL mp_bcast( angle2, ionode_id )
|
|
|
|
CALL mp_bcast( report, ionode_id )
|
2005-03-24 01:20:26 +08:00
|
|
|
CALL mp_bcast( constrained_magnetization, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( B_field, ionode_id )
|
|
|
|
CALL mp_bcast( fixed_magnetization, ionode_id )
|
2006-12-11 22:51:54 +08:00
|
|
|
CALL mp_bcast( lambda, ionode_id )
|
|
|
|
!
|
2009-09-14 00:48:24 +08:00
|
|
|
CALL mp_bcast( assume_isolated, ionode_id )
|
|
|
|
CALL mp_bcast( one_atom_occupations, ionode_id )
|
|
|
|
CALL mp_bcast( spline_ps, ionode_id )
|
2006-12-11 22:51:54 +08:00
|
|
|
!
|
2009-09-14 00:48:24 +08:00
|
|
|
CALL mp_bcast( london, ionode_id )
|
|
|
|
CALL mp_bcast( london_s6, ionode_id )
|
|
|
|
CALL mp_bcast( london_rcut, ionode_id )
|
2009-06-15 03:05:55 +08:00
|
|
|
!
|
2010-08-19 17:33:14 +08:00
|
|
|
CALL mp_bcast( no_t_rev, ionode_id )
|
2012-02-15 00:18:50 +08:00
|
|
|
#ifdef __ENVIRON
|
|
|
|
CALL mp_bcast( do_environ, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
#endif
|
2011-04-16 03:17:09 +08:00
|
|
|
!
|
|
|
|
! ... ESM method broadcast
|
|
|
|
!
|
|
|
|
CALL mp_bcast( esm_bc, ionode_id )
|
|
|
|
CALL mp_bcast( esm_efield, ionode_id )
|
|
|
|
CALL mp_bcast( esm_w, ionode_id )
|
|
|
|
CALL mp_bcast( esm_nfit, ionode_id )
|
|
|
|
CALL mp_bcast( esm_debug, ionode_id )
|
|
|
|
CALL mp_bcast( esm_debug_gpmax, ionode_id )
|
2008-06-11 18:47:40 +08:00
|
|
|
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
2012-02-15 00:18:50 +08:00
|
|
|
#ifdef __ENVIRON
|
2011-04-22 00:12:36 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
! Broadcast variables values for Namelist ENVIRON
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
2012-02-15 00:18:50 +08:00
|
|
|
SUBROUTINE environ_bcast()
|
2011-04-22 00:12:36 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
USE io_global, ONLY : ionode_id
|
|
|
|
USE mp, ONLY : mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CALL mp_bcast( verbose, ionode_id )
|
2012-02-15 00:18:50 +08:00
|
|
|
CALL mp_bcast( environ_thr, ionode_id )
|
2012-02-21 23:36:54 +08:00
|
|
|
CALL mp_bcast( environ_type, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
|
|
|
CALL mp_bcast( stype, ionode_id )
|
2012-02-15 00:18:50 +08:00
|
|
|
CALL mp_bcast( rhomax, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
CALL mp_bcast( rhomin, ionode_id )
|
|
|
|
CALL mp_bcast( tbeta, ionode_id )
|
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
CALL mp_bcast( env_static_permittivity, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
CALL mp_bcast( eps_mode, ionode_id )
|
|
|
|
CALL mp_bcast( solvationrad, ionode_id )
|
|
|
|
CALL mp_bcast( atomicspread, ionode_id )
|
2012-05-01 02:21:54 +08:00
|
|
|
CALL mp_bcast( add_jellium, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2011-10-03 23:36:45 +08:00
|
|
|
CALL mp_bcast( ifdtype, ionode_id )
|
|
|
|
CALL mp_bcast( nfdpoint, ionode_id )
|
|
|
|
!
|
2012-04-16 21:24:39 +08:00
|
|
|
CALL mp_bcast( mixtype, ionode_id )
|
|
|
|
CALL mp_bcast( ndiis, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
CALL mp_bcast( mixrhopol, ionode_id )
|
|
|
|
CALL mp_bcast( tolrhopol, ionode_id )
|
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
CALL mp_bcast( env_surface_tension, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
CALL mp_bcast( delta, ionode_id )
|
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
CALL mp_bcast( env_pressure, ionode_id )
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2012-04-16 21:24:39 +08:00
|
|
|
CALL mp_bcast( cion, ionode_id )
|
|
|
|
CALL mp_bcast( zion, ionode_id )
|
|
|
|
CALL mp_bcast( rhopb, ionode_id )
|
|
|
|
CALL mp_bcast( solvent_temperature, ionode_id )
|
|
|
|
!
|
2011-04-22 00:12:36 +08:00
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
!
|
|
|
|
#endif
|
2008-06-11 18:47:40 +08:00
|
|
|
! DCC
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Broadcast variables values for Namelist EE
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
SUBROUTINE ee_bcast()
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
USE io_global, ONLY : ionode_id
|
|
|
|
USE mp, ONLY : mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CALL mp_bcast( ecutcoarse, ionode_id )
|
|
|
|
CALL mp_bcast( mixing_charge_compensation, ionode_id )
|
|
|
|
CALL mp_bcast( errtol, ionode_id )
|
|
|
|
CALL mp_bcast( comp_thr, ionode_id )
|
|
|
|
CALL mp_bcast( nlev, ionode_id )
|
|
|
|
CALL mp_bcast( itmax, ionode_id )
|
|
|
|
CALL mp_bcast( whichbc, ionode_id )
|
|
|
|
CALL mp_bcast( n_charge_compensation, ionode_id )
|
|
|
|
CALL mp_bcast( ncompx, ionode_id )
|
|
|
|
CALL mp_bcast( ncompy, ionode_id )
|
|
|
|
CALL mp_bcast( ncompz, ionode_id )
|
|
|
|
CALL mp_bcast( mr1, ionode_id )
|
|
|
|
CALL mp_bcast( mr2, ionode_id )
|
|
|
|
CALL mp_bcast( mr3, ionode_id )
|
|
|
|
CALL mp_bcast( cellmin, ionode_id )
|
|
|
|
CALL mp_bcast( cellmax, ionode_id )
|
|
|
|
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Broadcast variables values for Namelist ELECTRONS
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE electrons_bcast()
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY : ionode_id
|
|
|
|
USE mp, ONLY : mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( emass, ionode_id )
|
|
|
|
CALL mp_bcast( emass_cutoff, ionode_id )
|
|
|
|
CALL mp_bcast( orthogonalization, ionode_id )
|
|
|
|
CALL mp_bcast( electron_maxstep, ionode_id )
|
|
|
|
CALL mp_bcast( ortho_eps, ionode_id )
|
|
|
|
CALL mp_bcast( ortho_max, ionode_id )
|
|
|
|
CALL mp_bcast( electron_dynamics, ionode_id )
|
|
|
|
CALL mp_bcast( electron_damping, ionode_id )
|
|
|
|
CALL mp_bcast( electron_velocities, ionode_id )
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL mp_bcast( electron_temperature, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( conv_thr, ionode_id )
|
|
|
|
CALL mp_bcast( ekincw, ionode_id )
|
|
|
|
CALL mp_bcast( fnosee, ionode_id )
|
|
|
|
CALL mp_bcast( startingwfc, ionode_id )
|
|
|
|
CALL mp_bcast( ampre, ionode_id )
|
|
|
|
CALL mp_bcast( grease, ionode_id )
|
|
|
|
CALL mp_bcast( startingpot, ionode_id )
|
|
|
|
CALL mp_bcast( diis_size, ionode_id )
|
|
|
|
CALL mp_bcast( diis_nreset, ionode_id )
|
|
|
|
CALL mp_bcast( diis_hcut, ionode_id )
|
|
|
|
CALL mp_bcast( diis_wthr, ionode_id )
|
|
|
|
CALL mp_bcast( diis_delt, ionode_id )
|
|
|
|
CALL mp_bcast( diis_maxstep, ionode_id )
|
|
|
|
CALL mp_bcast( diis_rot, ionode_id )
|
|
|
|
CALL mp_bcast( diis_fthr, ionode_id )
|
|
|
|
CALL mp_bcast( diis_temp, ionode_id )
|
|
|
|
CALL mp_bcast( diis_achmix, ionode_id )
|
|
|
|
CALL mp_bcast( diis_g0chmix, ionode_id )
|
|
|
|
CALL mp_bcast( diis_g1chmix, ionode_id )
|
|
|
|
CALL mp_bcast( diis_nchmix, ionode_id )
|
|
|
|
CALL mp_bcast( diis_nrot, ionode_id )
|
|
|
|
CALL mp_bcast( diis_rothr, ionode_id )
|
|
|
|
CALL mp_bcast( diis_ethr, ionode_id )
|
|
|
|
CALL mp_bcast( diis_chguess, ionode_id )
|
|
|
|
CALL mp_bcast( mixing_fixed_ns, ionode_id )
|
|
|
|
CALL mp_bcast( mixing_mode, ionode_id )
|
|
|
|
CALL mp_bcast( mixing_beta, ionode_id )
|
|
|
|
CALL mp_bcast( mixing_ndim, ionode_id )
|
2006-09-16 02:35:28 +08:00
|
|
|
CALL mp_bcast( tqr, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( diagonalization, ionode_id )
|
|
|
|
CALL mp_bcast( diago_thr_init, ionode_id )
|
|
|
|
CALL mp_bcast( diago_cg_maxiter, ionode_id )
|
|
|
|
CALL mp_bcast( diago_david_ndim, ionode_id )
|
2006-04-25 01:32:08 +08:00
|
|
|
CALL mp_bcast( diago_full_acc, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( sic, ionode_id )
|
|
|
|
CALL mp_bcast( sic_epsilon , ionode_id )
|
2005-12-06 22:55:23 +08:00
|
|
|
CALL mp_bcast( sic_alpha , ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( force_pairing , ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
! ... ensemble-DFT
|
|
|
|
!
|
2004-12-21 23:48:19 +08:00
|
|
|
CALL mp_bcast( fermi_energy, ionode_id )
|
|
|
|
CALL mp_bcast( n_inner, ionode_id )
|
2007-03-22 07:45:24 +08:00
|
|
|
CALL mp_bcast( niter_cold_restart, ionode_id )
|
|
|
|
CALL mp_bcast( lambda_cold, ionode_id )
|
2004-12-21 23:48:19 +08:00
|
|
|
CALL mp_bcast( rotation_dynamics, ionode_id )
|
|
|
|
CALL mp_bcast( occupation_dynamics,ionode_id )
|
|
|
|
CALL mp_bcast( rotmass, ionode_id )
|
|
|
|
CALL mp_bcast( occmass, ionode_id )
|
|
|
|
CALL mp_bcast( rotation_damping, ionode_id )
|
|
|
|
CALL mp_bcast( occupation_damping, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
!
|
|
|
|
! ... conjugate gradient
|
|
|
|
!
|
|
|
|
CALL mp_bcast( tcg, ionode_id )
|
|
|
|
CALL mp_bcast( maxiter, ionode_id )
|
|
|
|
CALL mp_bcast( etresh, ionode_id )
|
|
|
|
CALL mp_bcast( passop, ionode_id )
|
2007-03-22 07:45:24 +08:00
|
|
|
CALL mp_bcast( niter_cg_restart, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
!
|
|
|
|
! ... electric field
|
|
|
|
!
|
|
|
|
CALL mp_bcast( epol, ionode_id )
|
|
|
|
CALL mp_bcast( efield, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2006-02-02 01:57:32 +08:00
|
|
|
CALL mp_bcast( epol2, ionode_id )
|
|
|
|
CALL mp_bcast( efield2, ionode_id )
|
2009-06-02 21:27:51 +08:00
|
|
|
CALL mp_bcast( efield_cart, ionode_id )
|
2006-06-07 10:01:57 +08:00
|
|
|
!
|
2006-12-16 01:58:58 +08:00
|
|
|
! ... occupation constraints ...
|
|
|
|
!
|
|
|
|
CALL mp_bcast( occupation_constraints, ionode_id )
|
|
|
|
!
|
2009-04-03 00:05:09 +08:00
|
|
|
! ... real space ...
|
|
|
|
CALL mp_bcast( real_space, ionode_id)
|
2011-03-11 21:45:04 +08:00
|
|
|
CALL mp_bcast( adaptive_thr, ionode_id )
|
|
|
|
CALL mp_bcast( conv_thr_init, ionode_id )
|
|
|
|
CALL mp_bcast( conv_thr_multi, ionode_id )
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2009-02-12 16:07:11 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Broadcast variables values for Namelist IONS
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE ions_bcast()
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY: ionode_id
|
|
|
|
USE mp, ONLY: mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( phase_space, ionode_id )
|
|
|
|
CALL mp_bcast( ion_dynamics, ionode_id )
|
|
|
|
CALL mp_bcast( ion_radius, ionode_id )
|
|
|
|
CALL mp_bcast( ion_damping, ionode_id )
|
|
|
|
CALL mp_bcast( ion_positions, ionode_id )
|
|
|
|
CALL mp_bcast( ion_velocities, ionode_id )
|
|
|
|
CALL mp_bcast( ion_temperature, ionode_id )
|
|
|
|
CALL mp_bcast( tempw, ionode_id )
|
|
|
|
CALL mp_bcast( fnosep, ionode_id )
|
2006-06-30 06:43:58 +08:00
|
|
|
CALL mp_bcast( nhgrp, ionode_id )
|
2008-04-10 07:39:53 +08:00
|
|
|
CALL mp_bcast( fnhscl, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( nhpcl, ionode_id )
|
|
|
|
CALL mp_bcast( nhptyp, ionode_id )
|
|
|
|
CALL mp_bcast( ndega, ionode_id )
|
|
|
|
CALL mp_bcast( tranp, ionode_id )
|
|
|
|
CALL mp_bcast( amprp, ionode_id )
|
|
|
|
CALL mp_bcast( greasp, ionode_id )
|
|
|
|
CALL mp_bcast( tolp, ionode_id )
|
|
|
|
CALL mp_bcast( ion_nstepe, ionode_id )
|
|
|
|
CALL mp_bcast( ion_maxstep, ionode_id )
|
|
|
|
CALL mp_bcast( delta_t, ionode_id )
|
|
|
|
CALL mp_bcast( nraise, ionode_id )
|
2006-03-28 05:40:10 +08:00
|
|
|
CALL mp_bcast( refold_pos, ionode_id )
|
2006-06-15 22:27:14 +08:00
|
|
|
CALL mp_bcast( remove_rigid_rot, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( upscale, ionode_id )
|
2005-01-10 14:56:14 +08:00
|
|
|
CALL mp_bcast( pot_extrapolation, ionode_id )
|
|
|
|
CALL mp_bcast( wfc_extrapolation, ionode_id )
|
2004-05-03 16:48:48 +08:00
|
|
|
!
|
|
|
|
! ... BFGS
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( bfgs_ndim, ionode_id )
|
2004-05-03 16:48:48 +08:00
|
|
|
CALL mp_bcast( trust_radius_max, ionode_id )
|
|
|
|
CALL mp_bcast( trust_radius_min, ionode_id )
|
|
|
|
CALL mp_bcast( trust_radius_ini, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( w_1, ionode_id )
|
|
|
|
CALL mp_bcast( w_2, ionode_id )
|
2004-05-03 16:48:48 +08:00
|
|
|
!
|
2004-04-02 00:30:59 +08:00
|
|
|
CALL mp_bcast( sic_rloc, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( fe_step, ionode_id )
|
|
|
|
CALL mp_bcast( fe_nstep, ionode_id )
|
2006-09-13 02:28:57 +08:00
|
|
|
CALL mp_bcast( sw_nstep, ionode_id )
|
2006-12-11 22:51:54 +08:00
|
|
|
CALL mp_bcast( eq_nstep, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( g_amplitude, ionode_id )
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! Broadcast variables values for Namelist CELL
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE cell_bcast()
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY: ionode_id
|
|
|
|
USE mp, ONLY: mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( cell_parameters, ionode_id )
|
|
|
|
CALL mp_bcast( cell_dynamics, ionode_id )
|
|
|
|
CALL mp_bcast( cell_velocities, ionode_id )
|
|
|
|
CALL mp_bcast( cell_dofree, ionode_id )
|
|
|
|
CALL mp_bcast( press, ionode_id )
|
|
|
|
CALL mp_bcast( wmass, ionode_id )
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL mp_bcast( cell_temperature, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( temph, ionode_id )
|
|
|
|
CALL mp_bcast( fnoseh, ionode_id )
|
|
|
|
CALL mp_bcast( greash, ionode_id )
|
|
|
|
CALL mp_bcast( cell_factor, ionode_id )
|
|
|
|
CALL mp_bcast( cell_nstepe, ionode_id )
|
|
|
|
CALL mp_bcast( cell_damping, ionode_id )
|
2006-01-25 21:33:56 +08:00
|
|
|
CALL mp_bcast( press_conv_thr, ionode_id )
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2006-11-21 04:11:43 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Broadcast variables values for Namelist PRESS_AI
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
SUBROUTINE press_ai_bcast()
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
USE io_global, ONLY: ionode_id
|
|
|
|
USE mp, ONLY: mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
!
|
|
|
|
CALL mp_bcast( abivol, ionode_id )
|
|
|
|
CALL mp_bcast( abisur, ionode_id )
|
|
|
|
CALL mp_bcast( t_gauss, ionode_id )
|
|
|
|
CALL mp_bcast( cntr, ionode_id )
|
|
|
|
CALL mp_bcast( P_ext, ionode_id )
|
|
|
|
CALL mp_bcast( Surf_t, ionode_id )
|
|
|
|
CALL mp_bcast( pvar, ionode_id )
|
|
|
|
CALL mp_bcast( P_in, ionode_id )
|
|
|
|
CALL mp_bcast( P_fin, ionode_id )
|
|
|
|
CALL mp_bcast( delta_eps, ionode_id )
|
|
|
|
CALL mp_bcast( delta_sigma, ionode_id )
|
|
|
|
CALL mp_bcast( fill_vac, ionode_id )
|
|
|
|
CALL mp_bcast( scale_at, ionode_id )
|
|
|
|
CALL mp_bcast( n_cntr, ionode_id )
|
|
|
|
CALL mp_bcast( axis, ionode_id )
|
|
|
|
CALL mp_bcast( rho_thr, ionode_id )
|
|
|
|
CALL mp_bcast( dthr, ionode_id )
|
|
|
|
CALL mp_bcast( step_rad, ionode_id )
|
|
|
|
CALL mp_bcast( jellium, ionode_id )
|
|
|
|
CALL mp_bcast( R_j, ionode_id )
|
|
|
|
CALL mp_bcast( h_j, ionode_id )
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!=----------------------------------------------------------------------------=!
|
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
! Broadcast variables values for Namelist WANNIER
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2005-07-01 22:26:10 +08:00
|
|
|
SUBROUTINE wannier_bcast()
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
|
|
|
USE io_global, ONLY: ionode_id
|
|
|
|
USE mp, ONLY: mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( wf_efield, ionode_id )
|
|
|
|
CALL mp_bcast( wf_switch, ionode_id )
|
|
|
|
CALL mp_bcast( sw_len, ionode_id )
|
|
|
|
CALL mp_bcast( efx0, ionode_id )
|
|
|
|
CALL mp_bcast( efy0, ionode_id )
|
|
|
|
CALL mp_bcast( efz0, ionode_id )
|
|
|
|
CALL mp_bcast( efx1, ionode_id )
|
|
|
|
CALL mp_bcast( efy1, ionode_id )
|
|
|
|
CALL mp_bcast( efz1, ionode_id )
|
|
|
|
CALL mp_bcast( wfsd, ionode_id )
|
|
|
|
CALL mp_bcast( wfdt, ionode_id )
|
2011-08-01 17:57:39 +08:00
|
|
|
!=================================================================
|
|
|
|
!Lingzhu Kong
|
|
|
|
CALL mp_bcast( neigh, ionode_id )
|
|
|
|
CALL mp_bcast( poisson_eps, ionode_id )
|
|
|
|
CALL mp_bcast( dis_cutoff, ionode_id )
|
|
|
|
CALL mp_bcast( exx_ps_rcut, ionode_id )
|
|
|
|
CALL mp_bcast( exx_me_rcut, ionode_id )
|
|
|
|
CALL mp_bcast( vnbsp, ionode_id )
|
|
|
|
!=================================================================
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( maxwfdt, ionode_id )
|
|
|
|
CALL mp_bcast( wf_q, ionode_id )
|
2005-07-01 22:26:10 +08:00
|
|
|
CALL mp_bcast( wf_friction, ionode_id )
|
2005-09-20 23:17:18 +08:00
|
|
|
CALL mp_bcast( nit, ionode_id )
|
|
|
|
CALL mp_bcast( nsd, ionode_id )
|
|
|
|
CALL mp_bcast( nsteps, ionode_id )
|
|
|
|
CALL mp_bcast( tolw, ionode_id )
|
|
|
|
CALL mp_bcast( adapt, ionode_id )
|
|
|
|
CALL mp_bcast( calwf, ionode_id )
|
|
|
|
CALL mp_bcast( nwf, ionode_id )
|
|
|
|
CALL mp_bcast( wffort, ionode_id )
|
|
|
|
CALL mp_bcast( writev, ionode_id )
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
END SUBROUTINE
|
2009-02-12 16:07:11 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Broadcast variables values for Namelist WANNIER_NEW
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
SUBROUTINE wannier_ac_bcast()
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
USE io_global, ONLY: ionode_id
|
|
|
|
USE mp, ONLY: mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
!
|
|
|
|
CALL mp_bcast( plot_wannier,ionode_id )
|
|
|
|
CALL mp_bcast( use_energy_int,ionode_id )
|
|
|
|
CALL mp_bcast( print_wannier_coeff,ionode_id )
|
|
|
|
CALL mp_bcast( nwan, ionode_id )
|
|
|
|
CALL mp_bcast( plot_wan_num,ionode_id )
|
|
|
|
CALL mp_bcast( plot_wan_spin,ionode_id )
|
|
|
|
! CALL mp_bcast( wan_data,ionode_id )
|
|
|
|
CALL mp_bcast( constrain_pot, ionode_id )
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Check input values for Namelist CONTROL
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE control_checkin( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = ' control_checkin '
|
|
|
|
INTEGER :: i
|
|
|
|
LOGICAL :: allowed = .FALSE.
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
DO i = 1, SIZE( calculation_allowed )
|
2003-07-07 05:47:17 +08:00
|
|
|
IF( TRIM(calculation) == calculation_allowed(i) ) allowed = .TRUE.
|
2003-11-17 19:52:06 +08:00
|
|
|
END DO
|
|
|
|
IF( .NOT. allowed ) &
|
|
|
|
CALL errore( sub_name, ' calculation '''// &
|
2004-04-14 06:30:02 +08:00
|
|
|
& TRIM(calculation)//''' not allowed ',1)
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( ndr < 50 ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' ndr out of range ', 1 )
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( ndw > 0 .AND. ndw < 50 ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' ndw out of range ', 1 )
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( nstep < 0 ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' nstep out of range ', 1 )
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( iprint < 1 ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' iprint out of range ', 1 )
|
2004-07-20 01:19:16 +08:00
|
|
|
|
|
|
|
IF( prog == 'PW' ) THEN
|
|
|
|
IF( isave > 0 ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' isave not used in PW ' )
|
2004-07-20 01:19:16 +08:00
|
|
|
ELSE
|
|
|
|
IF( isave < 1 ) &
|
|
|
|
CALL errore( sub_name,' isave out of range ', 1 )
|
|
|
|
END IF
|
2007-03-12 21:50:45 +08:00
|
|
|
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( dt < 0.0_DP ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' dt out of range ', 1 )
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( max_seconds < 0.0_DP ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' max_seconds out of range ', 1 )
|
2004-07-20 01:19:16 +08:00
|
|
|
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ekin_conv_thr < 0.0_DP ) THEN
|
2004-07-20 01:19:16 +08:00
|
|
|
IF( prog == 'PW' ) THEN
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' ekin_conv_thr not used in PW ')
|
2007-03-12 21:50:45 +08:00
|
|
|
ELSE
|
2004-07-20 01:19:16 +08:00
|
|
|
CALL errore( sub_name,' ekin_conv_thr out of range ', 1 )
|
|
|
|
END IF
|
|
|
|
END IF
|
|
|
|
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( etot_conv_thr < 0.0_DP ) &
|
2003-08-25 06:25:53 +08:00
|
|
|
CALL errore( sub_name,' etot_conv_thr out of range ', 1 )
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( forc_conv_thr < 0.0_DP ) &
|
2009-01-13 16:29:45 +08:00
|
|
|
CALL errore( sub_name,' forc_conv_thr out of range ', 1 )
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2007-03-12 21:50:45 +08:00
|
|
|
IF( dipfield ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' dipfield not yet implemented ')
|
2007-03-12 21:50:45 +08:00
|
|
|
IF( lberry ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' lberry not implemented yet ')
|
2004-07-20 01:19:16 +08:00
|
|
|
IF( gdir /= 0 ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' gdir not used ')
|
2004-07-20 01:19:16 +08:00
|
|
|
IF( nppstr /= 0 ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' nppstr not used ')
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
2004-06-29 01:42:18 +08:00
|
|
|
!
|
2004-07-20 01:19:16 +08:00
|
|
|
IF( prog == 'PW' .AND. TRIM( restart_mode ) == 'reset_counters' ) THEN
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg ( sub_name, ' restart_mode == reset_counters' // &
|
|
|
|
& ' not implemented in PW ' )
|
2004-06-29 01:42:18 +08:00
|
|
|
END IF
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2005-09-06 17:27:34 +08:00
|
|
|
IF( refg < 0 ) &
|
|
|
|
CALL errore( sub_name, ' wrong table interval refg ', 1 )
|
- FPMD: pseudopotential variable wsg, wnl, fnl substituted with
dion, beta, bec everyware.
- subroutines formfn, compute_beta, nlsm1, nlsm2, ecc ... now are common
between FPMD and CPV, a lot of clean ups!
- Changes in stdout: relevant physical quantities ( positions velocities an cell )
are now printed with the seme format of the corresponding input card,
like in PW, as was suggested by SdG.
- exemple23 updated to reflect the new input namelist "wannier"
- Subroutine init_run now is used in FPMD too.
- WARNING in the stress computed with CP, for a pseudo with core-corrections,
a contribution is missing! Not yet fixed, I need to talk with PG for the
box staff.
- WARNING the examples reference are not updated, I'm on the IBM sp, and
I prefer to update them from a linux machine.
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2110 c92efa57-630b-4861-b058-cf58834340f0
2005-08-22 22:14:13 +08:00
|
|
|
!
|
2011-12-08 18:48:18 +08:00
|
|
|
#ifdef __LOWMEM
|
|
|
|
IF( wf_collect .EQ. .true. ) &
|
|
|
|
CALL errore( sub_name, ' wf_collect = .true. is not allowed with LOWMEM build ', 1 )
|
2011-12-11 17:08:26 +08:00
|
|
|
IF( prog /= 'CP' ) &
|
|
|
|
CALL errore( sub_name, ' LOWMEM not available in '//prog//' yet ', 1 )
|
2011-12-08 18:48:18 +08:00
|
|
|
#endif
|
|
|
|
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Check input values for Namelist SYSTEM
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE system_checkin( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = ' system_checkin '
|
2009-09-14 00:48:24 +08:00
|
|
|
INTEGER :: i
|
|
|
|
LOGICAL :: allowed
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!
|
2011-08-01 17:57:39 +08:00
|
|
|
IF( ( ibrav /= 0 ) .AND. (celldm(1) == 0.0_DP) .AND. ( a == 0.0_DP ) ) &
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL errore( ' iosys ', &
|
|
|
|
& ' invalid lattice parameters ( celldm or a )', 1 )
|
|
|
|
!
|
2006-09-14 17:55:56 +08:00
|
|
|
IF( nat < 0 ) &
|
|
|
|
CALL errore( sub_name ,' nat less than zero ', MAX( nat, 1) )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2006-09-14 17:55:56 +08:00
|
|
|
IF( ntyp < 0 ) &
|
|
|
|
CALL errore( sub_name ,' ntyp less than zero ', MAX( ntyp, 1) )
|
|
|
|
IF( ntyp < 0 .OR. ntyp > nsx ) &
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL errore( sub_name , &
|
|
|
|
& ' ntyp too large, increase NSX ', MAX( ntyp, 1) )
|
|
|
|
!
|
2009-03-11 21:21:35 +08:00
|
|
|
IF( nspin < 1 .OR. nspin > 4 .OR. nspin == 3 ) &
|
2003-08-21 00:16:26 +08:00
|
|
|
CALL errore( sub_name ,' nspin out of range ', MAX(nspin, 1 ) )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ecutwfc <= 0.0_DP ) &
|
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 errore( sub_name ,' ecutwfc out of range ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ecutrho < 0.0_DP ) &
|
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 errore( sub_name ,' ecutrho out of range ',1)
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( degauss /= 0.0_DP ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name ,' degauss is not used in CP ')
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ecfixed < 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name ,' ecfixed out of range ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( qcutz < 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name ,' qcutz out of range ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( q2sigma < 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name ,' q2sigma out of range ',1)
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
|
|
|
IF( ANY(starting_magnetization /= SM_NOT_SET ) ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL infomsg( sub_name ,&
|
2007-06-27 00:46:01 +08:00
|
|
|
& ' starting_magnetization is not used in CP ')
|
2010-05-08 17:58:35 +08:00
|
|
|
! IF( lda_plus_U ) &
|
|
|
|
! CALL infomsg( sub_name ,' lda_plus_U is not used in CP ')
|
2005-12-28 18:42:31 +08:00
|
|
|
IF( la2F ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name ,' la2F is not used in CP ')
|
2010-05-08 17:58:35 +08:00
|
|
|
! IF( ANY(Hubbard_U /= 0.0_DP) ) &
|
|
|
|
! CALL infomsg( sub_name ,' Hubbard_U is not used in CP ')
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ANY(Hubbard_alpha /= 0.0_DP) ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name ,' Hubbard_alpha is not used in CP ')
|
2003-07-07 05:47:17 +08:00
|
|
|
IF( nosym ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name ,' nosym not implemented in CP ')
|
2008-09-16 21:53:01 +08:00
|
|
|
IF( nosym_evc ) &
|
|
|
|
CALL infomsg( sub_name ,' nosym_evc not implemented in CP ')
|
2008-03-11 23:38:23 +08:00
|
|
|
IF( noinv ) &
|
|
|
|
CALL infomsg( sub_name ,' noinv not implemented in CP ')
|
2003-10-30 21:58:42 +08:00
|
|
|
END IF
|
2004-11-03 17:53:12 +08:00
|
|
|
!
|
|
|
|
! ... control on SIC variables
|
|
|
|
!
|
|
|
|
IF ( sic /= 'none' ) THEN
|
|
|
|
!
|
2007-06-12 01:13:15 +08:00
|
|
|
IF (sic_epsilon > 1.0_DP ) &
|
2004-11-03 17:53:12 +08:00
|
|
|
CALL errore( sub_name, &
|
|
|
|
& ' invalid sic_epsilon, greater than 1.',1 )
|
2007-06-12 01:13:15 +08:00
|
|
|
IF (sic_epsilon < 0.0_DP ) &
|
2004-11-03 17:53:12 +08:00
|
|
|
CALL errore( sub_name, &
|
|
|
|
& ' invalid sic_epsilon, less than 0 ',1 )
|
2007-06-12 01:13:15 +08:00
|
|
|
IF (sic_alpha > 1.0_DP ) &
|
2005-12-06 22:55:23 +08:00
|
|
|
CALL errore( sub_name, &
|
|
|
|
& ' invalid sic_alpha, greater than 1.',1 )
|
2007-06-12 01:13:15 +08:00
|
|
|
IF (sic_alpha < 0.0_DP ) &
|
2005-12-06 22:55:23 +08:00
|
|
|
CALL errore( sub_name, &
|
|
|
|
& ' invalid sic_alpha, less than 0 ',1 )
|
2004-11-03 17:53:12 +08:00
|
|
|
!
|
2004-04-03 22:42:28 +08:00
|
|
|
IF ( .NOT. force_pairing ) &
|
2004-11-03 17:53:12 +08:00
|
|
|
CALL errore( sub_name, &
|
|
|
|
& ' invalid force_pairing with sic activated', 1 )
|
|
|
|
IF ( nspin /= 2 ) &
|
|
|
|
CALL errore( sub_name, &
|
|
|
|
& ' invalid nspin with sic activated', 1 )
|
2009-09-18 21:30:27 +08:00
|
|
|
IF ( tot_magnetization /= 1._DP ) &
|
2004-11-03 17:53:12 +08:00
|
|
|
CALL errore( sub_name, &
|
2009-09-18 21:30:27 +08:00
|
|
|
& ' invalid tot_magnetization_ with sic activated', 1 )
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2004-11-03 17:53:12 +08:00
|
|
|
ENDIF
|
|
|
|
!
|
2009-09-14 00:48:24 +08:00
|
|
|
! ... control on EXX variables
|
|
|
|
!
|
|
|
|
DO i = 1, SIZE( exxdiv_treatment_allowed )
|
|
|
|
IF( TRIM(exxdiv_treatment) == exxdiv_treatment_allowed(i) ) allowed = .TRUE.
|
|
|
|
END DO
|
|
|
|
IF( .NOT. allowed ) CALL errore(sub_name, &
|
|
|
|
' invalid exxdiv_treatment: '//TRIM(exxdiv_treatment), 1 )
|
|
|
|
!
|
|
|
|
IF ( TRIM(exxdiv_treatment) == "yukawa" .AND. yukawa <= 0.0 ) &
|
|
|
|
CALL errore(sub_name, ' invalid value for yukawa', 1 )
|
|
|
|
!
|
|
|
|
IF ( TRIM(exxdiv_treatment) == "vcut_ws" .AND. ecutvcut <= 0.0 ) &
|
|
|
|
CALL errore(sub_name, ' invalid value for ecutvcut', 1 )
|
|
|
|
!
|
|
|
|
IF ( x_gamma_extrapolation .AND. ( TRIM(exxdiv_treatment) == "vcut_ws" .OR. &
|
|
|
|
TRIM(exxdiv_treatment) == "vcut_spherical" ) ) &
|
|
|
|
CALL errore(sub_name, ' x_gamma_extrapolation cannot be used with vcut', 1 )
|
|
|
|
!
|
2003-10-30 21:58:42 +08:00
|
|
|
RETURN
|
2004-11-03 17:53:12 +08:00
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Check input values for Namelist ELECTRONS
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE electrons_checkin( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = ' electrons_checkin '
|
|
|
|
INTEGER :: i
|
|
|
|
LOGICAL :: allowed = .FALSE.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
DO i = 1, SIZE(electron_dynamics_allowed)
|
|
|
|
IF( TRIM(electron_dynamics) == &
|
|
|
|
electron_dynamics_allowed(i) ) allowed = .TRUE.
|
|
|
|
END DO
|
|
|
|
IF( .NOT. allowed ) &
|
|
|
|
CALL errore( sub_name, ' electron_dynamics '''//&
|
2004-04-14 06:30:02 +08:00
|
|
|
& TRIM(electron_dynamics)//''' not allowed ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( emass <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name, ' emass less or equal 0 ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( emass_cutoff <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name, ' emass_cutoff less or equal 0 ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ortho_eps <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name, ' ortho_eps less or equal 0 ',1)
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( ortho_max < 1 ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name, ' ortho_max less than 1 ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( fnosee <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name, ' fnosee less or equal 0 ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( ekincw <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name, ' ekincw less or equal 0 ',1)
|
2006-12-16 01:58:58 +08:00
|
|
|
IF( occupation_constraints ) &
|
|
|
|
CALL errore( sub_name, ' occupation_constraints not yet implemented ',1)
|
2004-04-02 00:30:59 +08:00
|
|
|
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Check input values for Namelist IONS
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE ions_checkin( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = ' ions_checkin '
|
|
|
|
INTEGER :: i
|
|
|
|
LOGICAL :: allowed = .FALSE.
|
|
|
|
!
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
DO i = 1, SIZE( phase_space_allowed )
|
|
|
|
IF( TRIM( phase_space ) == phase_space_allowed(i) ) allowed = .TRUE.
|
|
|
|
END DO
|
|
|
|
IF ( .NOT. allowed ) &
|
|
|
|
CALL errore( sub_name, ' phase_space '''// &
|
|
|
|
& TRIM( phase_space )// ''' not allowed ', 1 )
|
|
|
|
!
|
2006-12-13 01:33:36 +08:00
|
|
|
allowed = .FALSE.
|
2003-11-17 19:52:06 +08:00
|
|
|
DO i = 1, SIZE(ion_dynamics_allowed)
|
2003-07-07 05:47:17 +08:00
|
|
|
IF( TRIM(ion_dynamics) == ion_dynamics_allowed(i) ) allowed = .TRUE.
|
2003-11-17 19:52:06 +08:00
|
|
|
END DO
|
|
|
|
IF( .NOT. allowed ) &
|
|
|
|
CALL errore( sub_name, ' ion_dynamics '''// &
|
2004-04-14 06:30:02 +08:00
|
|
|
& TRIM(ion_dynamics)//''' not allowed ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( tempw <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name,' tempw out of range ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( fnosep( 1 ) <= 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name,' fnosep out of range ',1)
|
2005-03-15 22:35:47 +08:00
|
|
|
IF( nhpcl > nhclm ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg ( sub_name,' nhpcl should be less than nhclm')
|
2005-03-15 22:35:47 +08:00
|
|
|
IF( nhpcl < 0 ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg ( sub_name,' nhpcl out of range ')
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( ion_nstepe <= 0 ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name,' ion_nstepe out of range ',1)
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( ion_maxstep < 0 ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name,' ion_maxstep out of range ',1)
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
IF (sic /= 'none' .and. sic_rloc == 0.0_DP) &
|
|
|
|
CALL errore( sub_name, ' invalid sic_rloc with sic activated ', 1 )
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Check input values for Namelist CELL
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE cell_checkin( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = ' cell_checkin '
|
|
|
|
INTEGER :: i
|
|
|
|
LOGICAL :: allowed = .FALSE.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
DO i = 1, SIZE(cell_dynamics_allowed)
|
|
|
|
IF( TRIM(cell_dynamics) == &
|
|
|
|
cell_dynamics_allowed(i) ) allowed = .TRUE.
|
|
|
|
END DO
|
|
|
|
IF( .NOT. allowed ) &
|
|
|
|
CALL errore( sub_name, ' cell_dynamics '''// &
|
2004-04-14 06:30:02 +08:00
|
|
|
TRIM(cell_dynamics)//''' not allowed ',1)
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( wmass < 0.0_DP ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name,' wmass out of range ',1)
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2007-06-12 01:13:15 +08:00
|
|
|
IF( cell_factor /= 0.0_DP ) &
|
2007-06-27 00:46:01 +08:00
|
|
|
CALL infomsg( sub_name,' cell_factor not used in CP ')
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
IF( cell_nstepe <= 0 ) &
|
2004-04-14 06:30:02 +08:00
|
|
|
CALL errore( sub_name,' cell_nstepe out of range ',1)
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
! Check input values for Namelist WANNIER
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2005-07-01 22:26:10 +08:00
|
|
|
SUBROUTINE wannier_checkin( prog )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2005-08-24 02:42:25 +08:00
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = 'wannier_checkin'
|
|
|
|
!
|
|
|
|
IF ( calwf < 1 .OR. calwf > 5 ) &
|
|
|
|
CALL errore( sub_name, ' calwf out of range ', 1 )
|
2005-07-01 22:26:10 +08:00
|
|
|
!
|
2009-03-04 18:43:48 +08:00
|
|
|
IF ( wfsd < 1 .OR. wfsd > 3 ) &
|
|
|
|
CALL errore( sub_name, ' wfsd out of range ', 1 ) !
|
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2009-02-12 16:07:11 +08:00
|
|
|
! Check input values for Namelist WANNIER_NEW
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
|
|
|
SUBROUTINE wannier_ac_checkin( prog )
|
|
|
|
!--------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = 'wannier_new_checkin'
|
|
|
|
!
|
|
|
|
!
|
|
|
|
IF ( nwan > nwanx ) &
|
|
|
|
CALL errore( sub_name, ' nwan out of range ', 1 )
|
|
|
|
|
|
|
|
IF ( plot_wan_num < 0 .OR. plot_wan_num > nwan ) &
|
|
|
|
CALL errore( sub_name, ' plot_wan_num out of range ', 1 )
|
|
|
|
|
|
|
|
IF ( plot_wan_spin < 0 .OR. plot_wan_spin > 2 ) &
|
|
|
|
CALL errore( sub_name, ' plot_wan_spin out of range ', 1 )
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! Set values according to the "calculation" variable
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
SUBROUTINE fixval( prog )
|
2007-03-12 21:50:45 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2006-06-23 21:46:13 +08:00
|
|
|
USE constants, ONLY : e2
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
CHARACTER(LEN=20) :: sub_name = ' fixval '
|
|
|
|
!
|
|
|
|
!
|
2005-09-17 10:14:39 +08:00
|
|
|
SELECT CASE( TRIM( calculation ) )
|
2003-07-07 05:47:17 +08:00
|
|
|
CASE ('scf')
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
|
|
|
electron_dynamics = 'damp'
|
2005-10-02 00:04:41 +08:00
|
|
|
ion_dynamics = 'none'
|
|
|
|
cell_dynamics = 'none'
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
2006-03-21 06:44:35 +08:00
|
|
|
CASE ('nscf', 'bands')
|
2003-11-17 19:52:06 +08:00
|
|
|
IF( prog == 'CP' ) occupations = 'bogus'
|
|
|
|
IF( prog == 'CP' ) electron_dynamics = 'damp'
|
2004-07-07 23:23:47 +08:00
|
|
|
CASE ( 'cp-wf' )
|
|
|
|
IF( prog == 'CP' ) THEN
|
|
|
|
electron_dynamics = 'damp'
|
2005-07-01 22:26:10 +08:00
|
|
|
ion_dynamics = 'damp'
|
2004-07-07 23:23:47 +08:00
|
|
|
END IF
|
2005-09-19 07:49:24 +08:00
|
|
|
IF ( prog == 'PW' ) &
|
2005-07-01 22:26:10 +08:00
|
|
|
CALL errore( sub_name, ' calculation ' // &
|
|
|
|
& TRIM( calculation ) // ' not implemented ', 1 )
|
2011-08-01 17:57:39 +08:00
|
|
|
!=========================================================================
|
|
|
|
!Lingzhu Kong
|
|
|
|
CASE ( 'cp-wf-nscf','cp-wf-pbe0','pbe0-nscf' )
|
|
|
|
IF( prog == 'CP' ) THEN
|
2012-03-09 17:58:43 +08:00
|
|
|
occupations = 'fixed'
|
2011-08-01 17:57:39 +08:00
|
|
|
electron_dynamics = 'damp'
|
|
|
|
ion_dynamics = 'damp'
|
|
|
|
END IF
|
|
|
|
IF ( prog == 'PW' ) &
|
|
|
|
CALL errore( sub_name, ' calculation ' // &
|
|
|
|
& TRIM( calculation ) // ' not implemented ', 1 )
|
|
|
|
!=========================================================================
|
2003-07-07 05:47:17 +08:00
|
|
|
CASE ('relax')
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2003-11-17 19:52:06 +08:00
|
|
|
electron_dynamics = 'damp'
|
2005-10-02 00:04:41 +08:00
|
|
|
ion_dynamics = 'damp'
|
2003-11-17 19:52:06 +08:00
|
|
|
ELSE IF( prog == 'PW' ) THEN
|
|
|
|
ion_dynamics = 'bfgs'
|
|
|
|
END IF
|
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
|
|
|
CASE ( 'md', 'cp' )
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2003-11-17 19:52:06 +08:00
|
|
|
electron_dynamics = 'verlet'
|
2005-10-02 00:04:41 +08:00
|
|
|
ion_dynamics = 'verlet'
|
2003-11-17 19:52:06 +08:00
|
|
|
ELSE IF( prog == 'PW' ) THEN
|
|
|
|
ion_dynamics = 'verlet'
|
|
|
|
END IF
|
2003-07-07 05:47:17 +08:00
|
|
|
CASE ('vc-relax')
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2003-11-17 19:52:06 +08:00
|
|
|
electron_dynamics = 'damp'
|
2005-10-02 00:04:41 +08:00
|
|
|
ion_dynamics = 'damp'
|
|
|
|
cell_dynamics = 'damp-pr'
|
2003-11-17 19:52:06 +08:00
|
|
|
ELSE IF( prog == 'PW' ) THEN
|
2009-06-12 22:48:23 +08:00
|
|
|
ion_dynamics = 'bfgs'
|
|
|
|
cell_dynamics= 'bfgs'
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
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
|
|
|
CASE ( 'vc-md', 'vc-cp' )
|
2005-09-19 07:49:24 +08:00
|
|
|
IF( prog == 'CP' ) THEN
|
2003-11-17 19:52:06 +08:00
|
|
|
electron_dynamics = 'verlet'
|
2005-10-02 00:04:41 +08:00
|
|
|
ion_dynamics = 'verlet'
|
|
|
|
cell_dynamics = 'pr'
|
2003-11-17 19:52:06 +08:00
|
|
|
ELSE IF( prog == 'PW' ) THEN
|
|
|
|
ion_dynamics = 'beeman'
|
|
|
|
END IF
|
2005-07-29 00:30:19 +08:00
|
|
|
!
|
2003-07-07 05:47:17 +08:00
|
|
|
CASE DEFAULT
|
2005-09-17 10:14:39 +08:00
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
CALL errore( sub_name,' calculation '// &
|
2005-12-29 22:27:31 +08:00
|
|
|
& TRIM(calculation)//' not implemented ', 1 )
|
2005-09-17 10:14:39 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
END SELECT
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2004-03-20 23:43:47 +08:00
|
|
|
IF ( prog == 'PW' ) THEN
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2009-09-25 03:44:04 +08:00
|
|
|
IF ( calculation == 'nscf' .OR. calculation == 'bands' ) THEN
|
2006-03-23 02:58:38 +08:00
|
|
|
!
|
|
|
|
startingpot = 'file'
|
2011-12-24 01:23:56 +08:00
|
|
|
startingwfc = 'atomic+random'
|
2006-03-23 02:58:38 +08:00
|
|
|
!
|
|
|
|
ELSE IF ( restart_mode == "from_scratch" ) THEN
|
2004-03-20 23:43:47 +08:00
|
|
|
!
|
2011-12-24 01:23:56 +08:00
|
|
|
startingwfc = 'atomic+random'
|
2004-03-20 23:43:47 +08:00
|
|
|
startingpot = 'atomic'
|
|
|
|
!
|
|
|
|
ELSE
|
|
|
|
!
|
|
|
|
startingwfc = 'file'
|
2007-03-12 21:50:45 +08:00
|
|
|
startingpot = 'file'
|
2004-03-20 23:43:47 +08:00
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
2011-12-24 01:23:56 +08:00
|
|
|
ELSE IF ( prog == 'CP' ) THEN
|
|
|
|
!
|
|
|
|
startingwfc = 'random'
|
|
|
|
startingpot = ' '
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
END IF
|
2004-08-18 23:53:01 +08:00
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
IF ( TRIM( sic ) /= 'none' ) THEN
|
2009-11-24 17:10:27 +08:00
|
|
|
force_pairing = ( nspin == 2 .AND. ( tot_magnetization==0._dp .OR. &
|
|
|
|
tot_magnetization==1._dp ) )
|
2004-04-02 00:30:59 +08:00
|
|
|
END IF
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
2005-05-27 06:42:05 +08:00
|
|
|
END SUBROUTINE
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
|
|
|
! Namelist parsing main routine
|
|
|
|
!
|
|
|
|
!=----------------------------------------------------------------------=!
|
|
|
|
!
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2011-01-11 01:12:55 +08:00
|
|
|
SUBROUTINE read_namelists( prog, unit )
|
2006-06-07 10:01:57 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
! this routine reads data from standard input and puts them into
|
|
|
|
! module-scope variables (accessible from other routines by including
|
|
|
|
! this module, or the one that contains them)
|
|
|
|
! ----------------------------------------------
|
|
|
|
!
|
|
|
|
! ... declare modules
|
|
|
|
!
|
|
|
|
USE io_global, ONLY : ionode, ionode_id
|
|
|
|
USE mp, ONLY : mp_bcast
|
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
! ... declare variables
|
|
|
|
!
|
|
|
|
CHARACTER(LEN=2) :: prog ! ... specify the calling program
|
|
|
|
! prog = 'PW' pwscf
|
|
|
|
! prog = 'CP' cpr
|
|
|
|
!
|
2011-01-11 01:12:55 +08:00
|
|
|
!
|
|
|
|
INTEGER, INTENT(IN), optional :: unit
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! ... declare other variables
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
INTEGER :: ios
|
|
|
|
!
|
2011-01-11 01:12:55 +08:00
|
|
|
INTEGER :: unit_loc=5
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
! ... end of declarations
|
|
|
|
!
|
|
|
|
! ----------------------------------------------
|
|
|
|
!
|
2011-01-11 01:12:55 +08:00
|
|
|
IF(PRESENT(unit)) unit_loc = unit
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2011-08-01 17:57:39 +08:00
|
|
|
IF( prog /= 'PW' .AND. prog /= 'CP' ) &
|
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 errore( ' read_namelists ', ' unknown calling program ', 1 )
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2004-03-20 23:43:47 +08:00
|
|
|
! ... default settings for all namelists
|
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
IF( prog == 'PW' .OR. prog == 'CP') THEN
|
|
|
|
CALL control_defaults( prog )
|
|
|
|
CALL system_defaults( prog )
|
|
|
|
CALL electrons_defaults( prog )
|
|
|
|
CALL ions_defaults( prog )
|
|
|
|
CALL cell_defaults( prog )
|
2012-02-15 00:18:50 +08:00
|
|
|
#ifdef __ENVIRON
|
|
|
|
CALL environ_defaults( prog )
|
2011-04-22 00:12:36 +08:00
|
|
|
#endif
|
2010-09-29 22:16:57 +08:00
|
|
|
CALL ee_defaults( prog )
|
|
|
|
ENDIF
|
2004-03-20 23:43:47 +08:00
|
|
|
!
|
2004-03-12 18:51:49 +08:00
|
|
|
! ... Here start reading standard input file
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2004-03-12 18:51:49 +08:00
|
|
|
! ... CONTROL namelist
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
IF(prog == 'PW' .OR. prog == 'CP' ) THEN
|
2003-11-17 19:52:06 +08:00
|
|
|
ios = 0
|
|
|
|
IF( ionode ) THEN
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, control, iostat = ios )
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist control ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL control_bcast( )
|
|
|
|
CALL control_checkin( prog )
|
|
|
|
!
|
2006-12-13 01:33:36 +08:00
|
|
|
! ... fixval changes some default values according to the value
|
|
|
|
! ... of "calculation" read in CONTROL namelist
|
2004-03-20 23:43:47 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL fixval( prog )
|
|
|
|
!
|
2004-03-12 18:51:49 +08:00
|
|
|
! ... SYSTEM namelist
|
2007-03-12 21:50:45 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
ios = 0
|
|
|
|
IF( ionode ) THEN
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, system, iostat = ios )
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist system ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL system_bcast( )
|
2006-03-24 02:00:12 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
CALL system_checkin( prog )
|
|
|
|
!
|
2010-12-21 06:14:44 +08:00
|
|
|
! CALL allocate_input_ions( ntyp, nat )
|
2006-03-24 02:00:12 +08:00
|
|
|
!
|
2004-03-12 18:51:49 +08:00
|
|
|
! ... ELECTRONS namelist
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
ios = 0
|
|
|
|
IF( ionode ) THEN
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, electrons, iostat = ios )
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist electrons ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL electrons_bcast( )
|
|
|
|
CALL electrons_checkin( prog )
|
|
|
|
!
|
2007-03-12 21:50:45 +08:00
|
|
|
! ... IONS namelist
|
2003-11-17 19:52:06 +08:00
|
|
|
!
|
|
|
|
ios = 0
|
2005-09-17 10:14:39 +08:00
|
|
|
IF ( ionode ) THEN
|
|
|
|
!
|
|
|
|
IF ( TRIM( calculation ) == 'relax' .OR. &
|
|
|
|
TRIM( calculation ) == 'md' .OR. &
|
|
|
|
TRIM( calculation ) == 'vc-relax' .OR. &
|
|
|
|
TRIM( calculation ) == 'vc-md' .OR. &
|
|
|
|
TRIM( calculation ) == 'cp' .OR. &
|
|
|
|
TRIM( calculation ) == 'vc-cp' .OR. &
|
|
|
|
TRIM( calculation ) == 'smd' .OR. &
|
2011-08-01 17:57:39 +08:00
|
|
|
TRIM( calculation ) == 'cp-wf-nscf' .OR. & !Lingzhu Kong
|
|
|
|
TRIM( calculation ) == 'cp-wf-pbe0' .OR. & !Lingzhu Kong
|
|
|
|
TRIM( calculation ) == 'pbe0-nscf' .OR. & !Lingzhu Kong
|
|
|
|
TRIM( calculation ) == 'cp-wf' ) READ( 5, ions, iostat = ios )
|
2010-11-02 05:08:30 +08:00
|
|
|
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist ions ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL ions_bcast( )
|
|
|
|
CALL ions_checkin( prog )
|
|
|
|
!
|
2004-03-12 18:51:49 +08:00
|
|
|
! ... CELL namelist
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
ios = 0
|
|
|
|
IF( ionode ) THEN
|
|
|
|
IF( TRIM( calculation ) == 'vc-relax' .OR. &
|
|
|
|
TRIM( calculation ) == 'vc-cp' .OR. &
|
|
|
|
TRIM( calculation ) == 'vc-md' .OR. &
|
|
|
|
TRIM( calculation ) == 'vc-md' ) THEN
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, cell, iostat = ios )
|
2003-07-07 05:47:17 +08:00
|
|
|
END IF
|
2003-11-17 19:52:06 +08:00
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist cell ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL cell_bcast()
|
|
|
|
CALL cell_checkin( prog )
|
|
|
|
!
|
2007-04-07 12:59:19 +08:00
|
|
|
ios = 0
|
2006-11-21 04:11:43 +08:00
|
|
|
IF( ionode ) THEN
|
|
|
|
if (tabps) then
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, press_ai, iostat = ios )
|
2006-11-21 04:11:43 +08:00
|
|
|
end if
|
|
|
|
END IF
|
2008-04-26 17:01:11 +08:00
|
|
|
CALL mp_bcast( ios, ionode_id )
|
2007-04-07 12:59:19 +08:00
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist press_ai ', ABS(ios) )
|
|
|
|
END IF
|
2006-11-21 04:11:43 +08:00
|
|
|
!
|
|
|
|
CALL press_ai_bcast()
|
2012-02-15 00:18:50 +08:00
|
|
|
#ifdef __ENVIRON
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
! ... ENVIRON namelist
|
2011-04-22 00:12:36 +08:00
|
|
|
!
|
2012-02-15 00:18:50 +08:00
|
|
|
IF ( do_environ ) THEN
|
2011-04-22 00:12:36 +08:00
|
|
|
ios = 0
|
2012-02-15 00:18:50 +08:00
|
|
|
IF( ionode ) READ( 5, environ, iostat = ios )
|
2011-04-22 00:12:36 +08:00
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) CALL errore( ' read_namelists ', &
|
2012-02-15 00:18:50 +08:00
|
|
|
& ' reading namelist environ ', ABS(ios) )
|
2011-04-22 00:12:36 +08:00
|
|
|
END IF
|
2012-02-15 00:18:50 +08:00
|
|
|
CALL environ_bcast()
|
2011-04-22 00:12:36 +08:00
|
|
|
#endif
|
2006-11-21 04:11:43 +08:00
|
|
|
!
|
2008-06-11 18:47:40 +08:00
|
|
|
! ... EE namelist
|
|
|
|
!
|
2009-09-16 04:29:07 +08:00
|
|
|
IF ( TRIM( assume_isolated ) == 'dcc' ) THEN
|
2008-06-11 18:47:40 +08:00
|
|
|
ios = 0
|
2011-01-11 01:12:55 +08:00
|
|
|
IF( ionode ) READ( unit_loc, ee, iostat = ios )
|
2008-06-11 18:47:40 +08:00
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist ee ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
CALL ee_bcast()
|
|
|
|
!
|
2005-07-01 22:26:10 +08:00
|
|
|
! ... WANNIER NAMELIST
|
|
|
|
!
|
|
|
|
CALL wannier_defaults( prog )
|
|
|
|
ios = 0
|
|
|
|
IF( ionode ) THEN
|
2011-08-01 17:57:39 +08:00
|
|
|
IF( TRIM( calculation ) == 'cp-wf' .OR. & ! Lingzhu Kong
|
|
|
|
TRIM( calculation ) == 'cp-wf-nscf' .OR. & ! Lingzhu Kong
|
|
|
|
TRIM( calculation ) == 'cp-wf-pbe0' .OR. & ! Lingzhu Kong
|
|
|
|
TRIM( calculation ) == 'pbe0-nscf' ) THEN ! Lingzhu Kong
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, wannier, iostat = ios )
|
2005-07-01 22:26:10 +08:00
|
|
|
END IF
|
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist wannier ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL wannier_bcast()
|
|
|
|
CALL wannier_checkin( prog )
|
|
|
|
!
|
2009-02-12 16:07:11 +08:00
|
|
|
! ... WANNIER_NEW NAMELIST
|
|
|
|
!
|
|
|
|
CALL wannier_ac_defaults( prog )
|
|
|
|
ios = 0
|
|
|
|
IF( ionode ) THEN
|
|
|
|
IF( use_wannier ) THEN
|
2011-01-11 01:12:55 +08:00
|
|
|
READ( unit_loc, wannier_ac, iostat = ios )
|
2009-02-12 16:07:11 +08:00
|
|
|
END IF
|
|
|
|
END IF
|
|
|
|
CALL mp_bcast( ios, ionode_id )
|
|
|
|
IF( ios /= 0 ) THEN
|
|
|
|
CALL errore( ' read_namelists ', &
|
|
|
|
& ' reading namelist wannier_new ', ABS(ios) )
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL wannier_ac_bcast()
|
|
|
|
CALL wannier_ac_checkin( prog )
|
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
ENDIF
|
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE read_namelists
|
|
|
|
!
|
2010-09-29 22:16:57 +08:00
|
|
|
!
|
2003-11-17 19:52:06 +08:00
|
|
|
END MODULE read_namelists_module
|