Newer
Older
/*----------------------------------------------------------------------
PuReMD - Purdue ReaxFF Molecular Dynamics Program
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 "reax_types.h"
Kurt A. O'Hearn
committed
#include "init_md.h"
#include "allocate.h"
#include "box.h"
#include "comm_tools.h"
#include "forces.h"
#include "grid.h"
#include "integrate.h"
#include "io_tools.h"
#include "list.h"
#include "lookup.h"
#include "neighbors.h"
#include "random.h"
#include "reset_tools.h"
#include "system_props.h"
#include "tool_box.h"
#include "vector.h"
Kurt A. O'Hearn
committed
#include "reax_init_md.h"
#include "reax_allocate.h"
#include "reax_forces.h"
#include "reax_io_tools.h"
#include "reax_list.h"
#include "reax_lookup.h"
#include "reax_reset_tools.h"
#include "reax_system_props.h"
#include "reax_tool_box.h"
#include "reax_vector.h"
#endif
#if defined(PURE_REAX)
/************************ initialize system ************************/
int Reposition_Atoms( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, mpi_datatypes *mpi_data, char *msg )
Kurt A. O'Hearn
committed
int i;
rvec dx;
Kurt A. O'Hearn
committed
/* fit atoms to periodic box */
if ( control->reposition_atoms == 0 )
Kurt A. O'Hearn
committed
/* put center of mass to center */
else if ( control->reposition_atoms == 1 )
{
rvec_Scale( dx, 0.5, system->big_box.box_norms );
rvec_ScaledAdd( dx, -1., data->xcm );
}
Kurt A. O'Hearn
committed
/* put center of mass to origin */
else if ( control->reposition_atoms == 2 )
{
rvec_Scale( dx, -1., data->xcm );
}
else
{
Kurt A. O'Hearn
committed
strcpy( msg, "[ERROR] reposition_atoms: invalid option" );
return FAILURE;
}
for ( i = 0; i < system->n; ++i )
Kurt A. O'Hearn
committed
{
// Inc_on_T3_Gen( system->my_atoms[i].x, dx, &(system->big_box) );
rvec_Add( system->my_atoms[i].x, dx );
Kurt A. O'Hearn
committed
}
}
void Generate_Initial_Velocities( reax_system *system, real T )
{
int i;
real m, scale, norm;
if ( T <= 0.1 )
{
for ( i = 0; i < system->n; i++ )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Randomize( );
for ( i = 0; i < system->n; i++ )
{
rvec_Random( system->my_atoms[i].v );
norm = rvec_Norm_Sqr( system->my_atoms[i].v );
m = system->reax_param.sbp[ system->my_atoms[i].type ].mass;
scale = SQRT( m * norm / (3.0 * K_B * T) );
Kurt A. O'Hearn
committed
rvec_Scale( system->my_atoms[i].v, 1.0 / scale, system->my_atoms[i].v );
Kurt A. O'Hearn
committed
void Init_System( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, storage *workspace,
Kurt A. O'Hearn
committed
mpi_datatypes *mpi_data )
int i;
reax_atom *atom;
int nrecv[MAX_NBRS];
Setup_New_Grid( system, control, MPI_COMM_WORLD );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Print_Grid( &system->my_grid, stderr );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Bin_My_Atoms( system, &workspace->realloc );
Reorder_My_Atoms( system, workspace );
/* estimate N and total capacity */
Kurt A. O'Hearn
committed
for ( i = 0; i < MAX_NBRS; ++i )
{
nrecv[i] = 0;
}
MPI_Barrier( MPI_COMM_WORLD );
system->N = SendRecv( system, mpi_data, mpi_data->boundary_atom_type, nrecv,
Estimate_Boundary_Atoms, Unpack_Estimate_Message, TRUE );
Kurt A. O'Hearn
committed
system->total_cap = MAX( (int)(system->N * SAFE_ZONE), MIN_CAP );
Bin_Boundary_Atoms( system );
/* estimate numH and Hcap */
system->numH = 0;
Kurt A. O'Hearn
committed
if ( control->hbond_cut > 0.0 )
{
Kurt A. O'Hearn
committed
for ( i = 0; i < system->N; ++i )
Kurt A. O'Hearn
committed
atom = &system->my_atoms[i];
Kurt A. O'Hearn
committed
if ( system->reax_param.sbp[ atom->type ].p_hbond == H_ATOM )
Kurt A. O'Hearn
committed
{
atom->Hindex = system->numH++;
Kurt A. O'Hearn
committed
}
else
{
atom->Hindex = -1;
}
Kurt A. O'Hearn
committed
}
//Tried fix
//system->Hcap = MAX( system->numH * SAFER_ZONE, MIN_CAP );
system->Hcap = MAX( system->n * SAFER_ZONE, MIN_CAP );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
/* list management */
system->far_nbrs = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->far_nbrs" );
system->max_far_nbrs = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->max_far_nbrs" );
system->bonds = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->bonds" );
system->max_bonds = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->max_bonds" );
system->hbonds = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->hbonds" );
system->max_hbonds = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->max_hbonds" );
system->cm_entries = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::system->cm_entries" );
system->max_cm_entries = (int *) smalloc( sizeof(int) * system->total_cap,
"ReAllocate_System::max_cm_entries->max_hbonds" );
fprintf( stderr, "p%d: n=%d local_cap=%d\n",
system->my_rank, system->n, system->local_cap );
fprintf( stderr, "p%d: N=%d total_cap=%d\n",
system->my_rank, system->N, system->total_cap );
fprintf( stderr, "p%d: numH=%d H_cap=%d\n",
system->my_rank, system->numH, system->Hcap );
Compute_Total_Mass( system, data, mpi_data->comm_mesh3D );
Kurt A. O'Hearn
committed
Compute_Center_of_Mass( system, data, mpi_data, mpi_data->comm_mesh3D );
Kurt A. O'Hearn
committed
// if( Reposition_Atoms( system, control, data, mpi_data ) == FAILURE )
Kurt A. O'Hearn
committed
// {
// return FAILURE;
// }
/* initialize velocities so that desired init T can be attained */
if ( !control->restart || (control->restart && control->random_vel) )
Kurt A. O'Hearn
committed
{
Generate_Initial_Velocities( system, control->T_init );
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
Compute_Kinetic_Energy( system, data, mpi_data->comm_mesh3D );
}
/************************ initialize simulation data ************************/
Kurt A. O'Hearn
committed
void Init_Simulation_Data( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data )
Reset_Simulation_Data( data );
if ( !control->restart )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
}
switch ( control->ensemble )
{
case NVE:
data->N_f = 3 * system->bigN;
Kurt A. O'Hearn
committed
control->Evolve = Velocity_Verlet_NVE;
control->virial = 0;
break;
case bNVT:
data->N_f = 3 * system->bigN + 1;
Kurt A. O'Hearn
committed
control->Evolve = Velocity_Verlet_Berendsen_NVT;
Kurt A. O'Hearn
committed
fprintf( stderr, "[WARNING] Nose-Hoover NVT is still under testing.\n" );
Kurt A. O'Hearn
committed
control->Evolve = Velocity_Verlet_Nose_Hoover_NVT_Klein;
control->virial = 0;
if ( !control->restart || (control->restart && control->random_vel) )
{
data->therm.G_xi = control->Tau_T *
(2.0 * data->sys_en.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;
}
break;
Kurt A. O'Hearn
committed
/* Semi-Isotropic NPT */
case sNPT:
Kurt A. O'Hearn
committed
control->Evolve = Velocity_Verlet_Berendsen_NPT;
Kurt A. O'Hearn
committed
/* Isotropic NPT */
case iNPT:
Kurt A. O'Hearn
committed
control->Evolve = Velocity_Verlet_Berendsen_NPT;
Kurt A. O'Hearn
committed
/* Anisotropic NPT */
case NPT:
fprintf( stderr, "[ERROR] p%d: init_simulation_data: option not yet implemented\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
Kurt A. O'Hearn
committed
control->Evolve = Velocity_Verlet_Berendsen_NPT;
control->virial = 1;
/*if( !control->restart ) {
data->therm.G_xi = control->Tau_T *
(2.0 * data->my_en.e_Kin - data->N_f * K_B * control->T );
data->therm.v_xi = data->therm.G_xi * control->dt;
Kurt A. O'Hearn
committed
data->iso_bar.eps = (1.0 / 3.0) * LOG(system->box.volume);
data->inv_W = 1.0 /
( data->N_f * K_B * control->T * SQR(control->Tau_P) );
Compute_Pressure( system, control, data, out_control );
}*/
break;
default:
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] p%d: init_simulation_data: ensemble not recognized\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
Kurt A. O'Hearn
committed
MPI_Barrier( MPI_COMM_WORLD );
if ( system->my_rank == MASTER_NODE )
{
data->timing.start = Get_Time( );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
void Init_System( reax_system *system )
system->big_box.V = 0;
system->big_box.box_norms[0] = 0;
system->big_box.box_norms[1] = 0;
system->big_box.box_norms[2] = 0;
system->local_cap = (int)(system->n * SAFE_ZONE);
system->total_cap = (int)(system->N * SAFE_ZONE);
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
system->far_nbrs = NULL;
system->max_far_nbrs = NULL;
system->bonds = NULL;
system->max_bonds = NULL;
system->hbonds = NULL;
system->max_hbonds = NULL;
system->cm_entries = NULL;
system->max_cm_entries = NULL;
fprintf( stderr, "p%d: local_cap=%d total_cap=%d\n",
system->my_rank, system->local_cap, system->total_cap );
Kurt A. O'Hearn
committed
ReAllocate_System( system, system->local_cap, system->total_cap );
Kurt A. O'Hearn
committed
void Init_Simulation_Data( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data )
Kurt A. O'Hearn
committed
//if( !control->restart )
data->step = data->prev_steps = 0;
}
#endif
/************************ initialize workspace ************************/
/* Initialize Taper params */
void Init_Taper( control_params *control, storage *workspace )
{
real d1, d7;
real swa, swa2, swa3;
real swb, swb2, swb3;
swa = control->nonb_low;
swb = control->nonb_cut;
Kurt A. O'Hearn
committed
if ( FABS( swa ) > 0.01 )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
fprintf( stderr, "[WARNING] non-zero lower Taper-radius cutoff in force field parameters\n" );
Kurt A. O'Hearn
committed
}
Kurt A. O'Hearn
committed
fprintf( stderr, "[ERROR] negative upper Taper-radius cutoff in force field parameters\n" );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
}
else if ( swb < 5 )
Kurt A. O'Hearn
committed
{
Kurt A. O'Hearn
committed
fprintf( stderr, "[WARNING] very low Taper-radius cutoff in force field parameters (%f)\n", swb );
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 );
workspace->Tap[7] = 20.0 / d7;
workspace->Tap[6] = -70.0 * (swa + swb) / d7;
workspace->Tap[5] = 84.0 * (swa2 + 3.0 * swa * swb + swb2) / d7;
workspace->Tap[4] = -35.0 * (swa3 + 9.0 * swa2 * swb + 9.0 * swa * swb2 + swb3 ) / d7;
workspace->Tap[3] = 140.0 * (swa3 * swb + 3.0 * swa2 * swb2 + swa * swb3 ) / d7;
workspace->Tap[2] = -210.0 * (swa3 * swb2 + swa2 * swb3) / d7;
workspace->Tap[1] = 140.0 * swa3 * swb3 / d7;
workspace->Tap[0] = (-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
storage *workspace )
Kurt A. O'Hearn
committed
Allocate_Workspace( system, control, workspace, system->local_cap,
Kurt A. O'Hearn
committed
system->total_cap );
workspace->realloc.far_nbrs = FALSE;
workspace->realloc.cm = FALSE;
workspace->realloc.hbonds = FALSE;
workspace->realloc.bonds = FALSE;
workspace->realloc.thbody = FALSE;
workspace->realloc.gcell_atoms = 0;
Kurt A. O'Hearn
committed
/************** setup communication data structures **************/
Kurt A. O'Hearn
committed
void Init_MPI_Datatypes( reax_system *system, storage *workspace,
Kurt A. O'Hearn
committed
mpi_datatypes *mpi_data )
Kurt A. O'Hearn
committed
int block[11];
int i;
Kurt A. O'Hearn
committed
MPI_Aint disp[11];
MPI_Aint base, size_entry;
MPI_Datatype type[11], temp_type;
mpi_atom sample[2];
boundary_atom b_sample[2];
restart_atom r_sample[2];
rvec rvec_sample[2];
rvec2 rvec2_sample[2];
/* mpi_atom */
Kurt A. O'Hearn
committed
block[0] = 1;
block[1] = 1;
block[2] = 1;
block[3] = 1;
block[4] = 1;
block[5] = MAX_ATOM_NAME_LEN;
Kurt A. O'Hearn
committed
block[6] = 3;
block[7] = 3;
block[8] = 3;
block[9] = 4;
block[10] = 4;
MPI_Get_address( sample, disp );
MPI_Get_address( &sample[0].imprt_id, disp + 1 );
MPI_Get_address( &sample[0].type, disp + 2 );
MPI_Get_address( &sample[0].num_bonds, disp + 3 );
MPI_Get_address( &sample[0].num_hbonds, disp + 4 );
MPI_Get_address( sample[0].name, disp + 5 );
MPI_Get_address( sample[0].x, disp + 6 );
MPI_Get_address( sample[0].v, disp + 7 );
MPI_Get_address( sample[0].f_old, disp + 8 );
MPI_Get_address( sample[0].s, disp + 9 );
MPI_Get_address( sample[0].t, disp + 10 );
base = disp[0];
for ( i = 0; i < 11; ++i )
{
disp[i] -= base;
}
Kurt A. O'Hearn
committed
type[0] = MPI_INT;
type[1] = MPI_INT;
type[2] = MPI_INT;
type[3] = MPI_INT;
type[4] = MPI_INT;
Kurt A. O'Hearn
committed
type[6] = MPI_DOUBLE;
type[7] = MPI_DOUBLE;
type[8] = MPI_DOUBLE;
type[9] = MPI_DOUBLE;
type[10] = MPI_DOUBLE;
MPI_Type_create_struct( 11, block, disp, type, &temp_type );
/* in case of compiler padding, compute difference
* between 2 consecutive struct elements */
MPI_Get_address( sample + 1, &size_entry );
size_entry = MPI_Aint_diff( size_entry, base );
MPI_Type_create_resized( temp_type, 0, size_entry, &mpi_data->mpi_atom_type );
Kurt A. O'Hearn
committed
MPI_Type_commit( &mpi_data->mpi_atom_type );
/* boundary_atom */
Kurt A. O'Hearn
committed
block[0] = 1;
block[1] = 1;
block[2] = 1;
block[3] = 1;
block[4] = 1;
MPI_Get_address( b_sample, disp );
MPI_Get_address( &b_sample[0].imprt_id, disp + 1 );
MPI_Get_address( &b_sample[0].type, disp + 2 );
MPI_Get_address( &b_sample[0].num_bonds, disp + 3 );
MPI_Get_address( &b_sample[0].num_hbonds, disp + 4 );
MPI_Get_address( b_sample[0].x, disp + 5 );
base = disp[0];
for ( i = 0; i < 6; ++i )
{
disp[i] -= base;
}
Kurt A. O'Hearn
committed
type[0] = MPI_INT;
type[1] = MPI_INT;
type[2] = MPI_INT;
type[3] = MPI_INT;
type[4] = MPI_INT;
MPI_Type_create_struct( 6, block, disp, type, &temp_type );
/* in case of compiler padding, compute difference
* between 2 consecutive struct elements */
MPI_Get_address( b_sample + 1, &size_entry );
size_entry = MPI_Aint_diff( size_entry, base );
MPI_Type_create_resized( temp_type, 0, size_entry, &mpi_data->boundary_atom_type );
Kurt A. O'Hearn
committed
MPI_Type_commit( &mpi_data->boundary_atom_type );
MPI_Get_address( &rvec_sample, disp );
base = disp[0];
for ( i = 0; i < 1; ++i )
{
disp[i] -= base;
}
MPI_Type_create_struct( 1, block, disp, type, &temp_type );
/* in case of compiler padding, compute difference
* between 2 consecutive struct elements */
MPI_Get_address( rvec_sample + 1, &size_entry );
size_entry = MPI_Aint_diff( size_entry, base );
MPI_Type_create_resized( temp_type, 0, size_entry, &mpi_data->mpi_rvec );
Kurt A. O'Hearn
committed
MPI_Type_commit( &mpi_data->mpi_rvec );
MPI_Get_address( &rvec2_sample, disp );
base = disp[0];
for ( i = 0; i < 1; ++i )
{
disp[i] -= base;
}
MPI_Type_create_struct( 1, block, disp, type, &temp_type );
/* in case of compiler padding, compute difference
* between 2 consecutive struct elements */
MPI_Get_address( rvec2_sample + 1, &size_entry );
size_entry = MPI_Aint_diff( size_entry, base );
MPI_Type_create_resized( temp_type, 0, size_entry, &mpi_data->mpi_rvec2 );
Kurt A. O'Hearn
committed
MPI_Type_commit( &mpi_data->mpi_rvec2 );
/* restart_atom */
Kurt A. O'Hearn
committed
block[0] = 1;
block[1] = 1 ;
block[2] = MAX_ATOM_NAME_LEN;
Kurt A. O'Hearn
committed
block[3] = 3;
block[4] = 3;
MPI_Get_address( &r_sample, disp );
MPI_Get_address( &r_sample[0].type, disp + 1 );
MPI_Get_address( r_sample[0].name, disp + 2 );
MPI_Get_address( r_sample[0].x, disp + 3 );
MPI_Get_address( r_sample[0].v, disp + 4 );
base = disp[0];
for ( i = 0; i < 5; ++i )
{
disp[i] -= base;
}
Kurt A. O'Hearn
committed
type[0] = MPI_INT;
type[1] = MPI_INT;
Kurt A. O'Hearn
committed
type[3] = MPI_DOUBLE;
type[4] = MPI_DOUBLE;
MPI_Type_create_struct( 5, block, disp, type, &temp_type );
/* in case of compiler padding, compute difference
* between 2 consecutive struct elements */
MPI_Get_address( r_sample + 1, &size_entry );
size_entry = MPI_Aint_diff( size_entry, base );
MPI_Type_create_resized( temp_type, 0, size_entry, &mpi_data->restart_atom_type );
Kurt A. O'Hearn
committed
MPI_Type_commit( &mpi_data->restart_atom_type );
mpi_data->in1_buffer = NULL;
mpi_data->in2_buffer = NULL;
}
/********************** allocate lists *************************/
Kurt A. O'Hearn
committed
void Init_Lists( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, storage *workspace, reax_list **lists,
Kurt A. O'Hearn
committed
mpi_datatypes *mpi_data )
Kurt A. O'Hearn
committed
int ret;
Kurt A. O'Hearn
committed
Estimate_Num_Neighbors( system );
Kurt A. O'Hearn
committed
Make_List( system->total_cap, system->total_far_nbrs, TYP_FAR_NEIGHBOR, lists[FAR_NBRS] );
Init_List_Indices( lists[FAR_NBRS], system->max_far_nbrs );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
ret = Generate_Neighbor_Lists( system, data, workspace, lists );
if ( ret != SUCCESS )
{
fprintf( stderr, "[ERROR] p%d: failed to generate neighbor lists. Terminating...\n", system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
Estimate_Storages( system, control, lists );
Kurt A. O'Hearn
committed
Allocate_Matrix( &workspace->H, system->local_cap, system->total_cm_entries );
Kurt A. O'Hearn
committed
Init_Matrix_Row_Indices( &workspace->H, system->max_cm_entries );
Kurt A. O'Hearn
committed
if ( control->hbond_cut > 0.0 )
Kurt A. O'Hearn
committed
Make_List( system->total_cap, system->total_hbonds, TYP_HBOND, lists[HBONDS] );
Init_List_Indices( lists[HBONDS], system->max_hbonds );
Kurt A. O'Hearn
committed
Make_List( system->total_cap, system->total_bonds, TYP_BOND, lists[BONDS] );
Init_List_Indices( lists[BONDS], system->max_bonds );
Kurt A. O'Hearn
committed
Make_List( system->total_bonds, system->total_thbodies, TYP_THREE_BODY, lists[THREE_BODIES] );
Kurt A. O'Hearn
committed
Make_List( system->total_cap, system->total_bonds * 8, TYP_DDELTA, lists[DDELTAS] );
Make_List( system->total_bonds, system->total_bonds * 50, TYP_DBO, lists[DBOS] );
Kurt A. O'Hearn
committed
void Initialize( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, storage *workspace,
reax_list **lists, output_controls *out_control,
mpi_datatypes *mpi_data )
Kurt A. O'Hearn
committed
Init_MPI_Datatypes( system, workspace, mpi_data );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_System( system, control, data, workspace, mpi_data );
Kurt A. O'Hearn
committed
Init_Simulation_Data( system, control, data );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_Workspace( system, control, workspace );
Kurt A. O'Hearn
committed
Init_Lists( system, control, data, workspace, lists, mpi_data );
Kurt A. O'Hearn
committed
Init_Output_Files( system, control, out_control, mpi_data );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_Lookup_Tables( system, control, workspace->Tap, mpi_data );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
#ifdef TEST_FORCES
// Init_Force_Test_Functions( );
// fprintf( stderr, "p%d: initialized force test functions\n", system->my_rank );
Kurt A. O'Hearn
committed
void Pure_Initialize( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, storage *workspace,
reax_list **lists, output_controls *out_control,
mpi_datatypes *mpi_data )
Kurt A. O'Hearn
committed
Init_Simulation_Data( system, control, data );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_Workspace( system, control, workspace );
Kurt A. O'Hearn
committed
Init_Lists( system, control, data, workspace, lists, mpi_data );
void Initialize( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, storage *workspace,
reax_list **lists, output_controls *out_control,
mpi_datatypes *mpi_data )
Kurt A. O'Hearn
committed
Init_System( system );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_Simulation_Data( system, control, data );
Kurt A. O'Hearn
committed
Init_Workspace( system, control, workspace );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_MPI_Datatypes( system, workspace, mpi_data );
Kurt A. O'Hearn
committed
Init_Lists( system, control, workspace, lists );
Kurt A. O'Hearn
committed
Init_Output_Files( system, control, out_control, mpi_data );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
Init_Lookup_Tables( system, control, workspace->Tap, mpi_data );