Newer
Older
/*----------------------------------------------------------------------
SerialReax - Reax Force Field Simulator
Copyright (2010) Purdue University
Hasan Metin Aktulga, haktulga@cs.purdue.edu
Joseph Fogarty, jcfogart@mail.usf.edu
Sagar Pandit, pandit@usf.edu
Ananth Y Grama, ayg@cs.purdue.edu
This program is free software; you can redistribute it and/or
modify it under the terms of the GNU General Public License as
published by the Free Software Foundation; either version 2 of
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
See the GNU General Public License for more details:
<http://www.gnu.org/licenses/>.
----------------------------------------------------------------------*/
#include "forces.h"
Kurt A. O'Hearn
committed
#include "bonds.h"
#include "charges.h"
Kurt A. O'Hearn
committed
#include "hydrogen_bonds.h"
#if defined(TEST_FORCES)
#include "io_tools.h"
#endif
Kurt A. O'Hearn
committed
#include "list.h"
#include "multi_body.h"
#include "nonbonded.h"
Kurt A. O'Hearn
committed
#include "torsion_angles.h"
#include "tool_box.h"
Kurt A. O'Hearn
committed
#include "valence_angles.h"
Kurt A. O'Hearn
committed
typedef enum
{
DIAGONAL = 0,
OFF_DIAGONAL = 1,
} MATRIX_ENTRY_POSITION;
Kurt A. O'Hearn
committed
#if defined(TEST_FORCES)
static int compare_bonds( const void *p1, const void *p2 )
{
return ((bond_data *)p1)->nbr - ((bond_data *)p2)->nbr;
}
#endif
void Init_Bonded_Force_Functions( control_params *control )
Kurt A. O'Hearn
committed
control->intr_funcs[0] = &BO;
control->intr_funcs[1] = &Bonds;
control->intr_funcs[2] = &Atom_Energy;
control->intr_funcs[3] = &Valence_Angles;
control->intr_funcs[4] = &Torsion_Angles;
Kurt A. O'Hearn
committed
if ( control->hbond_cut > 0.0 )
Kurt A. O'Hearn
committed
{
control->intr_funcs[5] = &Hydrogen_Bonds;
Kurt A. O'Hearn
committed
}
else
{
Kurt A. O'Hearn
committed
control->intr_funcs[5] = NULL;
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
control->intr_funcs[6] = NULL;
control->intr_funcs[7] = NULL;
control->intr_funcs[8] = NULL;
control->intr_funcs[9] = NULL;
Kurt A. O'Hearn
committed
static void Compute_Bonded_Forces( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, static_storage *workspace, reax_list **lists,
output_controls *out_control )
Kurt A. O'Hearn
committed
#if defined(TEST_ENERGY)
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
/* Mark beginning of a new timestep in each energy file */
fprintf( out_control->ebond, "step: %d\n%6s%6s%12s%12s%12s\n",
data->step, "atom1", "atom2", "bo", "ebond", "total" );
fprintf( out_control->elp, "step: %d\n%6s%12s%12s%12s\n",
data->step, "atom", "nlp", "elp", "total" );
fprintf( out_control->eov, "step: %d\n%6s%12s%12s\n",
data->step, "atom", "eov", "total" );
fprintf( out_control->eun, "step: %d\n%6s%12s%12s\n",
data->step, "atom", "eun", "total" );
fprintf( out_control->eval, "step: %d\n%6s%6s%6s%12s%12s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "atom3",
"angle", "bo(12)", "bo(23)", "eval", "epen", "total" );
fprintf( out_control->epen, "step: %d\n%6s%6s%6s%12s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "atom3",
"angle", "bo(12)", "bo(23)", "epen", "total" );
fprintf( out_control->ecoa, "step: %d\n%6s%6s%6s%12s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "atom3",
"angle", "bo(12)", "bo(23)", "ecoa", "total" );
fprintf( out_control->ehb, "step: %d\n%6s%6s%6s%12s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "atom3",
"r(23)", "angle", "bo(12)", "ehb", "total" );
fprintf( out_control->etor, "step: %d\n%6s%6s%6s%6s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "atom3", "atom4",
"phi", "bo(23)", "etor", "total" );
fprintf( out_control->econ, "step:%d\n%6s%6s%6s%6s%12s%12s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "atom3", "atom4",
"phi", "bo(12)", "bo(23)", "bo(34)", "econ", "total" );
#endif
Kurt A. O'Hearn
committed
/* function calls for bonded interactions */
Kurt A. O'Hearn
committed
for ( i = 0; i < NUM_INTRS; i++ )
Kurt A. O'Hearn
committed
if ( control->intr_funcs[i] != NULL )
{
(control->intr_funcs[i])( system, control, data, workspace,
lists, out_control );
}
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
#if defined(TEST_FORCES)
Kurt A. O'Hearn
committed
/* function calls for printing bonded interactions */
Kurt A. O'Hearn
committed
for ( i = 0; i < NUM_INTRS; i++ )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
if ( control->print_intr_funcs[i] != NULL )
{
(control->print_intr_funcs[i])( system, control, data, workspace,
lists, out_control );
}
Kurt A. O'Hearn
committed
#endif
Kurt A. O'Hearn
committed
static void Compute_NonBonded_Forces( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, static_storage *workspace,
reax_list** lists, output_controls *out_control, int realloc )
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
#if defined(TEST_ENERGY)
fprintf( out_control->evdw, "step: %d\n%6s%6s%12s%12s%12s\n",
data->step, "atom1", "atom2", "r12", "evdw", "total" );
fprintf( out_control->ecou, "step: %d\n%6s%6s%12s%12s%12s%12s%12s\n",
data->step, "atom1", "atom2", "r12", "q1", "q2", "ecou", "total" );
Compute_Charges( system, control, data, workspace, out_control, realloc );
Kurt A. O'Hearn
committed
data->timing.cm += t_elapsed;
if ( control->cm_solver_pre_comp_refactor == -1 )
{
if ( data->step <= 4 || is_refactoring_step( control, data ) )
{
if ( is_refactoring_step( control, data ) )
{
data->timing.cm_last_pre_comp = data->timing.cm_solver_pre_comp;
}
Kurt A. O'Hearn
committed
data->timing.cm_optimum = data->timing.cm_solver_pre_app
+ data->timing.cm_solver_spmv
+ data->timing.cm_solver_vector_ops
+ data->timing.cm_solver_orthog
+ data->timing.cm_solver_tri_solve;
Kurt A. O'Hearn
committed
data->timing.cm_total_loss = 0.0;
Kurt A. O'Hearn
committed
data->timing.cm_total_loss += data->timing.cm_solver_pre_app
+ data->timing.cm_solver_spmv
+ data->timing.cm_solver_vector_ops
+ data->timing.cm_solver_orthog
+ data->timing.cm_solver_tri_solve
- data->timing.cm_optimum;
Kurt A. O'Hearn
committed
if ( control->tabulate <= 0 )
{
vdW_Coulomb_Energy( system, control, data, workspace, lists, out_control );
}
{
Tabulated_vdW_Coulomb_Energy( system, control, data, workspace,
Kurt A. O'Hearn
committed
lists, out_control );
}
Kurt A. O'Hearn
committed
#if defined(TEST_FORCES)
Print_vdW_Coulomb_Forces( system, control, data, workspace,
Kurt A. O'Hearn
committed
lists, out_control );
/* This version of Compute_Total_Force computes forces from coefficients
accumulated by all interaction functions. Saves enormous time & space! */
Kurt A. O'Hearn
committed
static void Compute_Total_Force( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, static_storage *workspace, reax_list **lists )
Kurt A. O'Hearn
committed
int i;
Kurt A. O'Hearn
committed
reax_list *bonds;
bonds = lists[BONDS];
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
#if defined(_OPENMP)
Kurt A. O'Hearn
committed
#pragma omp parallel default(shared)
Kurt A. O'Hearn
committed
#endif
Kurt A. O'Hearn
committed
{
int pj;
Kurt A. O'Hearn
committed
#if defined(_OPENMP)
Kurt A. O'Hearn
committed
int j;
#endif
if ( control->compute_pressure == FALSE
&& (control->ensemble == NVE || control->ensemble == nhNVT
|| control->ensemble == bNVT) )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
#if defined(_OPENMP)
Kurt A. O'Hearn
committed
#pragma omp for schedule(static)
Kurt A. O'Hearn
committed
#endif
Kurt A. O'Hearn
committed
for ( i = 0; i < system->N; ++i )
Kurt A. O'Hearn
committed
for ( pj = Start_Index(i, bonds); pj < End_Index(i, bonds); ++pj )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
if ( i <= bonds->bond_list[pj].nbr )
Kurt A. O'Hearn
committed
{
Add_dBond_to_Forces( i, pj, system, data, workspace, lists );
}
Kurt A. O'Hearn
committed
}
}
}
Kurt A. O'Hearn
committed
else if ( control->ensemble == sNPT || control->ensemble == iNPT
|| control->ensemble == aNPT || control->compute_pressure == TRUE )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
#if defined(_OPENMP)
Kurt A. O'Hearn
committed
#pragma omp for schedule(static)
#endif
for ( i = 0; i < system->N; ++i )
{
for ( pj = Start_Index(i, bonds); pj < End_Index(i, bonds); ++pj )
{
Kurt A. O'Hearn
committed
if ( i <= bonds->bond_list[pj].nbr )
Kurt A. O'Hearn
committed
{
Add_dBond_to_Forces_NPT( i, pj, system, data, workspace, lists );
}
}
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
#if defined(_OPENMP)
Kurt A. O'Hearn
committed
/* reduction (sum) on thread-local force vectors */
Kurt A. O'Hearn
committed
#pragma omp for schedule(static)
for ( i = 0; i < system->N; ++i )
{
for ( j = 0; j < control->num_threads; ++j )
{
rvec_Add( system->atoms[i].f, workspace->f_local[j * system->N + i] );
}
}
#endif
}
Kurt A. O'Hearn
committed
static void Validate_Lists( static_storage *workspace, reax_list **lists, int step, int n,
Kurt A. O'Hearn
committed
int Hmax, int Htop, int num_bonds, int num_hbonds )
Kurt A. O'Hearn
committed
reax_list *bonds, *hbonds;
bonds = lists[BONDS];
hbonds = lists[HBONDS];
/* far neighbors */
if ( Htop > Hmax * DANGER_ZONE )
{
workspace->realloc.Htop = Htop;
if ( Htop > Hmax )
{
fprintf( stderr,
Kurt A. O'Hearn
committed
"[ERROR] step%d - ran out of space on H matrix: Htop=%d, max = %d",
Kurt A. O'Hearn
committed
exit( INSUFFICIENT_MEMORY );
workspace->realloc.num_bonds = num_bonds;
for ( i = 0; i < n - 1; ++i )
Kurt A. O'Hearn
committed
{
if ( End_Index(i, bonds) >= Start_Index(i + 1, bonds) - 2 )
{
workspace->realloc.bonds = 1;
if ( End_Index(i, bonds) > Start_Index(i + 1, bonds) )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] step%d-bondchk failed: i=%d end(i)=%d str(i+1)=%d\n",
step, flag, End_Index(flag, bonds), Start_Index(flag + 1, bonds) );
Kurt A. O'Hearn
committed
exit( INSUFFICIENT_MEMORY );
if ( End_Index(i, bonds) >= bonds->total_intrs - 2 )
if ( End_Index(i, bonds) > bonds->total_intrs )
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] step%d-bondchk failed: i=%d end(i)=%d bond_end=%d\n",
step, flag, End_Index(i, bonds), bonds->total_intrs );
Kurt A. O'Hearn
committed
exit( INSUFFICIENT_MEMORY );
/* hbonds list */
if ( workspace->num_H > 0 )
{
flag = -1;
workspace->realloc.num_hbonds = num_hbonds;
for ( i = 0; i < workspace->num_H - 1; ++i )
Kurt A. O'Hearn
committed
{
if ( Num_Entries(i, hbonds) >=
(Start_Index(i + 1, hbonds) - Start_Index(i, hbonds)) * DANGER_ZONE )
{
workspace->realloc.hbonds = 1;
if ( End_Index(i, hbonds) > Start_Index(i + 1, hbonds) )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] step%d-hbondchk failed: i=%d end(i)=%d str(i+1)=%d\n",
step, flag, End_Index(flag, hbonds), Start_Index(flag + 1, hbonds) );
Kurt A. O'Hearn
committed
exit( INSUFFICIENT_MEMORY );
(hbonds->total_intrs - Start_Index(i, hbonds)) * DANGER_ZONE )
if ( End_Index(i, hbonds) > hbonds->total_intrs )
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] step%d-hbondchk failed: i=%d end(i)=%d hbondend=%d\n",
step, flag, End_Index(i, hbonds), hbonds->total_intrs );
Kurt A. O'Hearn
committed
exit( INSUFFICIENT_MEMORY );
Kurt A. O'Hearn
committed
static inline real Init_Charge_Matrix_Entry_Tab( reax_system *system,
Kurt A. O'Hearn
committed
control_params *control, LR_lookup_table **LR, int i, int j,
Kurt A. O'Hearn
committed
real r_ij, MATRIX_ENTRY_POSITION pos )
{
int r;
Kurt A. O'Hearn
committed
real base, dif, val, ret;
Kurt A. O'Hearn
committed
LR_lookup_table *t;
Kurt A. O'Hearn
committed
ret = 0.0;
Kurt A. O'Hearn
committed
switch ( control->charge_method )
{
case QEQ_CM:
//TODO: tabulate other portions of matrices for EE, ACKS2?
case EE_CM:
case ACKS2_CM:
Kurt A. O'Hearn
committed
switch ( pos )
{
case OFF_DIAGONAL:
Kurt A. O'Hearn
committed
t = &LR[MIN( system->atoms[i].type, system->atoms[j].type )]
[MAX( system->atoms[i].type, system->atoms[j].type )];
Kurt A. O'Hearn
committed
/* cubic spline interpolation */
r = (int)(r_ij * t->inv_dx);
if ( r == 0 )
{
++r;
}
Kurt A. O'Hearn
committed
base = (real)(r + 1) * t->dx;
dif = r_ij - base;
Kurt A. O'Hearn
committed
val = ((t->ele[r].d * dif + t->ele[r].c) * dif + t->ele[r].b)
* dif + t->ele[r].a;
Kurt A. O'Hearn
committed
val *= EV_to_KCALpMOL / C_ELE;
Kurt A. O'Hearn
committed
ret = ((i == j) ? 0.5 : 1.0) * val;
break;
case DIAGONAL:
Kurt A. O'Hearn
committed
ret = system->reax_param.sbp[system->atoms[i].type].eta;
Kurt A. O'Hearn
committed
break;
default:
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] Invalid matrix position. Terminating...\n" );
Kurt A. O'Hearn
committed
exit( INVALID_INPUT );
break;
}
break;
default:
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] Invalid charge method. Terminating...\n" );
Kurt A. O'Hearn
committed
exit( INVALID_INPUT );
break;
}
return ret;
}
static inline real Init_Charge_Matrix_Entry( reax_system *system,
control_params *control, static_storage *workspace,
int i, int j, real r_ij, MATRIX_ENTRY_POSITION pos )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
real Tap, dr3gamij_1, dr3gamij_3, ret;
Kurt A. O'Hearn
committed
ret = 0.0;
Kurt A. O'Hearn
committed
switch ( control->charge_method )
{
case QEQ_CM:
case EE_CM:
Kurt A. O'Hearn
committed
case ACKS2_CM:
Kurt A. O'Hearn
committed
switch ( pos )
{
case OFF_DIAGONAL:
Tap = workspace->Tap[7] * r_ij + workspace->Tap[6];
Tap = Tap * r_ij + workspace->Tap[5];
Tap = Tap * r_ij + workspace->Tap[4];
Tap = Tap * r_ij + workspace->Tap[3];
Tap = Tap * r_ij + workspace->Tap[2];
Tap = Tap * r_ij + workspace->Tap[1];
Tap = Tap * r_ij + workspace->Tap[0];
/* shielding */
Kurt A. O'Hearn
committed
dr3gamij_1 = r_ij * r_ij * r_ij
Kurt A. O'Hearn
committed
+ POW( system->reax_param.tbp[system->atoms[i].type][system->atoms[j].type].gamma, -3.0 );
Kurt A. O'Hearn
committed
dr3gamij_3 = POW( dr3gamij_1 , 1.0 / 3.0 );
Kurt A. O'Hearn
committed
/* i == j: periodic self-interaction term
* i != j: general interaction term */
Kurt A. O'Hearn
committed
ret = ((i == j) ? 0.5 : 1.0) * Tap * EV_to_KCALpMOL / dr3gamij_3;
break;
case DIAGONAL:
Kurt A. O'Hearn
committed
ret = system->reax_param.sbp[system->atoms[i].type].eta;
Kurt A. O'Hearn
committed
break;
default:
fprintf( stderr, "[ERROR] Invalid matrix position. Terminating...\n" );
exit( INVALID_INPUT );
break;
}
break;
Kurt A. O'Hearn
committed
default:
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] Invalid charge method. Terminating...\n" );
Kurt A. O'Hearn
committed
exit( INVALID_INPUT );
break;
}
return ret;
}
Kurt A. O'Hearn
committed
static void Init_Charge_Matrix_Remaining_Entries( reax_system *system,
Kurt A. O'Hearn
committed
control_params *control, reax_list *far_nbr_list,
Kurt A. O'Hearn
committed
sparse_matrix * H, sparse_matrix * H_sp,
Kurt A. O'Hearn
committed
int * Htop, int * H_sp_top )
{
Kurt A. O'Hearn
committed
int i, j, pj, target, val_flag;
Kurt A. O'Hearn
committed
real d, xcut, bond_softness, * X_diag;
Kurt A. O'Hearn
committed
switch ( control->charge_method )
{
case QEQ_CM:
break;
Kurt A. O'Hearn
committed
if ( system->num_molec_charge_constraints == 0 )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
H->start[system->N_cm - 1] = *Htop;
H_sp->start[system->N_cm - 1] = *H_sp_top;
for ( i = 0; i < system->N_cm - 1; ++i )
Kaymak, Cagri
committed
{
Kurt A. O'Hearn
committed
/* total charge constraint on QM atoms */
#if defined(QMMM)
Kurt A. O'Hearn
committed
if ( system->atoms[i].qmmm_mask == TRUE )
{
Kaymak, Cagri
committed
H->j[*Htop] = i;
H->val[*Htop] = 1.0;
H_sp->j[*H_sp_top] = i;
H_sp->val[*H_sp_top] = 1.0;
#if defined(QMMM)
}
Kurt A. O'Hearn
committed
#endif
Kurt A. O'Hearn
committed
*Htop = *Htop + 1;
*H_sp_top = *H_sp_top + 1;
Kaymak, Cagri
committed
}
Kurt A. O'Hearn
committed
H->j[*Htop] = system->N_cm - 1;
H->val[*Htop] = 0.0;
Kaymak, Cagri
committed
*Htop = *Htop + 1;
Kurt A. O'Hearn
committed
H_sp->j[*H_sp_top] = system->N_cm - 1;
H_sp->val[*H_sp_top] = 0.0;
Kaymak, Cagri
committed
*H_sp_top = *H_sp_top + 1;
Kurt A. O'Hearn
committed
}
else
{
for ( i = 0; i < system->num_molec_charge_constraints; ++i )
{
H->start[system->N + i] = *Htop;
H_sp->start[system->N + i] = *H_sp_top;
for ( j = system->molec_charge_constraint_ranges[2 * i];
j <= system->molec_charge_constraint_ranges[2 * i + 1]; ++j )
{
/* molecule charge constraint on QM atoms */
#if defined(QMMM)
Kurt A. O'Hearn
committed
if ( system->atoms[j - 1].qmmm_mask == TRUE )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
H->j[*Htop] = j - 1;
H->val[*Htop] = 1.0;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
H_sp->j[*H_sp_top] = j - 1;
H_sp->val[*H_sp_top] = 1.0;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
*Htop = *Htop + 1;
*H_sp_top = *H_sp_top + 1;
#if defined(QMMM)
}
#endif
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
/* explicit zeros on diagonals */
H->j[*Htop] = system->N + i;
H->val[*Htop] = 0.0;
*Htop = *Htop + 1;
H_sp->j[*H_sp_top] = system->N + i;
H_sp->val[*H_sp_top] = 0.0;
*H_sp_top = *H_sp_top + 1;
}
}
Kurt A. O'Hearn
committed
break;
case ACKS2_CM:
Kurt A. O'Hearn
committed
X_diag = smalloc( sizeof(real) * system->N, __FILE__, __LINE__ );
Kurt A. O'Hearn
committed
for ( i = 0; i < system->N; ++i )
{
Kurt A. O'Hearn
committed
X_diag[i] = 0.0;
Kurt A. O'Hearn
committed
}
for ( i = 0; i < system->N; ++i )
{
Kurt A. O'Hearn
committed
H->start[system->N + i] = *Htop;
H_sp->start[system->N + i] = *H_sp_top;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
/* constraint on ref. value for kinetic energy potential */
Kurt A. O'Hearn
committed
H->j[*Htop] = i;
Kurt A. O'Hearn
committed
*Htop = *Htop + 1;
H_sp->j[*H_sp_top] = i;
Kurt A. O'Hearn
committed
*H_sp_top = *H_sp_top + 1;
Kurt A. O'Hearn
committed
/* kinetic energy terms */
Kurt A. O'Hearn
committed
for ( pj = Start_Index(i, far_nbr_list); pj < End_Index(i, far_nbr_list); ++pj )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
/* exclude self-periodic images of atoms for
* kinetic energy term because the effective
* potential is the same on an atom and its periodic image */
Kurt A. O'Hearn
committed
if ( far_nbr_list->far_nbr_list[pj].d <= control->nonb_cut )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
j = far_nbr_list->far_nbr_list[pj].nbr;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
xcut = 0.5 * ( system->reax_param.sbp[ system->atoms[i].type ].b_s_acks2
+ system->reax_param.sbp[ system->atoms[j].type ].b_s_acks2 );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
if ( far_nbr_list->far_nbr_list[pj].d < xcut )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
d = far_nbr_list->far_nbr_list[pj].d / xcut;
Kurt A. O'Hearn
committed
bond_softness = system->reax_param.gp.l[34] * POW( d, 3.0 )
Kurt A. O'Hearn
committed
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
* POW( 1.0 - d, 6.0 );
if ( bond_softness > 0.0 )
{
val_flag = FALSE;
for ( target = H->start[system->N + i]; target < *Htop; ++target )
{
if ( H->j[target] == system->N + j )
{
H->val[target] += bond_softness;
val_flag = TRUE;
break;
}
}
if ( val_flag == FALSE )
{
H->j[*Htop] = system->N + j;
H->val[*Htop] = bond_softness;
++(*Htop);
}
val_flag = FALSE;
for ( target = H_sp->start[system->N + i]; target < *H_sp_top; ++target )
{
if ( H_sp->j[target] == system->N + j )
{
H_sp->val[target] += bond_softness;
val_flag = TRUE;
break;
}
}
if ( val_flag == FALSE )
{
H_sp->j[*H_sp_top] = system->N + j;
H_sp->val[*H_sp_top] = bond_softness;
++(*H_sp_top);
}
X_diag[i] -= bond_softness;
X_diag[j] -= bond_softness;
}
Kurt A. O'Hearn
committed
}
}
}
Kurt A. O'Hearn
committed
/* placeholders for diagonal entries, to be replaced below */
H->j[*Htop] = system->N + i;
Kurt A. O'Hearn
committed
H->val[*Htop] = 0.0;
*Htop = *Htop + 1;
Kurt A. O'Hearn
committed
H_sp->j[*H_sp_top] = system->N + i;
Kurt A. O'Hearn
committed
H_sp->val[*H_sp_top] = 0.0;
*H_sp_top = *H_sp_top + 1;
}
Kurt A. O'Hearn
committed
/* second to last row */
H->start[system->N_cm - 2] = *Htop;
H_sp->start[system->N_cm - 2] = *H_sp_top;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
/* place accumulated diagonal entries (needed second to last row marker above before this code) */
for ( i = system->N; i < 2 * system->N; ++i )
Kurt A. O'Hearn
committed
{
for ( pj = H->start[i]; pj < H->start[i + 1]; ++pj )
{
if ( H->j[pj] == i )
{
Kurt A. O'Hearn
committed
H->val[pj] = X_diag[i - system->N];
Kurt A. O'Hearn
committed
break;
}
}
for ( pj = H_sp->start[i]; pj < H_sp->start[i + 1]; ++pj )
{
if ( H_sp->j[pj] == i )
{
Kurt A. O'Hearn
committed
H_sp->val[pj] = X_diag[i - system->N];
Kurt A. O'Hearn
committed
break;
}
}
}
Kurt A. O'Hearn
committed
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
/* coupling with the kinetic energy potential */
for ( i = 0; i < system->N; ++i )
{
H->j[*Htop] = system->N + i;
H->val[*Htop] = 1.0;
*Htop = *Htop + 1;
H_sp->j[*H_sp_top] = system->N + i;
H_sp->val[*H_sp_top] = 1.0;
*H_sp_top = *H_sp_top + 1;
}
/* explicitly store zero on diagonal */
H->j[*Htop] = system->N_cm - 2;
H->val[*Htop] = 0.0;
*Htop = *Htop + 1;
H_sp->j[*H_sp_top] = system->N_cm - 2;
H_sp->val[*H_sp_top] = 0.0;
*H_sp_top = *H_sp_top + 1;
/* last row */
H->start[system->N_cm - 1] = *Htop;
H_sp->start[system->N_cm - 1] = *H_sp_top;
for ( i = 0; i < system->N; ++i )
Kurt A. O'Hearn
committed
{
H->j[*Htop] = i;
Kurt A. O'Hearn
committed
*Htop = *Htop + 1;
H_sp->j[*H_sp_top] = i;
Kurt A. O'Hearn
committed
*H_sp_top = *H_sp_top + 1;
}
Kurt A. O'Hearn
committed
/* explicitly store zero on diagonal */
Kurt A. O'Hearn
committed
H->j[*Htop] = system->N_cm - 1;
H->val[*Htop] = 0.0;
*Htop = *Htop + 1;
H_sp->j[*H_sp_top] = system->N_cm - 1;
H_sp->val[*H_sp_top] = 0.0;
*H_sp_top = *H_sp_top + 1;
Kurt A. O'Hearn
committed
sfree( X_diag, __FILE__, __LINE__ );
Kurt A. O'Hearn
committed
break;
default:
break;
}
}
Kurt A. O'Hearn
committed
/* Generate bond list (full format), hydrogen bond list (full format),
* and charge matrix (half symmetric format)
* from the far neighbors list (with distance updates, if necessary) */
Kurt A. O'Hearn
committed
static void Init_Forces( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, static_storage *workspace,
Kurt A. O'Hearn
committed
reax_list **lists, output_controls *out_control )
Kurt A. O'Hearn
committed
int i, j, pj, target;
int Htop, H_sp_top, btop_i, btop_j, num_bonds, num_hbonds;
Kurt A. O'Hearn
committed
int flag, flag_sp, val_flag, renbr;
Kurt A. O'Hearn
committed
real r_ij, r2, val;
real C12, C34, C56;
real Cln_BOp_s, Cln_BOp_pi, Cln_BOp_pi2;
real BO, BO_s, BO_pi, BO_pi2;
sparse_matrix *H, *H_sp;
Kurt A. O'Hearn
committed
reax_list *far_nbrs, *bonds, *hbonds;
single_body_parameters *sbp_i, *sbp_j;
two_body_parameters *twbp;
far_neighbor_data *nbr_pj;
reax_atom *atom_i, *atom_j;
bond_data *ibond, *jbond;
bond_order_data *bo_ij, *bo_ji;
far_nbrs = lists[FAR_NBRS];
bonds = lists[BONDS];
hbonds = lists[HBONDS];
Kurt A. O'Hearn
committed
H = &workspace->H;
H_sp = &workspace->H_sp;
Kurt A. O'Hearn
committed
btop_i = 0;
btop_j = 0;
Kurt A. O'Hearn
committed
renbr = ((data->step - data->prev_steps) % control->reneighbor) == 0 ? TRUE : FALSE;
atom_i = &system->atoms[i];
type_i = atom_i->type;
Kurt A. O'Hearn
committed
start_i = Start_Index( i, far_nbrs );
end_i = End_Index( i, far_nbrs );
H_sp->start[i] = H_sp_top;
Kurt A. O'Hearn
committed
sbp_i = &system->reax_param.sbp[type_i];
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
if ( control->hbond_cut > 0.0 )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
ihb = sbp_i->p_hbond;
Kurt A. O'Hearn
committed
if ( ihb == H_ATOM )
Kurt A. O'Hearn
committed
{
ihb_top = End_Index( workspace->hbond_index[i], hbonds );
}
else
{
ihb_top = -1;
}
Kurt A. O'Hearn
committed
}
else
{
ihb = NON_H_BONDING_ATOM;
Kurt A. O'Hearn
committed
nbr_pj = &far_nbrs->far_nbr_list[pj];
Kurt A. O'Hearn
committed
flag = FALSE;
flag_sp = FALSE;
#if defined(QMMM)
if ( system->atoms[i].qmmm_mask == TRUE
|| system->atoms[j].qmmm_mask == TRUE )
{
#endif
Kurt A. O'Hearn
committed
/* check if reneighboring step --
* atomic distances just computed via
* Verlet list, so use current distances */
Kurt A. O'Hearn
committed
if ( renbr == TRUE )
Kurt A. O'Hearn
committed
if ( nbr_pj->d <= control->nonb_cut )
Kurt A. O'Hearn
committed
flag = TRUE;
if ( nbr_pj->d <= control->nonb_sp_cut )
{
flag_sp = TRUE;
}
Kurt A. O'Hearn
committed
/* update atomic distances */
Kurt A. O'Hearn
committed
else
atom_j = &system->atoms[j];
Kurt A. O'Hearn
committed
nbr_pj->d = control->compute_atom_distance( &system->box,
Kurt A. O'Hearn
committed
atom_i->x, atom_j->x, atom_i->rel_map,
atom_j->rel_map, nbr_pj->rel_box,
Kurt A. O'Hearn
committed
nbr_pj->dvec );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
if ( nbr_pj->d <= control->nonb_cut )
Kurt A. O'Hearn
committed
flag = TRUE;
if ( nbr_pj->d <= control->nonb_sp_cut )
Kurt A. O'Hearn
committed
{
flag_sp = TRUE;
}
Kurt A. O'Hearn
committed
if ( flag == TRUE )
type_j = system->atoms[j].type;
Kurt A. O'Hearn
committed
sbp_j = &system->reax_param.sbp[type_j];
twbp = &system->reax_param.tbp[type_i][type_j];
Kurt A. O'Hearn
committed
val = Init_Charge_Matrix_Entry( system, control,
workspace, i, j, r_ij, OFF_DIAGONAL );
val_flag = FALSE;
for ( target = H->start[i]; target < Htop; ++target )
{
if ( H->j[target] == j )
{
H->val[target] += val;
val_flag = TRUE;
break;
}
}
if ( val_flag == FALSE )
{
H->j[Htop] = j;
H->val[Htop] = val;
++Htop;
}
/* H_sp matrix entry */
Kurt A. O'Hearn
committed
if ( flag_sp == TRUE )
Kurt A. O'Hearn
committed
val_flag = FALSE;
for ( target = H_sp->start[i]; target < H_sp_top; ++target )
{
if ( H_sp->j[target] == j )
{
H_sp->val[target] += val;
val_flag = TRUE;
break;
}
}
if ( val_flag == FALSE )
{
H_sp->j[H_sp_top] = j;
H_sp->val[H_sp_top] = val;
++H_sp_top;
}
#if defined(QMMM)
if ( system->atoms[i].qmmm_mask == TRUE
&& system->atoms[j].qmmm_mask == TRUE )
{
#endif
/* Only non-dummy atoms can form bonds */
if ( system->atoms[i].is_dummy == FALSE
&& system->atoms[j].is_dummy == FALSE )
Kurt A. O'Hearn
committed
/* hydrogen bond lists */
if ( control->hbond_cut > 0.0
&& (ihb == H_ATOM || ihb == H_BONDING_ATOM)
&& nbr_pj->d <= control->hbond_cut )
jhb = sbp_j->p_hbond;
if ( ihb == H_ATOM && jhb == H_BONDING_ATOM )
{
hbonds->hbond_list[ihb_top].nbr = j;
hbonds->hbond_list[ihb_top].scl = 1;
hbonds->hbond_list[ihb_top].ptr = nbr_pj;
++ihb_top;
++num_hbonds;
}
else if ( ihb == H_BONDING_ATOM && jhb == H_ATOM )
{
jhb_top = End_Index( workspace->hbond_index[j], hbonds );
hbonds->hbond_list[jhb_top].nbr = i;
hbonds->hbond_list[jhb_top].scl = -1;
hbonds->hbond_list[jhb_top].ptr = nbr_pj;
Set_End_Index( workspace->hbond_index[j], jhb_top + 1, hbonds );
++num_hbonds;
}
/* uncorrected bond orders */
if ( nbr_pj->d <= control->bond_cut )
if ( sbp_i->r_s > 0.0 && sbp_j->r_s > 0.0 )
{
C12 = twbp->p_bo1 * POW( r_ij / twbp->r_s, twbp->p_bo2 );
BO_s = (1.0 + control->bo_cut) * EXP( C12 );
}
else
{
C12 = 0.0;
BO_s = 0.0;
}
if ( sbp_i->r_pi > 0.0 && sbp_j->r_pi > 0.0 )
{
C34 = twbp->p_bo3 * POW( r_ij / twbp->r_p, twbp->p_bo4 );
BO_pi = EXP( C34 );
}
else
{
C34 = 0.0;
BO_pi = 0.0;
}
if ( sbp_i->r_pi_pi > 0.0 && sbp_j->r_pi_pi > 0.0 )
{
C56 = twbp->p_bo5 * POW( r_ij / twbp->r_pp, twbp->p_bo6 );
BO_pi2 = EXP( C56 );
}
else
{
C56 = 0.0;
BO_pi2 = 0.0;
}
/* Initially BO values are the uncorrected ones, page 1 */
BO = BO_s + BO_pi + BO_pi2;
if ( BO >= control->bo_cut )
{
num_bonds += 2;
/****** bonds i-j and j-i ******/
ibond = &bonds->bond_list[btop_i];
btop_j = End_Index( j, bonds );
jbond = &bonds->bond_list[btop_j];
ibond->nbr = j;
jbond->nbr = i;
ibond->d = r_ij;