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 "init_md.h"
#include "allocate.h"
#include "box.h"
#include "forces.h"
#include "grid.h"
#include "integrate.h"
#include "neighbors.h"
#include "list.h"
#include "lookup.h"
#include "reset_utils.h"
#include "system_props.h"
#include "tool_box.h"
#include "vector.h"
void Generate_Initial_Velocities( reax_system *system, real T )
{
Kurt A. O'Hearn
committed
for ( i = 0; i < system->N; i++ )
{
Kurt A. O'Hearn
committed
}
else
{
for ( i = 0; i < system->N; i++ )
{
rvec_Random( system->atoms[i].v );
norm = rvec_Norm_Sqr( system->atoms[i].v );
scale = SQRT( system->reaxprm.sbp[ system->atoms[i].type ].mass *
norm / (3.0 * K_B * T) );
rvec_Scale( system->atoms[i].v, 1.0 / scale, system->atoms[i].v );
/*fprintf( stderr, "v = %f %f %f\n",
system->atoms[i].v[0],system->atoms[i].v[1],system->atoms[i].v[2]);
fprintf( stderr, "scale = %f\n", scale );
fprintf( stderr, "v = %f %f %f\n",
system->atoms[i].v[0],system->atoms[i].v[1],system->atoms[i].v[2]);*/
}
}
void Init_System( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
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
115
116
117
118
119
120
121
122
Compute_Total_Mass( system, data );
Compute_Center_of_Mass( system, data, stderr );
/* reposition atoms */
// just fit the atoms to the periodic box
if ( control->reposition_atoms == 0 )
{
rvec_MakeZero( dx );
}
// put the center of mass to the center of the box
else if ( control->reposition_atoms == 1 )
{
rvec_Scale( dx, 0.5, system->box.box_norms );
rvec_ScaledAdd( dx, -1., data->xcm );
}
// put the center of mass to the origin
else if ( control->reposition_atoms == 2 )
{
rvec_Scale( dx, -1., data->xcm );
}
else
{
fprintf( stderr, "UNKNOWN OPTION: reposition_atoms. Terminating...\n" );
exit( UNKNOWN_OPTION );
}
for ( i = 0; i < system->N; ++i )
{
Inc_on_T3( system->atoms[i].x, dx, &(system->box) );
/*fprintf( stderr, "%6d%2d%8.3f%8.3f%8.3f\n",
i, system->atoms[i].type,
system->atoms[i].x[0], system->atoms[i].x[1], system->atoms[i].x[2] );*/
}
/* Initialize velocities so that desired init T can be attained */
if ( !control->restart || (control->restart && control->random_vel) )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
void Init_Simulation_Data( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, output_controls *out_control,
evolve_function *Evolve )
Reset_Simulation_Data( data );
if ( !control->restart )
Kurt A. O'Hearn
committed
{
data->step = 0;
data->prev_steps = 0;
}
switch ( control->ensemble )
{
case NVE:
data->N_f = 3 * system->N;
*Evolve = Velocity_Verlet_NVE;
break;
case NVT:
data->N_f = 3 * system->N + 1;
//control->Tau_T = 100 * data->N_f * K_B * control->T_final;
if ( !control->restart || (control->restart && control->random_vel) )
{
data->therm.G_xi = control->Tau_T * (2.0 * data->E_Kin -
data->N_f * K_B * control->T );
data->therm.v_xi = data->therm.G_xi * control->dt;
data->therm.v_xi_old = 0;
data->therm.xi = 0;
fprintf( stderr, "init_md: G_xi=%f Tau_T=%f E_kin=%f N_f=%f v_xi=%f\n",
data->therm.G_xi, control->Tau_T, data->E_Kin,
data->N_f, data->therm.v_xi );
}
*Evolve = Velocity_Verlet_Nose_Hoover_NVT_Klein;
break;
case NPT: // Anisotropic NPT
fprintf( stderr, "THIS OPTION IS NOT YET IMPLEMENTED! TERMINATING...\n" );
exit( UNKNOWN_OPTION );
data->N_f = 3 * system->N + 9;
if ( !control->restart )
{
data->therm.G_xi = control->Tau_T * (2.0 * data->E_Kin -
Kurt A. O'Hearn
committed
data->N_f * K_B * control->T);
Kurt A. O'Hearn
committed
data->iso_bar.eps = 1.0 / 3.0 * LOG( system->box.volume );
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
//data->inv_W = 1. / (data->N_f*K_B*control->T*SQR(control->Tau_P));
//Compute_Pressure( system, data, workspace );
}
*Evolve = Velocity_Verlet_Berendsen_Isotropic_NPT;
break;
case sNPT: // Semi-Isotropic NPT
data->N_f = 3 * system->N + 4;
*Evolve = Velocity_Verlet_Berendsen_SemiIsotropic_NPT;
break;
case iNPT: // Isotropic NPT
data->N_f = 3 * system->N + 2;
*Evolve = Velocity_Verlet_Berendsen_Isotropic_NPT;
break;
case bNVT:
data->N_f = 3 * system->N + 1;
*Evolve = Velocity_Verlet_Berendsen_NVT;
fprintf (stderr, " Initializing Velocity_Verlet_Berendsen_NVT .... \n");
break;
default:
break;
Compute_Kinetic_Energy( system, data );
/* init timing info */
data->timing.start = Get_Time( );
data->timing.total = data->timing.start;
data->timing.nbrs = 0;
data->timing.init_forces = 0;
data->timing.bonded = 0;
data->timing.nonb = 0;
Kurt A. O'Hearn
committed
data->timing.cm = ZERO;
data->timing.cm_sort_mat_rows = ZERO;
data->timing.cm_solver_pre_comp = ZERO;
data->timing.cm_solver_pre_app = ZERO;
data->timing.cm_solver_iters = 0;
data->timing.cm_solver_spmv = ZERO;
data->timing.cm_solver_vector_ops = ZERO;
data->timing.cm_solver_orthog = ZERO;
data->timing.cm_solver_tri_solve = ZERO;
Kurt A. O'Hearn
committed
/* Initialize Taper params */
void Init_Taper( control_params *control )
{
real d1, d7;
real swa, swa2, swa3;
real swb, swb2, swb3;
swa = control->r_low;
swb = control->r_cut;
Kurt A. O'Hearn
committed
if ( FABS( swa ) > 0.01 )
{
Kurt A. O'Hearn
committed
fprintf( stderr, "Warning: non-zero value for lower Taper-radius cutoff\n" );
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
if ( swb < 0.0 )
Kurt A. O'Hearn
committed
{
fprintf( stderr, "Negative value for upper Taper-radius cutoff\n" );
exit( INVALID_INPUT );
}
Kurt A. O'Hearn
committed
else if ( swb < 5.0 )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
fprintf( stderr, "Warning: low value for upper Taper-radius cutoff:%f\n", swb );
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
d1 = swb - swa;
d7 = POW( d1, 7.0 );
swa2 = SQR( swa );
swa3 = CUBE( swa );
swb2 = SQR( swb );
swb3 = CUBE( swb );
control->Tap7 = 20.0 / d7;
control->Tap6 = -70.0 * (swa + swb) / d7;
control->Tap5 = 84.0 * (swa2 + 3.0 * swa * swb + swb2) / d7;
control->Tap4 = -35.0 * (swa3 + 9.0 * swa2 * swb + 9.0 * swa * swb2 + swb3 ) / d7;
control->Tap3 = 140.0 * (swa3 * swb + 3.0 * swa2 * swb2 + swa * swb3 ) / d7;
control->Tap2 = -210.0 * (swa3 * swb2 + swa2 * swb3) / d7;
control->Tap1 = 140.0 * swa3 * swb3 / d7;
control->Tap0 = (-35.0 * swa3 * swb2 * swb2 + 21.0 * swa2 * swb3 * swb2 +
Kurt A. O'Hearn
committed
7.0 * swa * swb3 * swb3 + swb3 * swb3 * swb ) / d7;
Kurt A. O'Hearn
committed
}
void Init_Workspace( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
static_storage *workspace )
{
int i;
/* Allocate space for hydrogen bond list */
workspace->hbond_index = (int *) smalloc( system->N * sizeof( int ),
"Init_Workspace::workspace->hbond_index" );
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
workspace->total_bond_order = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->bond_order" );
workspace->Deltap = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Deltap" );
workspace->Deltap_boc = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Deltap_boc" );
workspace->dDeltap_self = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->dDeltap_self" );
workspace->Delta = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Delta" );
workspace->Delta_lp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Delta_lp" );
workspace->Delta_lp_temp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Delta_lp_temp" );
workspace->dDelta_lp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->dDelta_lp" );
workspace->dDelta_lp_temp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->dDelta_lp_temp" );
workspace->Delta_e = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Delta_e" );
workspace->Delta_boc = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Delta_boc" );
workspace->nlp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->nlp" );
workspace->nlp_temp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->nlp_temp" );
workspace->Clp = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->Clp" );
workspace->CdDelta = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->CdDelta" );
workspace->vlpex = (real *) smalloc( system->N * sizeof( real ),
"Init_Workspace::workspace->vlpex" );
Kurt A. O'Hearn
committed
/* charge method storage */
switch ( control->charge_method )
{
case QEQ_CM:
system->N_cm = system->N;
break;
Kurt A. O'Hearn
committed
system->N_cm = system->N + 1;
break;
case ACKS2_CM:
system->N_cm = 2 * system->N + 2;
break;
default:
fprintf( stderr, "Unknown charge method type. Terminating...\n" );
exit( INVALID_INPUT );
break;
}
Kurt A. O'Hearn
committed
/* sparse matrices */
Kurt A. O'Hearn
committed
workspace->H = NULL;
Kurt A. O'Hearn
committed
workspace->H_full = NULL;
Kurt A. O'Hearn
committed
workspace->H_sp = NULL;
Kurt A. O'Hearn
committed
workspace->H_p = NULL;
workspace->H_spar_patt = NULL;
workspace->H_app_inv = NULL;
Kurt A. O'Hearn
committed
workspace->L = NULL;
Kurt A. O'Hearn
committed
workspace->U = NULL;
workspace->Hdia_inv = NULL;
Kurt A. O'Hearn
committed
if ( control->cm_solver_pre_comp_type == ICHOLT_PC ||
control->cm_solver_pre_comp_type == ILUT_PAR_PC )
{
workspace->droptol = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->droptol" );
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
//TODO: check if unused
//workspace->w = (real *) scalloc( cm_lin_sys_size, sizeof( real ),
//"Init_Workspace::workspace->droptol" );
Kurt A. O'Hearn
committed
//TODO: check if unused
workspace->b = (real *) scalloc( system->N_cm * 2, sizeof( real ),
"Init_Workspace::workspace->b" );
workspace->b_s = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->b_s" );
workspace->b_t = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->b_t" );
workspace->b_prc = (real *) scalloc( system->N_cm * 2, sizeof( real ),
"Init_Workspace::workspace->b_prc" );
workspace->b_prm = (real *) scalloc( system->N_cm * 2, sizeof( real ),
"Init_Workspace::workspace->b_prm" );
workspace->s = (real**) scalloc( 5, sizeof( real* ),
"Init_Workspace::workspace->s" );
workspace->t = (real**) scalloc( 5, sizeof( real* ),
"Init_Workspace::workspace->t" );
workspace->s[i] = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->s[i]" );
workspace->t[i] = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->t[i]" );
Kurt A. O'Hearn
committed
switch ( control->charge_method )
Kurt A. O'Hearn
committed
case QEQ_CM:
for ( i = 0; i < system->N; ++i )
{
workspace->b_s[i] = -system->reaxprm.sbp[ system->atoms[i].type ].chi;
workspace->b_t[i] = -1.0;
//TODO: check if unused (redundant)
workspace->b[i] = -system->reaxprm.sbp[ system->atoms[i].type ].chi;
workspace->b[i + system->N] = -1.0;
}
break;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
for ( i = 0; i < system->N; ++i )
{
workspace->b_s[i] = -system->reaxprm.sbp[ system->atoms[i].type ].chi;
//TODO: check if unused (redundant)
workspace->b[i] = -system->reaxprm.sbp[ system->atoms[i].type ].chi;
}
workspace->b_s[system->N] = control->cm_q_net;
workspace->b[system->N] = control->cm_q_net;
break;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
case ACKS2_CM:
Kurt A. O'Hearn
committed
for ( i = 0; i < system->N; ++i )
{
workspace->b_s[i] = -system->reaxprm.sbp[ system->atoms[i].type ].chi;
//TODO: check if unused (redundant)
workspace->b[i] = -system->reaxprm.sbp[ system->atoms[i].type ].chi;
}
workspace->b_s[system->N] = control->cm_q_net;
workspace->b[system->N] = control->cm_q_net;
for ( i = system->N + 1; i < system->N_cm; ++i )
{
workspace->b_s[i] = 0.0;
//TODO: check if unused (redundant)
workspace->b[i] = 0.0;
}
Kurt A. O'Hearn
committed
break;
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
default:
fprintf( stderr, "Unknown charge method type. Terminating...\n" );
exit( INVALID_INPUT );
break;
Kurt A. O'Hearn
committed
switch ( control->cm_solver_type )
Kurt A. O'Hearn
committed
/* GMRES storage */
case GMRES_S:
case GMRES_H_S:
workspace->y = (real *) scalloc( control->cm_solver_restart + 1, sizeof( real ),
"Init_Workspace::workspace->y" );
workspace->z = (real *) scalloc( control->cm_solver_restart + 1, sizeof( real ),
"Init_Workspace::workspace->z" );
workspace->g = (real *) scalloc( control->cm_solver_restart + 1, sizeof( real ),
"Init_Workspace::workspace->g" );
workspace->h = (real **) scalloc( control->cm_solver_restart + 1, sizeof( real*),
"Init_Workspace::workspace->h" );
workspace->hs = (real *) scalloc( control->cm_solver_restart + 1, sizeof( real ),
"Init_Workspace::workspace->hs" );
workspace->hc = (real *) scalloc( control->cm_solver_restart + 1, sizeof( real ),
"Init_Workspace::workspace->hc" );
workspace->rn = (real **) scalloc( control->cm_solver_restart + 1, sizeof( real*),
"Init_Workspace::workspace->rn" );
workspace->v = (real **) scalloc( control->cm_solver_restart + 1, sizeof( real*),
"Init_Workspace::workspace->v" );
Kurt A. O'Hearn
committed
for ( i = 0; i < control->cm_solver_restart + 1; ++i )
Kurt A. O'Hearn
committed
{
workspace->h[i] = (real *) scalloc( control->cm_solver_restart + 1, sizeof( real ),
"Init_Workspace::workspace->h[i]" );
workspace->rn[i] = (real *) scalloc( system->N_cm * 2, sizeof( real ),
"Init_Workspace::workspace->rn[i]" );
workspace->v[i] = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->v[i]" );
Kurt A. O'Hearn
committed
}
workspace->r = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->r" );
workspace->d = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->d" );
workspace->q = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->q" );
workspace->p = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->p" );
Kurt A. O'Hearn
committed
break;
/* CG storage */
case CG_S:
workspace->r = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->r" );
workspace->d = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->d" );
workspace->q = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->q" );
workspace->p = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->p" );
Kurt A. O'Hearn
committed
break;
case SDM_S:
workspace->r = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->r" );
workspace->d = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->d" );
workspace->q = (real *) scalloc( system->N_cm, sizeof( real ),
"Init_Workspace::workspace->q" );
Kurt A. O'Hearn
committed
break;
default:
fprintf( stderr, "Unknown charge method linear solver type. Terminating...\n" );
exit( INVALID_INPUT );
break;
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
/* SpMV related */
#ifdef _OPENMP
workspace->b_local = (real*) smalloc( control->num_threads * system->N_cm * sizeof(real),
"Init_Workspace::b_local" );
#endif
/* level scheduling related */
workspace->levels_L = 1;
workspace->levels_U = 1;
if ( control->cm_solver_pre_app_type == TRI_SOLVE_LEVEL_SCHED_PA ||
control->cm_solver_pre_app_type == TRI_SOLVE_GC_PA )
{
workspace->row_levels_L = (unsigned int*) smalloc( system->N_cm * sizeof(unsigned int),
"Init_Workspace::row_levels_L" );
workspace->level_rows_L = (unsigned int*) smalloc( system->N_cm * sizeof(unsigned int),
"Init_Workspace::level_rows_L" );
workspace->level_rows_cnt_L = (unsigned int*) smalloc( (system->N_cm + 1) * sizeof(unsigned int),
"Init_Workspace::level_rows_cnt_L" );
workspace->row_levels_U = (unsigned int*) smalloc( system->N_cm * sizeof(unsigned int),
"Init_Workspace::row_levels_U" );
workspace->level_rows_U = (unsigned int*) smalloc( system->N_cm * sizeof(unsigned int),
"Init_Workspace::level_rows_U" );
workspace->level_rows_cnt_U = (unsigned int*) smalloc( (system->N_cm + 1) * sizeof(unsigned int),
"Init_Workspace::level_rows_cnt_U" );
workspace->top = (unsigned int*) smalloc( (system->N_cm + 1) * sizeof(unsigned int),
"Init_Workspace::top" );
}
else
{
workspace->row_levels_L = NULL;
workspace->level_rows_L = NULL;
workspace->level_rows_cnt_L = NULL;
workspace->row_levels_U = NULL;
workspace->level_rows_U = NULL;
workspace->level_rows_cnt_U = NULL;
workspace->top = NULL;
}
Kurt A. O'Hearn
committed
/* graph coloring related */
workspace->recolor_cnt = 0;
Kurt A. O'Hearn
committed
if ( control->cm_solver_pre_app_type == TRI_SOLVE_GC_PA )
{
workspace->color = (unsigned int*) smalloc( sizeof(unsigned int) * system->N_cm,
"Init_Workspace::color" );
Kurt A. O'Hearn
committed
workspace->to_color = (unsigned int*) smalloc( sizeof(unsigned int) * system->N_cm,
"Init_Workspace::to_color" );
Kurt A. O'Hearn
committed
workspace->conflict = (unsigned int*) smalloc( sizeof(unsigned int) * system->N_cm,
"setup_graph_coloring::conflict" );
Kurt A. O'Hearn
committed
workspace->conflict_cnt = (unsigned int*) smalloc( sizeof(unsigned int) * (control->num_threads + 1),
"Init_Workspace::conflict_cnt" );
Kurt A. O'Hearn
committed
workspace->recolor = (unsigned int*) smalloc( sizeof(unsigned int) * system->N_cm,
"Init_Workspace::recolor" );
Kurt A. O'Hearn
committed
workspace->color_top = (unsigned int*) smalloc( sizeof(unsigned int) * (system->N_cm + 1),
"Init_Workspace::color_top" );
Kurt A. O'Hearn
committed
workspace->permuted_row_col = (unsigned int*) smalloc( sizeof(unsigned int) * system->N_cm,
"Init_Workspace::premuted_row_col" );
Kurt A. O'Hearn
committed
workspace->permuted_row_col_inv = (unsigned int*) smalloc( sizeof(unsigned int) * system->N_cm,
"Init_Workspace::premuted_row_col_inv" );
Kurt A. O'Hearn
committed
workspace->y_p = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::y_p" );
Kurt A. O'Hearn
committed
workspace->x_p = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::x_p" );
}
else
{
workspace->color = NULL;
workspace->to_color = NULL;
workspace->conflict = NULL;
workspace->conflict_cnt = NULL;
workspace->recolor = NULL;
workspace->color_top = NULL;
workspace->permuted_row_col = NULL;
workspace->permuted_row_col_inv = NULL;
workspace->y_p = NULL;
workspace->x_p = NULL;
}
/* Jacobi iteration related */
if ( control->cm_solver_pre_app_type == JACOBI_ITER_PA )
{
workspace->Dinv_L = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::Dinv_L" );
workspace->Dinv_U = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::Dinv_U" );
workspace->Dinv_b = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::Dinv_b" );
workspace->rp = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::rp" );
workspace->rp2 = (real*) smalloc( sizeof(real) * system->N_cm,
"Init_Workspace::rp2" );
}
else
{
workspace->Dinv_L = NULL;
workspace->Dinv_U = NULL;
workspace->Dinv_b = NULL;
workspace->rp = NULL;
workspace->rp2 = NULL;
}
workspace->a = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->a" );
workspace->f_old = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_old" );
workspace->v_const = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->v_const" );
Kurt A. O'Hearn
committed
#ifdef _OPENMP
workspace->f_local = (rvec *) smalloc( control->num_threads * system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_local" );
Kurt A. O'Hearn
committed
#endif
/* storage for analysis */
if ( control->molec_anal || control->diffusion_coef )
workspace->mark = (int *) scalloc( system->N, sizeof(int),
"Init_Workspace::workspace->mark" );
workspace->old_mark = (int *) scalloc( system->N, sizeof(int),
"Init_Workspace::workspace->old_mark" );
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
{
workspace->x_old = (rvec *) scalloc( system->N, sizeof( rvec ),
"Init_Workspace::workspace->x_old" );
Kurt A. O'Hearn
committed
}
else
{
workspace->x_old = NULL;
}
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
workspace->dDelta = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->dDelta" );
workspace->f_ele = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_ele" );
workspace->f_vdw = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_vdw" );
workspace->f_bo = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_bo" );
workspace->f_be = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_be" );
workspace->f_lp = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_lp" );
workspace->f_ov = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_ov" );
workspace->f_un = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_un" );
workspace->f_ang = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_ang" );
workspace->f_coa = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_coa" );
workspace->f_pen = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_pen" );
workspace->f_hb = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_hb" );
workspace->f_tor = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_tor" );
workspace->f_con = (rvec *) smalloc( system->N * sizeof( rvec ),
"Init_Workspace::workspace->f_con" );
workspace->realloc.num_far = -1;
workspace->realloc.Htop = -1;
workspace->realloc.hbonds = -1;
workspace->realloc.bonds = -1;
workspace->realloc.num_3body = -1;
workspace->realloc.gcell_atoms = -1;
Kurt A. O'Hearn
committed
/* Initialize Taper function */
Init_Taper( control );
void Init_Lists( 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, num_nbrs, num_bonds, num_3body, Htop, max_nnz;
Kurt A. O'Hearn
committed
#if defined(DEBUG_FOCUS)
int num_hbonds;
#endif
num_nbrs = Estimate_NumNeighbors( system, control, workspace, lists );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Make_List( system->N, num_nbrs, TYP_FAR_NEIGHBOR, &(*lists)[FAR_NBRS] );
Kurt A. O'Hearn
committed
fprintf( stderr, "memory allocated: far_nbrs = %ldMB\n",
num_nbrs * sizeof(far_neighbor_data) / (1024 * 1024) );
Generate_Neighbor_Lists(system, control, data, workspace, lists, out_control);
Htop = 0;
hb_top = (int*) scalloc( system->N, sizeof(int),
"Init_Lists::hb_top" );
bond_top = (int*) scalloc( system->N, sizeof(int),
"Init_Lists::bond_top" );
Kurt A. O'Hearn
committed
Estimate_Storage_Sizes( system, control, lists, &Htop,
hb_top, bond_top, &num_3body );
Kurt A. O'Hearn
committed
num_3body = MAX( num_3body, MIN_BONDS );
Kurt A. O'Hearn
committed
switch ( control->charge_method )
{
case QEQ_CM:
max_nnz = Htop;
break;
Kurt A. O'Hearn
committed
max_nnz = Htop + system->N_cm;
break;
case ACKS2_CM:
max_nnz = 2 * Htop + 3 * system->N + 2;
break;
default:
max_nnz = Htop;
break;
}
if ( Allocate_Matrix( &(workspace->H), system->N_cm, max_nnz ) == FAILURE )
{
fprintf( stderr, "Not enough space for init matrices. Terminating...\n" );
exit( INSUFFICIENT_MEMORY );
}
/* TODO: better estimate for H_sp?
* If so, need to refactor Estimate_Storage_Sizes
* to use various cut-off distances as parameters
* (non-bonded, hydrogen, 3body, etc.) */
Kurt A. O'Hearn
committed
if ( Allocate_Matrix( &(workspace->H_sp), system->N_cm, max_nnz ) == FAILURE )
{
fprintf( stderr, "Not enough space for init matrices. Terminating...\n" );
exit( INSUFFICIENT_MEMORY );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "estimated storage - Htop: %d\n", Htop );
fprintf( stderr, "memory allocated: H = %ldMB\n",
Kurt A. O'Hearn
committed
Htop * sizeof(sparse_matrix_entry) / (1024 * 1024) );
Kurt A. O'Hearn
committed
if ( control->hb_cut > 0.0 )
{
/* init H indexes */
for ( i = 0; i < system->N; ++i )
Kurt A. O'Hearn
committed
{
// H atom
if ( system->reaxprm.sbp[ system->atoms[i].type ].p_hbond == 1 )
{
Kurt A. O'Hearn
committed
}
else
{
workspace->hbond_index[i] = -1;
}
}
Kurt A. O'Hearn
committed
if ( workspace->num_H == 0 )
{
control->hb_cut = 0.0;
}
else
{
Allocate_HBond_List( system->N, workspace->num_H, workspace->hbond_index,
Kurt A. O'Hearn
committed
hb_top, &(*lists)[HBONDS] );
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
num_hbonds = hb_top[system->N - 1];
fprintf( stderr, "estimated storage - num_hbonds: %d\n", num_hbonds );
fprintf( stderr, "memory allocated: hbonds = %ldMB\n",
num_hbonds * sizeof(hbond_data) / (1024 * 1024) );
Kurt A. O'Hearn
committed
Allocate_Bond_List( system->N, bond_top, &(*lists)[BONDS] );
Kurt A. O'Hearn
committed
fprintf( stderr, "estimated storage - num_bonds: %d\n", num_bonds );
fprintf( stderr, "memory allocated: bonds = %ldMB\n",
num_bonds * sizeof(bond_data) / (1024 * 1024) );
Kurt A. O'Hearn
committed
Make_List( num_bonds, num_3body, TYP_THREE_BODY, &(*lists)[THREE_BODIES] );
Kurt A. O'Hearn
committed
fprintf( stderr, "estimated storage - num_3body: %d\n", num_3body );
fprintf( stderr, "memory allocated: 3-body = %ldMB\n",
num_3body * sizeof(three_body_interaction_data) / (1024 * 1024) );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Make_List( system->N, num_bonds * 8, TYP_DDELTA, &(*lists)[DDELTA] );
Kurt A. O'Hearn
committed
Make_List( num_bonds, num_bonds * MAX_BONDS * 3, TYP_DBO, &(*lists)[DBO] );
sfree( hb_top, "Init_Lists::hb_top" );
sfree( bond_top, "Init_Lists::bond_top" );
Kurt A. O'Hearn
committed
void Init_Out_Controls( reax_system *system, control_params *control,
static_storage *workspace, output_controls *out_control )
#define TEMP_SIZE (1000)
char temp[TEMP_SIZE];
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
/* Init trajectory file */
if ( out_control->write_steps > 0 )
{
strcpy( temp, control->sim_name );
strcat( temp, ".trj" );
out_control->trj = fopen( temp, "w" );
out_control->write_header( system, control, workspace, out_control );
}
if ( out_control->energy_update_freq > 0 )
{
/* Init out file */
strcpy( temp, control->sim_name );
strcat( temp, ".out" );
out_control->out = fopen( temp, "w" );
fprintf( out_control->out, "%-6s%16s%16s%16s%11s%11s%13s%13s%13s\n",
"step", "total energy", "poten. energy", "kin. energy",
"temp.", "target", "volume", "press.", "target" );
fflush( out_control->out );
/* Init potentials file */
strcpy( temp, control->sim_name );
strcat( temp, ".pot" );
out_control->pot = fopen( temp, "w" );
fprintf( out_control->pot,
"%-6s%13s%13s%13s%13s%13s%13s%13s%13s%13s%13s%13s\n",
"step", "ebond", "eatom", "elp", "eang", "ecoa", "ehb",
"etor", "econj", "evdw", "ecoul", "epol" );
fflush( out_control->pot );
/* Init log file */
strcpy( temp, control->sim_name );
strcat( temp, ".log" );
out_control->log = fopen( temp, "w" );
Kurt A. O'Hearn
committed
fprintf( out_control->log, "%-6s %10s %10s %10s %10s %10s %10s %10s %10s %10s %10s %10s %10s %10s %10s\n",
Kurt A. O'Hearn
committed
"nonbonded", "CM", "CM Sort", "S iters", "Pre Comp", "Pre App",
Kurt A. O'Hearn
committed
"S spmv", "S vec ops", "S orthog", "S tsolve" );
}
/* Init pressure file */
if ( control->ensemble == NPT ||
control->ensemble == iNPT ||
control->ensemble == sNPT )
{
strcpy( temp, control->sim_name );
strcat( temp, ".prs" );
out_control->prs = fopen( temp, "w" );
fprintf( out_control->prs, "%-6s%13s%13s%13s%13s%13s%13s%13s%13s\n",
"step", "norm_x", "norm_y", "norm_z",
"press_x", "press_y", "press_z", "target_p", "volume" );
fflush( out_control->prs );
}
/* Init molecular analysis file */
if ( control->molec_anal )
{
Kurt A. O'Hearn
committed
snprintf( temp, TEMP_SIZE, "%.*s.mol", TEMP_SIZE - 5, control->sim_name );
out_control->mol = fopen( temp, "w" );
if ( control->num_ignored )
{
Kurt A. O'Hearn
committed
snprintf( temp, TEMP_SIZE, "%.*s.ign", TEMP_SIZE - 5, control->sim_name );
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
out_control->ign = fopen( temp, "w" );
}
}
/* Init electric dipole moment analysis file */
if ( control->dipole_anal )
{
strcpy( temp, control->sim_name );
strcat( temp, ".dpl" );
out_control->dpl = fopen( temp, "w" );
fprintf( out_control->dpl,
"Step Molecule Count Avg. Dipole Moment Norm\n" );
fflush( out_control->dpl );
}
/* Init diffusion coef analysis file */
if ( control->diffusion_coef )
{
strcpy( temp, control->sim_name );
strcat( temp, ".drft" );
out_control->drft = fopen( temp, "w" );
fprintf( out_control->drft, "Step Type Count Avg Squared Disp\n" );
fflush( out_control->drft );
}
#ifdef TEST_ENERGY
/* open bond energy file */
strcat( temp, ".ebond" );
out_control->ebond = fopen( temp, "w" );
strcat( temp, ".elp" );
out_control->elp = fopen( temp, "w" );
/* open overcoordination energy file */
strcpy( temp, control->sim_name );
strcat( temp, ".eov" );
out_control->eov = fopen( temp, "w" );
/* open undercoordination energy file */
strcat( temp, ".eun" );
out_control->eun = fopen( temp, "w" );
/* open angle energy file */
strcat( temp, ".eval" );
out_control->eval = fopen( temp, "w" );
/* open penalty energy file */
strcat( temp, ".epen" );
out_control->epen = fopen( temp, "w" );
/* open coalition energy file */
strcat( temp, ".ecoa" );
out_control->ecoa = fopen( temp, "w" );
/* open hydrogen bond energy file */
strcat( temp, ".ehb" );
out_control->ehb = fopen( temp, "w" );
/* open torsion energy file */
strcpy( temp, control->sim_name );
strcat( temp, ".etor" );
out_control->etor = fopen( temp, "w" );
/* open conjugation energy file */
strcpy( temp, control->sim_name );
strcat( temp, ".econ" );
out_control->econ = fopen( temp, "w" );
/* open vdWaals energy file */
strcpy( temp, control->sim_name );
strcat( temp, ".evdw" );
out_control->evdw = fopen( temp, "w" );
/* open coulomb energy file */
strcpy( temp, control->sim_name );
strcat( temp, ".ecou" );
out_control->ecou = fopen( temp, "w" );
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* open bond orders file */
strcpy( temp, control->sim_name );
strcat( temp, ".fbo" );
out_control->fbo = fopen( temp, "w" );
/* open bond orders derivatives file */
strcpy( temp, control->sim_name );
strcat( temp, ".fdbo" );
out_control->fdbo = fopen( temp, "w" );
/* open bond forces file */
strcpy( temp, control->sim_name );
strcat( temp, ".fbond" );
out_control->fbond = fopen( temp, "w" );
/* open lone-pair forces file */
strcpy( temp, control->sim_name );
strcat( temp, ".flp" );
out_control->flp = fopen( temp, "w" );
/* open overcoordination forces file */
strcpy( temp, control->sim_name );
strcat( temp, ".fatom" );
out_control->fatom = fopen( temp, "w" );
/* open angle forces file */