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;
/* reposition atoms */
if ( control->reposition_atoms == 0 ) //fit atoms to periodic box
{
rvec_MakeZero( dx );
}
else if ( control->reposition_atoms == 1 ) //put center of mass to center
{
rvec_Scale( dx, 0.5, system->big_box.box_norms );
rvec_ScaledAdd( dx, -1., data->xcm );
}
else if ( control->reposition_atoms == 2 ) //put center of mass to origin
{
rvec_Scale( dx, -1., data->xcm );
}
else
{
strcpy( msg, "reposition_atoms: invalid option" );
return FAILURE;
}
for ( i = 0; i < system->n; ++i )
// Inc_on_T3_Gen( system->my_atoms[i].x, dx, &(system->big_box) );
rvec_Add( system->my_atoms[i].x, dx );
return SUCCESS;
}
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) );
rvec_Scale( system->my_atoms[i].v, 1. / scale, system->my_atoms[i].v );
// fprintf( stderr, "v = %f %f %f\n",
// system->my_atoms[i].v[0],
// system->my_atoms[i].v[1],
// system->my_atoms[i].v[2] );
//
// fprintf( stderr, "scale = %f\n", scale );
// fprintf( stderr, "v = %f %f %f\n",
// system->my_atoms[i].v[0],
// system->my_atoms[i].v[1],
// system->my_atoms[i].v[2] );
int Init_System( reax_system *system, control_params *control,
Kurt A. O'Hearn
committed
simulation_data *data, storage *workspace,
mpi_datatypes *mpi_data, char *msg )
int i;
reax_atom *atom;
int nrecv[MAX_NBRS];
Setup_New_Grid( system, control, MPI_COMM_WORLD );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d GRID:\n", system->my_rank );
Print_Grid( &(system->my_grid), stderr );
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 );
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 )
{
for ( i = 0; i < system->n; ++i )
{
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 );
// Sudhir-style below
/*
Kurt A. O'Hearn
committed
if ( control->hbond_cut > 0.0 )
for ( i = 0; i < system->n; ++i )
{
atom = &(system->my_atoms[i]);
if ( system->reax_param.sbp[ atom->type ].p_hbond == H_ATOM )
atom->Hindex = system->numH++;
else atom->Hindex = -1;
}
system->Hcap = MAX( system->numH * SAFER_ZONE, MIN_CAP );
Kurt A. O'Hearn
committed
//Sync_System( system );
//Allocate_System( system, system->local_cap, system->total_cap, msg );
Kurt A. O'Hearn
committed
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 );
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, msg ) == FAILURE )
// {
// 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
}
Compute_Kinetic_Energy( system, data, mpi_data->comm_mesh3D );
return SUCCESS;
}
/************************ 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, char *msg )
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;
Evolve = Velocity_Verlet_NVE;
control->virial = 0;
break;
case bNVT:
data->N_f = 3 * system->bigN + 1;
Evolve = Velocity_Verlet_Berendsen_NVT;
control->virial = 0;
break;
case nhNVT:
Kurt A. O'Hearn
committed
fprintf( stderr, "[WARNING] Nose-Hoover NVT is still under testing.\n" );
data->N_f = 3 * system->bigN + 1;
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;
case sNPT: /* Semi-Isotropic NPT */
data->N_f = 3 * system->bigN + 4;
Evolve = Velocity_Verlet_Berendsen_NPT;
control->virial = 1;
if ( !control->restart )
break;
case iNPT: /* Isotropic NPT */
data->N_f = 3 * system->bigN + 2;
Evolve = Velocity_Verlet_Berendsen_NPT;
control->virial = 1;
if ( !control->restart )
fprintf( stderr, "p%d: init_simulation_data: option not yet implemented\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
data->N_f = 3 * system->bigN + 9;
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:
fprintf( stderr, "p%d: init_simulation_data: ensemble not recognized\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
/* initialize the timer(s) */
MPI_Barrier( MPI_COMM_WORLD ); // wait for everyone to come here
if ( system->my_rank == MASTER_NODE )
{
data->timing.start = Get_Time( );
Kurt A. O'Hearn
committed
Kurt A. O'Hearn
committed
#elif defined(LAMMPS_REAX)
int Init_System( reax_system *system, char *msg )
{
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
fprintf( stderr, "p%d: local_cap=%d total_cap=%d\n",
system->my_rank, system->local_cap, system->total_cap );
Allocate_System( system, system->local_cap, system->total_cap, msg );
return SUCCESS;
Kurt A. O'Hearn
committed
void Init_Simulation_Data( reax_system *system, control_params *control,
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
{
fprintf( stderr, "Warning: non-zero lower Taper-radius cutoff\n" );
Kurt A. O'Hearn
committed
}
if ( swb < 0 )
{
fprintf( stderr, "Negative upper Taper-radius cutoff\n" );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
}
else if ( swb < 5 )
Kurt A. O'Hearn
committed
{
fprintf( stderr, "Warning: very low Taper-radius cutoff: %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, char *msg )
Kurt A. O'Hearn
committed
Allocate_Workspace( system, control, workspace, system->local_cap,
system->total_cap, msg );
memset( &(workspace->realloc), 0, sizeof(reallocate_data) );
Reset_Workspace( system, workspace );
/* Initialize the Taper function */
Init_Taper( control, workspace );
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, char *msg )
Kurt A. O'Hearn
committed
int i, block[11];
MPI_Aint base, disp[11];
MPI_Datatype type[11];
mpi_atom sample;
Kurt A. O'Hearn
committed
restart_atom r_sample;
rvec rvec_sample;
rvec2 rvec2_sample;
/* setup the world */
mpi_data->world = MPI_COMM_WORLD;
Kurt A. O'Hearn
committed
/* mpi_atom: orig_id, imprt_id, type, num_bonds, num_hbonds, name,
* x, v, f_old, s, t */
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;
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
// MPI_Get_address( &sample, &base );
// MPI_Get_address( &(sample.orig_id), disp + 0 );
// MPI_Get_address( &(sample.imprt_id), disp + 1 );
// MPI_Get_address( &(sample.type), disp + 2 );
// MPI_Get_address( &(sample.num_bonds), disp + 3 );
// MPI_Get_address( &(sample.num_hbonds), disp + 4 );
// MPI_Get_address( &(sample.name), disp + 5 );
// MPI_Get_address( &(sample.x[0]), disp + 6 );
// MPI_Get_address( &(sample.v[0]), disp + 7 );
// MPI_Get_address( &(sample.f_old[0]), disp + 8 );
// MPI_Get_address( &(sample.s[0]), disp + 9 );
// MPI_Get_address( &(sample.t[0]), disp + 10 );
// for ( i = 0; i < 11; ++i )
// {
// disp[i] -= base;
// }
disp[0] = offsetof( mpi_atom, orig_id );
disp[1] = offsetof( mpi_atom, imprt_id );
disp[2] = offsetof( mpi_atom, type );
disp[3] = offsetof( mpi_atom, num_bonds );
disp[4] = offsetof( mpi_atom, num_hbonds );
disp[5] = offsetof( mpi_atom, name );
disp[6] = offsetof( mpi_atom, x );
disp[7] = offsetof( mpi_atom, v );
disp[8] = offsetof( mpi_atom, f_old );
disp[9] = offsetof( mpi_atom, s );
disp[10] = offsetof( mpi_atom, t );
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, &(mpi_data->mpi_atom_type) );
MPI_Type_commit( &(mpi_data->mpi_atom_type) );
/* boundary_atom - [orig_id, imprt_id, type, num_bonds, num_hbonds, x] */
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, &base );
// MPI_Get_address( &(b_sample.orig_id), disp + 0 );
// MPI_Get_address( &(b_sample.imprt_id), disp + 1 );
// MPI_Get_address( &(b_sample.type), disp + 2 );
// MPI_Get_address( &(b_sample.num_bonds), disp + 3 );
// MPI_Get_address( &(b_sample.num_hbonds), disp + 4 );
// MPI_Get_address( &(b_sample.x[0]), disp + 5 );
// for ( i = 0; i < 6; ++i )
// {
// disp[i] -= base;
// }
disp[0] = offsetof( boundary_atom, orig_id );
disp[1] = offsetof( boundary_atom, imprt_id );
disp[2] = offsetof( boundary_atom, type );
disp[3] = offsetof( boundary_atom, num_bonds );
disp[4] = offsetof( boundary_atom, num_hbonds );
disp[5] = offsetof( boundary_atom, x );
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, &(mpi_data->boundary_atom_type) );
MPI_Type_commit( &(mpi_data->boundary_atom_type) );
/* mpi_rvec */
block[0] = 3;
// MPI_Get_address( &rvec_sample, &base );
// MPI_Get_address( &(rvec_sample[0]), disp + 0 );
// for ( i = 0; i < 1; ++i )
// {
// disp[i] -= base;
// }
disp[0] = 0;
MPI_Type_create_struct( 1, block, disp, type, &(mpi_data->mpi_rvec) );
MPI_Type_commit( &(mpi_data->mpi_rvec) );
/* mpi_rvec2 */
block[0] = 2;
// MPI_Get_address( &rvec2_sample, &base );
// MPI_Get_address( &(rvec2_sample[0]), disp + 0 );
// for ( i = 0; i < 1; ++i )
// {
// disp[i] -= base;
// }
disp[0] = 0;
MPI_Type_create_struct( 1, block, disp, type, &(mpi_data->mpi_rvec2) );
/* restart_atom - [orig_id, type, name, x, v] */
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, &base );
// MPI_Get_address( &(r_sample.orig_id), disp + 0 );
// MPI_Get_address( &(r_sample.type), disp + 1 );
// MPI_Get_address( &(r_sample.name), disp + 2 );
// MPI_Get_address( &(r_sample.x[0]), disp + 3 );
// MPI_Get_address( &(r_sample.v[0]), disp + 4 );
// for ( i = 0; i < 5; ++i )
// {
// disp[i] -= base;
// }
disp[0] = offsetof( restart_atom, orig_id );
disp[1] = offsetof( restart_atom, type );
disp[2] = offsetof( restart_atom, name );
disp[3] = offsetof( restart_atom, x );
disp[4] = offsetof( restart_atom, v );
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, &(mpi_data->restart_atom_type) );
}
/********************** allocate lists *************************/
Kurt A. O'Hearn
committed
int Init_Lists( reax_system *system, control_params *control,
simulation_data *data, storage *workspace, reax_list **lists,
mpi_datatypes *mpi_data, char *msg )
int i, num_nbrs;
int total_hbonds, total_bonds, bond_cap, num_3body, cap_3body, Htop;
int *hb_top, *bond_top;
//for( i = 0; i < MAX_NBRS; ++i ) nrecv[i] = system->my_nbrs[i].est_recv;
//system->N = SendRecv( system, mpi_data, mpi_data->boundary_atom_type, nrecv,
// Sort_Boundary_Atoms, Unpack_Exchange_Message, TRUE );
Make_List( system->total_cap, num_nbrs, TYP_FAR_NEIGHBOR, lists[FAR_NBRS] );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: allocated far_nbrs: num_far=%d, space=%dMB\n",
Kurt A. O'Hearn
committed
system->my_rank, num_nbrs,
(int)(num_nbrs * sizeof(far_neighbor_data) / (1024 * 1024)) );
Generate_Neighbor_Lists( system, data, workspace, lists );
Kurt A. O'Hearn
committed
bond_top = (int*) scalloc( system->total_cap, sizeof(int), "Init_Lists::bond_top" );
hb_top = (int*) scalloc( system->local_cap, sizeof(int), "Init_Lists::hb_top" );
// hb_top = (int*) scalloc( system->Hcap, sizeof(int), "Init_Lists::hb_top" );
Kurt A. O'Hearn
committed
&Htop, hb_top, bond_top, &num_3body );
// Host_Estimate_Sparse_Matrix( system, control, lists, system->local_cap, system->total_cap,
// &Htop, hb_top, bond_top, &num_3body );
Kurt A. O'Hearn
committed
Allocate_Matrix( &(workspace->H), system->n, Htop );
//MATRIX CHANGES
//workspace->L = NULL;
//workspace->U = NULL;
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: allocated H matrix: Htop=%d, space=%dMB\n",
system->my_rank, Htop,
(int)(Htop * sizeof(sparse_matrix_entry) / (1024 * 1024)) );
Kurt A. O'Hearn
committed
if ( control->hbond_cut > 0.0 )
// init H indexes
total_hbonds = 0;
for ( i = 0; i < system->n; ++i )
{
system->my_atoms[i].num_hbonds = hb_top[i];
total_hbonds += hb_top[i];
}
total_hbonds = MAX( total_hbonds * SAFER_ZONE, MIN_CAP * MIN_HBONDS );
Kurt A. O'Hearn
committed
// DANIEL, to make Mpi_Not_Gpu_Validate_Lists() not complain that max_hbonds is 0
system->max_hbonds = total_hbonds * SAFER_ZONE;
Make_List( system->Hcap, total_hbonds, TYP_HBOND, lists[HBONDS] );
fprintf( stderr, "p%d: allocated hbonds: total_hbonds=%d, space=%dMB\n",
system->my_rank, total_hbonds,
(int)(total_hbonds * sizeof(hbond_data) / (1024 * 1024)) );
}
/* bonds list */
total_bonds = 0;
for ( i = 0; i < system->N; ++i )
{
system->my_atoms[i].num_bonds = bond_top[i];
total_bonds += bond_top[i];
Kurt A. O'Hearn
committed
// DANIEL, to make Mpi_Not_Gpu_Validate_Lists() not complain that max_bonds is 0
system->max_bonds[i] = MAX( bond_top[i], MIN_BONDS );
}
bond_cap = MAX( total_bonds * SAFE_ZONE, MIN_CAP * MIN_BONDS );
Make_List( system->total_cap, bond_cap, TYP_BOND, lists[BONDS] );
fprintf( stderr, "p%d: allocated bonds: total_bonds=%d, space=%dMB\n",
system->my_rank, bond_cap,
(int)(bond_cap * sizeof(bond_data) / (1024 * 1024)) );
/* 3bodies list */
cap_3body = MAX( num_3body * SAFE_ZONE, MIN_3BODIES );
Make_List( bond_cap, cap_3body, TYP_THREE_BODY, lists[THREE_BODIES] );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: allocated 3-body list: num_3body=%d, space=%dMB\n",
system->my_rank, cap_3body,
(int)(cap_3body * sizeof(three_body_interaction_data) / (1024 * 1024)) );
Make_List( system->total_cap, bond_cap * 8, TYP_DDELTA, lists[DDELTAS] );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: allocated dDelta list: num_ddelta=%d space=%ldMB\n",
system->my_rank, bond_cap * 30,
bond_cap * 8 * sizeof(dDelta_data) / (1024 * 1024) );
Make_List( bond_cap, bond_cap * 50, TYP_DBO, lists[DBOS] );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: allocated dbond list: num_dbonds=%d space=%ldMB\n",
system->my_rank, bond_cap * MAX_BONDS * 3,
bond_cap * MAX_BONDS * 3 * sizeof(dbond_data) / (1024 * 1024) );
Kurt A. O'Hearn
committed
sfree( hb_top, "Init_Lists::hb_top" );
sfree( bond_top, "Init_Lists::bond_top" );
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
host_scratch = (void *) smalloc( HOST_SCRATCH_SIZE, "Initialize::host_scratch" );
Kurt A. O'Hearn
committed
Init_MPI_Datatypes( system, workspace, mpi_data, msg );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized mpi datatypes\n", system->my_rank );
if ( Init_System(system, control, data, workspace, mpi_data, msg) == FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: system could not be initialized! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: system initialized\n", system->my_rank );
Kurt A. O'Hearn
committed
Init_Simulation_Data( system, control, data, msg );
fprintf( stderr, "p%d: initialized simulation data\n", system->my_rank );
Kurt A. O'Hearn
committed
Init_Workspace( system, control, workspace, msg );
fprintf( stderr, "p%d: initialized workspace\n", system->my_rank );
if ( Init_Lists( system, control, data, workspace, lists, mpi_data, msg ) ==
FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: system could not be initialized! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized lists\n", system->my_rank );
if ( Init_Output_Files(system, control, out_control, mpi_data, msg) == FAILURE)
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: could not open output files! terminating...\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: output files opened\n", system->my_rank );
Kurt A. O'Hearn
committed
if ( Init_Lookup_Tables(system, control, workspace->Tap, mpi_data, msg) == FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: couldn't create lookup table! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized lookup tables\n", system->my_rank );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized force functions\n", system->my_rank );
Kurt A. O'Hearn
committed
/*#ifdef TEST_FORCES
Init_Force_Test_Functions();
fprintf(stderr,"p%d: initialized force test functions\n",system->my_rank);
#endif */
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, msg );
fprintf( stderr, "p%d: initialized simulation data\n", system->my_rank );
fprintf( stderr, "p%d: pure initialized simulation data\n", system->my_rank );
Kurt A. O'Hearn
committed
Init_Workspace( system, control, workspace, msg );
fprintf( stderr, "p%d: initialized workspace\n", system->my_rank );
fprintf( stderr, "p%d: pure initialized workspace\n", system->my_rank );
if ( Init_Lists( system, control, data, workspace, lists, mpi_data, msg ) ==
FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: system could not be initialized! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
fprintf( stderr, "p%d: initialized lists\n", system->my_rank );
fprintf( stderr, "p%d: pure initialized lists done \n", system->my_rank );
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
Kurt A. O'Hearn
committed
host_scratch = (void *) smalloc( HOST_SCRATCH_SIZE, "Initialize::host_scratch" );
Kurt A. O'Hearn
committed
if ( Init_System(system, msg) == FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: system could not be initialized! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: system initialized\n", system->my_rank );
Kurt A. O'Hearn
committed
Init_Simulation_Data( system, control, data, msg );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized simulation data\n", system->my_rank );
Kurt A. O'Hearn
committed
Init_Workspace( system, control, workspace, msg );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized workspace\n", system->my_rank );
Kurt A. O'Hearn
committed
Init_MPI_Datatypes( system, workspace, mpi_data, msg );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized mpi datatypes\n", system->my_rank );
if ( Init_Lists( system, control, workspace, lists, msg ) == FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: system could not be initialized! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized lists\n", system->my_rank );
if ( Init_Output_Files(system, control, out_control, mpi_data, msg) == FAILURE)
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: could not open output files! terminating...\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: output files opened\n", system->my_rank );
Kurt A. O'Hearn
committed
if ( Init_Lookup_Tables( system, control, workspace->Tap, mpi_data, msg ) == FAILURE )
{
fprintf( stderr, "p%d: %s\n", system->my_rank, msg );
fprintf( stderr, "p%d: couldn't create lookup table! terminating.\n",
system->my_rank );
MPI_Abort( MPI_COMM_WORLD, CANNOT_INITIALIZE );
}
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized lookup tables\n", system->my_rank );
Kurt A. O'Hearn
committed
fprintf( stderr, "p%d: initialized force functions\n", system->my_rank );
Kurt A. O'Hearn
committed
/*#if defined(TEST_FORCES)
Init_Force_Test_Functions();
fprintf(stderr,"p%d: initialized force test functions\n",system->my_rank);
#endif*/
}