Skip to content
Snippets Groups Projects
Commit 58a9a0e0 authored by Kurt A. O'Hearn's avatar Kurt A. O'Hearn
Browse files

PG-PuReMD: port sPuReMD logic for seperating into library and driver.

parent dd042bce
No related branches found
No related tags found
No related merge requests found
......@@ -16,15 +16,15 @@ NVCCFLAGS += --compiler-options "$(DEFS) $(NVCC_OPT_CODE_DEFS) -O3 -funroll-loop
endif
bin_PROGRAMS = bin/pg-puremd
bin_pg_puremd_SOURCES = src/allocate.c src/basic_comm.c src/ffield.c src/grid.c src/list.c \
lib_LTLIBRARIES = lib/libpuremd.la
lib_libpuremd_la_SOURCES = src/allocate.c src/basic_comm.c src/ffield.c src/grid.c src/list.c \
src/lookup.c src/io_tools.c src/reset_tools.c src/restart.c src/random.c \
src/tool_box.c src/traj.c src/analyze.c src/box.c src/system_props.c \
src/control.c src/comm_tools.c src/geo_tools.c src/lin_alg.c src/neighbors.c \
src/charges.c src/bond_orders.c src/multi_body.c src/bonds.c src/valence_angles.c \
src/hydrogen_bonds.c src/torsion_angles.c src/nonbonded.c src/forces.c \
src/integrate.c src/init_md.c src/parallelreax.c
include_HEADERS = src/reax_types.h src/index_utils.h \
src/integrate.c src/init_md.c src/puremd.c
lib_libpuremd_la_SOURCES += src/reax_types.h src/index_utils.h \
src/allocate.h src/basic_comm.h src/ffield.h src/grid.h src/list.h \
src/lookup.h src/io_tools.h src/reset_tools.h src/restart.h src/random.h src/vector.h \
src/tool_box.h src/traj.h src/analyze.h src/box.h src/system_props.h \
......@@ -32,9 +32,10 @@ include_HEADERS = src/reax_types.h src/index_utils.h \
src/charges.h src/bond_orders.h src/multi_body.h src/bonds.h src/valence_angles.h \
src/hydrogen_bonds.h src/torsion_angles.h src/nonbonded.h src/forces.h \
src/integrate.h src/init_md.h
include_HEADERS = src/puremd.h
if USE_CUDA
bin_pg_puremd_SOURCES += src/cuda/cuda_utils.cu src/cuda/cuda_allocate.cu src/cuda/cuda_environment.cu \
lib_libpuremd_la_SOURCES += src/cuda/cuda_utils.cu src/cuda/cuda_allocate.cu src/cuda/cuda_environment.cu \
src/cuda/cuda_system_props.cu src/cuda/cuda_reduction.cu src/cuda/cuda_box.cu src/cuda/cuda_list.cu \
src/cuda/cuda_copy.cu src/cuda/cuda_reset_tools.cu src/cuda/cuda_random.cu \
src/cuda/cuda_neighbors.cu src/cuda/cuda_bond_orders.cu src/cuda/cuda_bonds.cu \
......@@ -43,7 +44,7 @@ bin_pg_puremd_SOURCES += src/cuda/cuda_utils.cu src/cuda/cuda_allocate.cu src/cu
src/cuda/cuda_charges.cu src/cuda/cuda_lin_alg.cu \
src/cuda/cuda_nonbonded.cu src/cuda/cuda_integrate.cu src/cuda/cuda_post_evolve.cu \
src/cuda/cuda_init_md.cu src/cuda/cuda_lookup.cu
include_HEADERS += src/cuda/cuda_helpers.h src/cuda/cuda_shuffle.h \
lib_libpuremd_la_SOURCES += src/cuda/cuda_helpers.h src/cuda/cuda_shuffle.h \
src/cuda/cuda_utils.h src/cuda/cuda_allocate.h src/cuda/cuda_environment.h \
src/cuda/cuda_system_props.h src/cuda/cuda_reduction.h src/cuda/cuda_box.h src/cuda/cuda_list.h \
src/cuda/cuda_copy.h src/cuda/cuda_reset_tools.h src/cuda/cuda_random.h src/cuda/cuda_vector.h \
......@@ -55,20 +56,28 @@ include_HEADERS += src/cuda/cuda_helpers.h src/cuda/cuda_shuffle.h \
src/cuda/cuda_init_md.h src/cuda/cuda_lookup.h
if DEBUG
bin_pg_puremd_SOURCES += src/cuda/cuda_validation.cu
include_HEADERS += src/cuda/cuda_validation.h
lib_libpuremd_la_SOURCES += src/cuda/cuda_validation.cu
lib_libpuremd_la_SOURCES += src/cuda/cuda_validation.h
endif
# dummy source to cause C linking
nodist_EXTRA_bin_pg_puremd_SOURCES = src/dummy.c
nodist_EXTRA_lib_libpuremd_la_SOURCES = src/dummy.c
endif
lib_libpuremd_la_CFLAGS = $(AM_CFLAGS) $(MPI_CFLAGS) $(CFLAGS)
lib_libpuremd_la_CPPFLAGS = $(AM_CPPFLAGS) $(CPPFLAGS)
lib_libpuremd_la_LIBADD = $(AM_LDADD) $(MPI_LIBS) $(LDADD) -lstdc++
lib_libpuremd_la_LDFLAGS = -version-info 1:0:0
if USE_CUDA
lib_libpuremd_la_CFLAGS += $(CUDA_CFLAGS)
lib_libpuremd_la_LIBADD += $(CUDA_LIBS)
endif
bin_PROGRAMS = bin/pg-puremd
bin_pg_puremd_SOURCES = src/driver.c
bin_pg_puremd_CFLAGS = $(AM_CFLAGS) $(MPI_CFLAGS) $(CFLAGS)
bin_pg_puremd_CPPFLAGS = $(AM_CPPFLAGS) $(CPPFLAGS)
bin_pg_puremd_LDADD = $(AM_LDADD) $(MPI_LIBS) $(LDADD) -lstdc++
bin_pg_puremd_CPPFLAGS = -I src $(AM_CPPFLAGS) $(CPPFLAGS)
bin_pg_puremd_LDADD = lib/libpuremd.la $(AM_LDADD) $(MPI_LIBS) $(LDADD) -lstdc++
if USE_CUDA
bin_pg_puremd_CFLAGS += $(CUDA_CFLAGS)
bin_pg_puremd_LDADD += $(CUDA_LIBS)
......
......@@ -56,6 +56,66 @@ m4_ifndef([AC_AUTOCONF_VERSION],
[m4_copy([m4_PACKAGE_VERSION], [AC_AUTOCONF_VERSION])])dnl
_AM_AUTOCONF_VERSION(m4_defn([AC_AUTOCONF_VERSION]))])
# Copyright (C) 2011-2017 Free Software Foundation, Inc.
#
# This file is free software; the Free Software Foundation
# gives unlimited permission to copy and/or distribute it,
# with or without modifications, as long as this notice is preserved.
# AM_PROG_AR([ACT-IF-FAIL])
# -------------------------
# Try to determine the archiver interface, and trigger the ar-lib wrapper
# if it is needed. If the detection of archiver interface fails, run
# ACT-IF-FAIL (default is to abort configure with a proper error message).
AC_DEFUN([AM_PROG_AR],
[AC_BEFORE([$0], [LT_INIT])dnl
AC_BEFORE([$0], [AC_PROG_LIBTOOL])dnl
AC_REQUIRE([AM_AUX_DIR_EXPAND])dnl
AC_REQUIRE_AUX_FILE([ar-lib])dnl
AC_CHECK_TOOLS([AR], [ar lib "link -lib"], [false])
: ${AR=ar}
AC_CACHE_CHECK([the archiver ($AR) interface], [am_cv_ar_interface],
[AC_LANG_PUSH([C])
am_cv_ar_interface=ar
AC_COMPILE_IFELSE([AC_LANG_SOURCE([[int some_variable = 0;]])],
[am_ar_try='$AR cru libconftest.a conftest.$ac_objext >&AS_MESSAGE_LOG_FD'
AC_TRY_EVAL([am_ar_try])
if test "$ac_status" -eq 0; then
am_cv_ar_interface=ar
else
am_ar_try='$AR -NOLOGO -OUT:conftest.lib conftest.$ac_objext >&AS_MESSAGE_LOG_FD'
AC_TRY_EVAL([am_ar_try])
if test "$ac_status" -eq 0; then
am_cv_ar_interface=lib
else
am_cv_ar_interface=unknown
fi
fi
rm -f conftest.lib libconftest.a
])
AC_LANG_POP([C])])
case $am_cv_ar_interface in
ar)
;;
lib)
# Microsoft lib, so override with the ar-lib wrapper script.
# FIXME: It is wrong to rewrite AR.
# But if we don't then we get into trouble of one sort or another.
# A longer-term fix would be to have automake use am__AR in this case,
# and then we could set am__AR="$am_aux_dir/ar-lib \$(AR)" or something
# similar.
AR="$am_aux_dir/ar-lib $AR"
;;
unknown)
m4_default([$1],
[AC_MSG_ERROR([could not determine $AR interface])])
;;
esac
AC_SUBST([AR])dnl
])
# AM_AUX_DIR_EXPAND -*- Autoconf -*-
# Copyright (C) 2001-2017 Free Software Foundation, Inc.
......@@ -1153,3 +1213,8 @@ AC_SUBST([am__untar])
m4_include([../m4/acx_mpi.m4])
m4_include([../m4/ax_compiler_vendor.m4])
m4_include([../m4/ax_cuda.m4])
m4_include([../m4/libtool.m4])
m4_include([../m4/ltoptions.m4])
m4_include([../m4/ltsugar.m4])
m4_include([../m4/ltversion.m4])
m4_include([../m4/lt~obsolete.m4])
......@@ -10,8 +10,9 @@ sav_CFLAGS="$CFLAGS"
AM_INIT_AUTOMAKE([1.15 subdir-objects -Wall -Werror foreign])
# Enable silent build rules by default.
m4_ifdef([AM_SILENT_RULES], [AM_SILENT_RULES([yes])], [AC_SUBST([AM_DEFAULT_VERBOSITY],[1])])
#LT_PREREQ([2.2])
#LT_INIT([dlopen])
AM_PROG_AR
LT_PREREQ([2.2])
LT_INIT([dlopen])
AC_CONFIG_MACRO_DIR([../m4])
......
/*----------------------------------------------------------------------
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
the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
See the GNU General Public License for more details:
<http://www.gnu.org/licenses/>.
----------------------------------------------------------------------*/
#include <stdio.h>
#include <stdlib.h>
#include "puremd.h"
#define INVALID_INPUT (-1)
static void usage( char * argv[] )
{
fprintf( stderr, "usage: ./%s geometry_file force_field_params_file control_file\n", argv[0] );
}
int main( int argc, char* argv[] )
{
void *handle;
int ret;
MPI_Init( &argc, &argv );
if ( argc != 4 )
{
usage( argv );
MPI_Abort( MPI_COMM_WORLD, INVALID_INPUT );
}
handle = setup( argv[1], argv[2], argv[3] );
ret = PUREMD_FAILURE;
if ( handle != NULL )
{
ret = simulate( handle );
}
MPI_Finalized( &ret );
if ( !ret )
{
MPI_Finalize( );
}
if ( ret == PUREMD_SUCCESS )
{
ret = cleanup( handle );
}
return (ret == PUREMD_SUCCESS) ? 0 : (-1);
}
......@@ -30,7 +30,7 @@
#endif
void Read_Force_Field( char *ffield_file, reax_interaction *reax,
void Read_Force_Field_File( char *ffield_file, reax_interaction *reax,
reax_system *system, control_params *control )
{
FILE *fp;
......
......@@ -29,7 +29,7 @@
extern "C" {
#endif
void Read_Force_Field( char*, reax_interaction*, reax_system *, control_params* );
void Read_Force_Field_File( char*, reax_interaction*, reax_system *, control_params* );
#ifdef __cplusplus
}
......
......@@ -67,7 +67,7 @@ void Count_Geo_Atoms( FILE *geo, reax_system *system )
}
void Read_Geo( char* geo_file, reax_system* system, control_params *control,
void Read_Geo_File( char* geo_file, reax_system* system, control_params *control,
simulation_data *data, storage *workspace, mpi_datatypes *mpi_data )
{
int i, j, serial, top;
......@@ -79,7 +79,7 @@ void Read_Geo( char* geo_file, reax_system* system, control_params *control,
reax_atom *atom;
/* open the geometry file */
geo = sfopen( geo_file, "r", "Read_Geo::geo" );
geo = sfopen( geo_file, "r", "Read_Geo_File::geo" );
/* read box information */
fscanf( geo, CUSTOM_BOXGEO_FORMAT,
......@@ -136,7 +136,7 @@ void Read_Geo( char* geo_file, reax_system* system, control_params *control,
}
}
sfclose( geo, "Read_Geo::geo" );
sfclose( geo, "Read_Geo_File::geo" );
#if defined(DEBUG_FOCUS)
fprintf( stderr, "p%d: finished reading the geo file\n", system->my_rank );
......@@ -247,7 +247,7 @@ void Count_PDB_Atoms( FILE *geo, reax_system *system )
}
void Read_PDB( char* pdb_file, reax_system* system, control_params *control,
void Read_PDB_File( char* pdb_file, reax_system* system, control_params *control,
simulation_data *data, storage *workspace, mpi_datatypes *mpi_data )
{
FILE *pdb;
......@@ -264,7 +264,7 @@ void Read_PDB( char* pdb_file, reax_system* system, control_params *control,
rvec x;
reax_atom *atom;
pdb = sfopen( pdb_file, "r", "Read_PDB::pdb" );
pdb = sfopen( pdb_file, "r", "Read_PDB_File::pdb" );
/* allocate memory for tokenizing pdb lines */
Allocate_Tokenizer_Space( &s, &s1, &tmp );
......@@ -451,14 +451,14 @@ void Read_PDB( char* pdb_file, reax_system* system, control_params *control,
}
}
sfclose( pdb, "Read_PDB::pdb" );
sfclose( pdb, "Read_PDB_File::pdb" );
}
/* PDB serials are written without regard to the order, we'll see if this
* cause trouble, if so we'll have to rethink this approach
* Also, we do not write connect lines yet. */
void Write_PDB( reax_system* system, reax_list* bonds, simulation_data *data,
void Write_PDB_File( reax_system* system, reax_list* bonds, simulation_data *data,
control_params *control, mpi_datatypes *mpi_data, output_controls *out_control )
{
int i, cnt, me, np, buffer_req, buffer_len;
......@@ -477,7 +477,7 @@ void Write_PDB( reax_system* system, reax_list* bonds, simulation_data *data,
np = control->nprocs;
/* Allocation*/
line = smalloc( sizeof(char) * PDB_ATOM_FORMAT_O_LENGTH, "Write_PDB::line" );
line = smalloc( sizeof(char) * PDB_ATOM_FORMAT_O_LENGTH, "Write_PDB_File::line" );
if ( me == MASTER_NODE )
{
buffer_req = system->bigN * PDB_ATOM_FORMAT_O_LENGTH;
......@@ -487,7 +487,7 @@ void Write_PDB( reax_system* system, reax_list* bonds, simulation_data *data,
buffer_req = system->n * PDB_ATOM_FORMAT_O_LENGTH;
}
buffer = smalloc( sizeof(char) * buffer_req, "Write_PDB::buffer" );
buffer = smalloc( sizeof(char) * buffer_req, "Write_PDB_File::buffer" );
pdb = NULL;
line[0] = 0;
......@@ -512,7 +512,7 @@ void Write_PDB( reax_system* system, reax_list* bonds, simulation_data *data,
sprintf( fname, "%s-%d.pdb", control->sim_name, data->step );
pdb = sfopen( fname, "w", "Write_PDB::pdb" );
pdb = sfopen( fname, "w", "Write_PDB_File::pdb" );
fprintf( pdb, PDB_CRYST1_FORMAT_O,
"CRYST1",
system->big_box.box_norms[0], system->big_box.box_norms[1],
......@@ -563,7 +563,7 @@ void Write_PDB( reax_system* system, reax_list* bonds, simulation_data *data,
if ( me == MASTER_NODE)
{
fprintf( pdb, "%s", buffer );
sfclose( pdb, "Write_PDB::pdb" );
sfclose( pdb, "Write_PDB_File::pdb" );
}
/* Writing connect information */
......@@ -585,6 +585,6 @@ void Write_PDB( reax_system* system, reax_list* bonds, simulation_data *data,
}
*/
sfree( buffer, "Write_PDB::buffer" );
sfree( line, "Write_PDB::line" );
sfree( buffer, "Write_PDB_File::buffer" );
sfree( line, "Write_PDB_File::line" );
}
......@@ -115,13 +115,13 @@ COLUMNS DATA TYPE FIELD DEFINITION
extern "C" {
#endif
void Read_Geo( char*, reax_system*, control_params*,
void Read_Geo_File( char*, reax_system*, control_params*,
simulation_data*, storage*, mpi_datatypes* );
void Read_PDB( char*, reax_system*, control_params*,
void Read_PDB_File( char*, reax_system*, control_params*,
simulation_data*, storage*, mpi_datatypes* );
void Write_PDB( reax_system*, reax_list*, simulation_data*,
void Write_PDB_File( reax_system*, reax_list*, simulation_data*,
control_params*, mpi_datatypes*, output_controls* );
#ifdef __cplusplus
......
This diff is collapsed.
/*----------------------------------------------------------------------
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
the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
See the GNU General Public License for more details:
<http://www.gnu.org/licenses/>.
----------------------------------------------------------------------*/
#ifndef __PUREMD_H_
#define __PUREMD_H_
#include "reax_types.h"
#define PUREMD_SUCCESS (0)
#define PUREMD_FAILURE (-1)
#ifdef __cplusplus
extern "C" {
#endif
void* setup( const char * const, const char * const,
const char * const );
int setup_callback( const void * const, const callback_function );
int simulate( const void * const );
int cleanup( const void * const );
reax_atom* get_atoms( const void * const );
int set_output_enabled( const void * const, const int );
#ifdef __cplusplus
}
#endif
#endif
......@@ -548,6 +548,7 @@ typedef struct molecule molecule;
typedef struct LR_data LR_data;
typedef struct cubic_spline_coef cubic_spline_coef;
typedef struct LR_lookup_table LR_lookup_table;
typedef struct puremd_handle puremd_handle;
/* function pointer definitions */
......@@ -563,6 +564,8 @@ typedef real (*lookup_function)( real );
typedef void (*message_sorter)( reax_system*, int, int, int, mpi_out_data* );
/**/
typedef void (*unpacker)( reax_system*, int, void*, int, neighbor_proc*, int );
/**/
typedef void (*callback_function)(reax_atom*, simulation_data*, reax_list*);
/* struct definitions */
......@@ -2458,6 +2461,30 @@ struct LR_lookup_table
};
/* Handle for working with an instance of the PuReMD library */
struct puremd_handle
{
/* System info. struct pointer */
reax_system *system;
/* System struct pointer */
control_params *control;
/* Control parameters struct pointer */
simulation_data *data;
/* Internal workspace struct pointer */
storage *workspace;
/* Reax interaction list struct pointer */
reax_list **lists;
/* Output controls struct pointer */
output_controls *out_control;
/* MPI datatypes struct pointer */
mpi_datatypes *mpi_data;
/* TRUE if file I/O for simulation output enabled, FALSE otherwise */
int output_enabled;
/* Callback for getting simulation state at the end of each time step */
callback_function callback;
};
/* CUDA-specific globals */
extern reax_list **dev_lists;
extern storage *dev_workspace;
......
......@@ -29,7 +29,7 @@
#include "vector.h"
void Write_Binary_Restart( reax_system *system, control_params *control,
void Write_Binary_Restart_File( reax_system *system, control_params *control,
simulation_data *data, output_controls *out_control,
mpi_datatypes *mpi_data )
{
......@@ -49,7 +49,7 @@ void Write_Binary_Restart( reax_system *system, control_params *control,
{
/* master handles the restart file */
sprintf( fname, "%s.res%d", control->sim_name, data->step );
fres = sfopen( fname, "wb", "Write_Binary_Restart::fres" );
fres = sfopen( fname, "wb", "Write_Binary_Restart_File::fres" );
/* master can write the header by itself */
res_header.step = data->step;
......@@ -64,12 +64,12 @@ void Write_Binary_Restart( reax_system *system, control_params *control,
/* master needs to allocate space for all atoms */
buffer = scalloc( system->bigN, sizeof(restart_atom),
"Write_Binary_Restart::buffer" );
"Write_Binary_Restart_File::buffer" );
}
else
{
buffer = scalloc( system->n, sizeof(restart_atom),
"Write_Binary_Restart::buffer" );
"Write_Binary_Restart_File::buffer" );
}
/* fill in the buffers */
......@@ -108,14 +108,14 @@ void Write_Binary_Restart( reax_system *system, control_params *control,
if ( me == MASTER_NODE )
{
fwrite( buffer, system->bigN, sizeof(restart_atom), fres );
sfclose( fres, "Write_Binary_Restart::fres" );
sfclose( fres, "Write_Binary_Restart_File::fres" );
}
sfree( buffer, "Write_Binary_Restart::buffer" );
sfree( buffer, "Write_Binary_Restart_File::buffer" );
}
void Write_Restart( reax_system *system, control_params *control,
void Write_Restart_File( reax_system *system, control_params *control,
simulation_data *data, output_controls *out_control,
mpi_datatypes *mpi_data )
{
......@@ -135,7 +135,7 @@ void Write_Restart( reax_system *system, control_params *control,
if ( me == MASTER_NODE )
{
sprintf( fname, "%s.res%d", control->sim_name, data->step );
fres = sfopen( fname, "w", "Write_Restart::fres" );
fres = sfopen( fname, "w", "Write_Restart_File::fres" );
/* write the header - only master writes it */
fprintf( fres, RESTART_HEADER,
......@@ -156,7 +156,7 @@ void Write_Restart( reax_system *system, control_params *control,
buffer_req = system->n * RESTART_LINE_LEN + 1;
}
buffer = smalloc( sizeof(char) * buffer_req, "Write_Restart::buffer" );
buffer = smalloc( sizeof(char) * buffer_req, "Write_Restart_File::buffer" );
line[0] = 0;
buffer[0] = 0;
......@@ -200,10 +200,10 @@ void Write_Restart( reax_system *system, control_params *control,
if ( me == MASTER_NODE )
{
fprintf( fres, "%s", buffer );
sfclose( fres, "Write_Restart::fres" );
sfclose( fres, "Write_Restart_File::fres" );
}
sfree( buffer, "Write_Restart::buffer" );
sfree( line, "Write_Restart::line" );
sfree( buffer, "Write_Restart_File::buffer" );
sfree( line, "Write_Restart_File::line" );
}
......@@ -235,7 +235,7 @@ void Count_Binary_Restart_Atoms( FILE *fres, reax_system *system )
}
void Read_Binary_Restart( char *res_file, reax_system *system,
void Read_Binary_Restart_File( char *res_file, reax_system *system,
control_params *control, simulation_data *data,
storage *workspace, mpi_datatypes *mpi_data )
{
......@@ -245,7 +245,7 @@ void Read_Binary_Restart( char *res_file, reax_system *system,
restart_atom res_atom;
reax_atom *p_atom;
fres = sfopen( res_file, "rb", "Read_Binary_Restart::fres" );
fres = sfopen( res_file, "rb", "Read_Binary_Restart_File::fres" );
/* first read the header lines */
fread( &res_header, sizeof(restart_header), 1, fres );
......@@ -305,7 +305,7 @@ void Read_Binary_Restart( char *res_file, reax_system *system,
}
}
sfclose( fres, "Read_Binary_Restart::fres" );
sfclose( fres, "Read_Binary_Restart_File::fres" );
data->step = data->prev_steps;
// nsteps is updated based on the number of steps in the previous run
......@@ -350,7 +350,7 @@ void Count_Restart_Atoms( FILE *fres, reax_system *system )
}
void Read_Restart( char *res_file, reax_system *system,
void Read_Restart_File( char *res_file, reax_system *system,
control_params *control, simulation_data *data,
storage *workspace, mpi_datatypes *mpi_data )
{
......@@ -362,13 +362,13 @@ void Read_Restart( char *res_file, reax_system *system,
rvec x_temp, v_temp;
rtensor box;
fres = sfopen( res_file, "r", "Read_Restart::fres" );
fres = sfopen( res_file, "r", "Read_Restart_File::fres" );
s = smalloc( sizeof(char) * MAX_LINE, "Read_Restart::s" );
tmp = smalloc( sizeof(char*) * MAX_TOKENS, "Read_Restart::tmp" );
s = smalloc( sizeof(char) * MAX_LINE, "Read_Restart_File::s" );
tmp = smalloc( sizeof(char*) * MAX_TOKENS, "Read_Restart_File::tmp" );
for (i = 0; i < MAX_TOKENS; i++)
{
tmp[i] = smalloc( sizeof(char) * MAX_LINE, "Read_Restart::tmp[i]" );
tmp[i] = smalloc( sizeof(char) * MAX_LINE, "Read_Restart_File::tmp[i]" );
}
//read first header lines
......@@ -486,15 +486,15 @@ void Read_Restart( char *res_file, reax_system *system,
top++;
}
}
sfclose( fres, "Read_Restart::fres" );
sfclose( fres, "Read_Restart_File::fres" );
/* free memory allocations at the top */
for ( i = 0; i < MAX_TOKENS; i++ )
{
sfree( tmp[i], "Read_Restart::tmp[i]" );
sfree( tmp[i], "Read_Restart_File::tmp[i]" );
}
sfree( tmp, "Read_Restart::tmp" );
sfree( s, "Read_Restart::s" );
sfree( tmp, "Read_Restart_File::tmp" );
sfree( s, "Read_Restart_File::s" );
data->step = data->prev_steps;
// nsteps is updated based on the number of steps in the previous run
......
......@@ -45,16 +45,16 @@
extern "C" {
#endif
void Write_Binary_Restart( reax_system*, control_params*,
void Write_Binary_Restart_File( reax_system*, control_params*,
simulation_data*, output_controls*, mpi_datatypes* );
void Write_Restart( reax_system*, control_params*,
void Write_Restart_File( reax_system*, control_params*,
simulation_data*, output_controls*, mpi_datatypes* );
void Read_Binary_Restart( char*, reax_system*, control_params*,
void Read_Binary_Restart_File( char*, reax_system*, control_params*,
simulation_data*, storage*, mpi_datatypes* );
void Read_Restart( char*, reax_system*, control_params*,
void Read_Restart_File( char*, reax_system*, control_params*,
simulation_data*, storage*, mpi_datatypes* );
#ifdef __cplusplus
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment