Files
lammps/examples/PACKAGES/atc/molecule/in.polarize
2021-06-29 11:23:47 -04:00

121 lines
3.1 KiB
Plaintext

###
echo both
log log.polarize
variable nstepsequil equal 1000000
variable nsteps equal 100
variable temp equal 300.0
variable tdamp equal 0.100
variable x_dim equal 50
variable y_dim equal 50
variable z_fluid equal 37.7919
variable efieldz equal -0.0723725329978652551934
variable x_max equal ${x_dim}/2
variable y_max equal ${y_dim}/2
variable fluid_half_thickness equal ${z_fluid}/2
variable buffer equal 2.45
variable zhiwall equal ${fluid_half_thickness}+${buffer}
variable zlowall equal -${fluid_half_thickness}-${buffer}
variable zmaxatc equal ${fluid_half_thickness}+10
variable zmaxsub equal ${zmaxatc}-0.1
variable zmaxsup equal ${zmaxatc}+0.1
variable cellatc equal 5
units metal
atom_style full
dimension 3
newton off
neighbor 2 bin
neigh_modify every 1 delay 1 check no
boundary p p f
kspace_style pppm 1e-5
kspace_modify slab 3.0
pair_style lj/cut/coul/long 13.0 13.0
bond_style harmonic
angle_style harmonic
read_data waterequil.init
#velocity all create 0.0 482748 dist uniform
lattice sc 0.05
#region ATC block EDGE EDGE EDGE EDGE -${zmaxatc} ${zmaxatc} units box
#variable nelemelec equal round(2*${zmaxatc}/0.05)
region ATC block EDGE EDGE EDGE EDGE EDGE EDGE units box
# region ATC block EDGE EDGE EDGE EDGE 0.0 0.05 units box
#region ATC block -${fem_atc} ${fem_atc} -${fem_atc} ${fem_atc} -${fem_atc} ${fem_atc} units box
variable nelem equal 2
# variable nelem equal 1
group water type 1 2
group hyd type 1
group oxy type 2
timestep 0.0005
dielectric 1.0
pair_coeff 1 * 0.0 0.0
pair_coeff 2 2 0.006596 3.1507
bond_coeff 1 19.513855 0.957200
angle_coeff 1 2.385027 104.519997
special_bonds amber
#fix 2 all nvt temp ${temp} ${temp} ${tdamp}
fix lowall oxy wall/lj1043 zlo ${zlowall} 0.1351 3.1507 11.0275 units box
fix hiwall oxy wall/lj1043 zhi ${zhiwall} 0.1351 3.1507 11.0275 units box
fix efield all efield 0. 0. v_efieldz
# assign SHAKE fixes
fix 1 all shake 0.00001 500 5000 b 1 a 1
thermo 100
thermo_style custom step temp ke pe etotal f_efield f_lowall f_hiwall
variable dumpfreq equal 1
#dump 1 all atom 500 waterequil.dump
#dump_modify 1 image yes scale no flush yes
#run ${nstepsequil}
#write_restart waterequil.restart
#unfix 2
fix 3 all nve
fix AtC all atc hardy
fix_modify AtC mesh create 1 1 ${nelem} ATC p p f
#fix_modify AtC_half mesh create 1 1 ${nelem} ATC p p f
fix_modify AtC kernel cell 25 25 ${cellatc}
fix_modify AtC atom_element_map eulerian ${dumpfreq}
fix_modify AtC fields none
fix_modify AtC fields add dipole_moment mass_density
fix_modify AtC add_molecule small Water water
fix_modify AtC output polarizeFE ${dumpfreq} text vector_components
run ${nsteps}
#fix AtCH hyd atc field
#fix_modify AtCH mesh create 1 1 ${nelemelec} ATC p p f
#fix_modify AtCH atom_element_map eulerian ${dumpfreq}
#fix_modify AtCH fields add number_density
#fix_modify AtCH output FEH ${dumpfreq} text binary
#fix AtCO oxy atc field
#fix_modify AtCO mesh create 1 1 ${nelemelec} ATC p p f
#fix_modify AtCO atom_element_map eulerian ${dumpfreq}
#fix_modify AtCO fields add number_density
#fix_modify AtCO output FEO ${dumpfreq} text binary
#run ${nsteps}
#write_restart waterendpolarize.restart