Newer
Older
#include "Angle.h"
#include "Dihedral.h"
#include "Restraint.h"
cmaffeo2
committed
#include <cmath>
cmaffeo2
committed
#include <stdlib.h> /* srand, rand */
#include <time.h> /* time */
#include <iostream>
using namespace std;
#define gpuErrchk(ans) { gpuAssert((ans), __FILE__, __LINE__); }
cmaffeo2
committed
inline void gpuAssert(cudaError_t code, const char *file, int line, bool abort=true) {
if (code != cudaSuccess) {
fprintf(stderr,"CUDA Error: %s %s %d\n", cudaGetErrorString(code), __FILE__, line);
if (abort) exit(code);
}
}
namespace
{
template<class T>
void convertString(const String& token, void* data)
{
exit(1);
}
template<>
void convertString<float>(const String& token, void* data)
{
float* tmp = (float*)data;
*tmp = atof(token);
}
template<>
void convertString<String>(const String& token, void* data)
{
String* tmp = (String*)data;
*tmp = token;
}
template<class T>
void stringToArray(String* str, int& size, T** array)
{
register int num;
String *token;
num = str->tokenCount();
size = num;
*array = new T[num];
token = new String[num];
str->tokenize(token);
for(int i = 0; i < num; ++i)
convertString<T>(token[i], (*array)+i);
delete [] token;
}
}
Configuration::Configuration(const char* config_file, int simNum, bool debug) :
simNum(simNum) {
// Read the parameters.
eguzman6
committed
//type_d = NULL;
eguzman6
committed
//bonds_d = NULL;
//bondMap_d = NULL;
//excludes_d = NULL;
//excludeMap_d = NULL;
//angles_d = NULL;
//dihedrals_d = NULL;
Maxim Belkin
committed
// Get the number of particles
// printf("\nCounting particles specified in the ");
if (restartCoordinates.length() > 0) {
Maxim Belkin
committed
// Read them from the restart file.
num = countRestart(restartCoordinates.val());
if (copyReplicaCoordinates <= 0) {
num /= simNum;
}
Maxim Belkin
committed
} else {
if (readPartsFromFile) readAtoms();
if (numPartsFromFile > 0) {
// Determine number of particles from input file (PDB-style)
Maxim Belkin
committed
num = numPartsFromFile;
} else {
// Sum up all particles in config file
Maxim Belkin
committed
//int num0 = 0;
num = 0;
for (int i = 0; i < numParts; i++) num += part[i].num;
//num = num0;
}
} // end result: variable "num" is set
Maxim Belkin
committed
// Set the number capacity
if (numCap <= 0) numCap = numCapFactor*num; // max number of particles
if (numCap <= 0) numCap = 20;
Maxim Belkin
committed
// Allocate particle variables.
pos = new Vector3[num * simNum];
cmaffeo2
committed
//Han-Yi Chou
if (ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
momentum = new Vector3[num * simNum];
type = new int[num * simNum];
Maxim Belkin
committed
serial = new int[num * simNum];
posLast = new Vector3[num * simNum];
cmaffeo2
committed
//Han-Yi Chou
if(ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
momLast = new Vector3[num * simNum];
Maxim Belkin
committed
name = new String[num * simNum];
currSerial = 0;
// Now, load the coordinates
loadedCoordinates = false;
cmaffeo2
committed
loadedMomentum = false; //Han-Yi Chou
Maxim Belkin
committed
cmaffeo2
committed
//I need kT here Han-Yi Chou
kT = temperature * 0.0019872065f; // `units "k K" "kcal_mol"`
//kT = temperature * 0.593f;
Maxim Belkin
committed
// If we have a restart file - use it
if (restartCoordinates.length() > 0) {
cmaffeo2
committed
loadRestart(restartCoordinates.val());
Maxim Belkin
committed
printf("Loaded %d restart coordinates from `%s'.\n", num, restartCoordinates.val());
printf("Particle numbers specified in the configuration file will be ignored.\n");
loadedCoordinates = true;
cmaffeo2
committed
//Han-Yi Chou Langevin dynamic
if (ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
{
if (restartMomentum.length() > 0)
{
loadRestartMomentum(restartMomentum.val());
printf("Loaded %d restart momentum from `%s'.\n", num, restartMomentum.val());
printf("Particle numbers specified in the configuration file will be ignored.\n");
loadedMomentum = true;
}
else
{
printf("Warning: There is no restart momentum file when using restart coordinates in Langevin Dynamics\n");
printf("Initialize with Boltzmann distribution\n");
loadedMomentum = Boltzmann(COM_Velocity, num * simNum);
}
}
}
else
{
Maxim Belkin
committed
// Load coordinates from a file?
if (numPartsFromFile > 0) {
loadedCoordinates = true;
for (int i = 0; i < num; i++) {
int numTokens = partsFromFile[i].tokenCount();
// Break the line down into pieces (tokens) so we can process them individually
String* tokenList = new String[numTokens];
partsFromFile[i].tokenize(tokenList);
int currType = 0;
for (int j = 0; j < numParts; j++)
if (tokenList[2] == part[j].name)
currType = j;
for (int s = 0; s < simNum; ++s)
cmaffeo2
committed
type[i + s*num] = currType;
Maxim Belkin
committed
serial[i] = currSerial++;
cmaffeo2
committed
pos[i] = Vector3(atof(tokenList[3].val()), atof(tokenList[4].val()), atof(tokenList[5].val()));
//Han-Yi Chou
if (ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
{
loadedMomentum = true;
if(numTokens == 9)
momentum[i] = Vector3(atof(tokenList[6].val()), atof(tokenList[7].val()), atof(tokenList[8].val()));
else
{
printf("Error occurs in %s at line %d. Please specify momentum\n", __FILE__, __LINE__);
assert(1==2);
}
}
Maxim Belkin
committed
}
delete[] partsFromFile;
partsFromFile = NULL;
cmaffeo2
committed
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
//Han-Yi Chou
for(int i = 1; i < simNum; ++i)
for(int j = 0; j < num; ++j)
serial[j + num * i] = currSerial++;
}
else
{
// Not loading coordinates from a file
populate();
if (inputCoordinates.length() > 0)
{
printf("Loading coordinates from %s ... ", inputCoordinates.val());
loadedCoordinates = loadCoordinates(inputCoordinates.val());
if (loadedCoordinates)
printf("done!\n");
}
if(ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
{
if (inputMomentum.length() > 0)
{
printf("Loading momentum from %s ... ", inputMomentum.val());
loadedMomentum = loadMomentum(inputMomentum.val());
if (loadedMomentum)
printf("done!\n");
}
else
loadedMomentum = Boltzmann(COM_Velocity, num * simNum);
}
}
}
//Check initialize momentum
//if(ParticleDynamicType == String("Langevin"))
//PrintMomentum();
cmaffeo2
committed
/* Initialize exclusions */
excludeCapacity = 256;
numExcludes = 0;
excludes = new Exclude[excludeCapacity];
cmaffeo2
committed
if (readBondsFromFile) readBonds();
if (readAnglesFromFile) readAngles();
if (readDihedralsFromFile) readDihedrals();
if (readRestraintsFromFile) readRestraints();
if (temperatureGridFile.length() != 0) {
printf("\nFound temperature grid file: %s\n", temperatureGridFile.val());
tGrid = new BaseGrid(temperatureGridFile.val());
printf("Loaded `%s'.\n", temperatureGridFile.val());
printf("Grid size %s.\n", tGrid->getExtent().toString().val());
// TODO: ask Max Belkin what this is about and how to remove hard-coded temps
Maxim Belkin
committed
float ToSo = 1.0f / (295.0f * 4.634248239f); // 1 / (To * sigma(To))
sigmaT = new BaseGrid(*tGrid);
sigmaT->shift(-122.8305f);
sigmaT->scale(0.0269167f);
sigmaT->mult(*tGrid);
sigmaT->scale(ToSo);
kTGrid = new BaseGrid(*tGrid);
float factor = 0.0019872065f; // `units "k K" "kcal_mol"`
kTGrid->scale(factor);
// char outFile[256];
// char comment[256]; sprintf(comment,"KTGrid");
// sprintf(outFile,"kTGrid.dx");
// kTGrid->write(outFile, comment);
}
printf("\nFound %d particle types.\n", numParts);
// Load the potential grids.
printf("Loading the potential grids...\n");
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
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
319
for (int i = 0; i < numParts; i++)
{
String map = partGridFile[i][0];
int len = map.length();
if (len >= 3 && map[len-3]=='.' && map[len-2]=='d' && map[len-1]=='x')
{
part[i].pmf = new BaseGrid[part[i].numPartGridFiles];
part[i].meanPmf = new float[part[i].numPartGridFiles];
for(int j = 0; j < part[i].numPartGridFiles; ++j)
{
map = partGridFile[i][j];
len = map.length();
if (!(len >= 3 && map[len-3]=='.' && map[len-2]=='d' && map[len-1]=='x'))
{
cout << "currently do not support different format " << endl;
exit(1);
}
// Decide which type of grid is given.
//String map = partGridFile[i];
// A dx file. Load the old-fashioned way.
//part[i].pmf[j] = new BaseGrid(map.val());
BaseGrid tmp(map.val());
part[i].pmf[j] = tmp;
if (partGridFileScale[i][j] != 1.0f)
part[i].pmf[j].scale(partGridFileScale[i][j]);
//part[i].meanPmf = part[i].pmf->mean();
part[i].meanPmf[j] = part[i].pmf[j].mean();
printf("Loaded dx potential grid `%s'.\n", map.val());
printf("Grid size %s.\n", part[i].pmf[j].getExtent().toString().val());
}
}
else if (len >= 4 && map[len-4]=='.' && map[len-3]=='d' && map[len-2]=='e' && map[len-1]=='f')
{
OverlordGrid* over = new OverlordGrid[part[i].numPartGridFiles];
part[i].meanPmf = new float[part[i].numPartGridFiles];
for(int j = 0; j < part[i].numPartGridFiles; ++j)
{
map = partGridFile[i][j];
len = map.length();
if (!(len >= 4 && map[len-4]=='.' && map[len-3]=='d' && map[len-2]=='e' && map[len-1]=='f'))
{
cout << "currently do not support different format " << endl;
exit(1);
}
String rootGrid = OverlordGrid::readDefFirst(map);
over[j] = OverlordGrid(rootGrid.val());
int count = over->readDef(map);
printf("Loaded system def file `%s'.\n", map.val());
printf("Found %d unique grids.\n", over->getUniqueGridNum());
printf("Linked %d subgrids.\n", count);
part[i].meanPmf[j] = part[i].pmf[j].mean();
}
part[i].pmf = static_cast<BaseGrid*>(over);
}
else
{
printf("WARNING: Unrecognized gridFile extension. Must be *.def or *.dx.\n");
exit(-1);
}
if (partForceXGridFile[i].length() != 0) {
part[i].forceXGrid = new BaseGrid(partForceXGridFile[i].val());
printf("Loaded `%s'.\n", partForceXGridFile[i].val());
printf("Grid size %s.\n", part[i].forceXGrid->getExtent().toString().val());
}
if (partForceYGridFile[i].length() != 0) {
part[i].forceYGrid = new BaseGrid(partForceYGridFile[i].val());
printf("Loaded `%s'.\n", partForceYGridFile[i].val());
printf("Grid size %s.\n", part[i].forceYGrid->getExtent().toString().val());
Maxim Belkin
committed
}
if (partForceZGridFile[i].length() != 0) {
part[i].forceZGrid = new BaseGrid(partForceZGridFile[i].val());
printf("Loaded `%s'.\n", partForceZGridFile[i].val());
printf("Grid size %s.\n", part[i].forceZGrid->getExtent().toString().val());
}
if (partDiffusionGridFile[i].length() != 0) {
part[i].diffusionGrid = new BaseGrid(partDiffusionGridFile[i].val());
printf("Loaded `%s'.\n", partDiffusionGridFile[i].val());
printf("Grid size %s.\n", part[i].diffusionGrid->getExtent().toString().val());
if (temperatureGridFile.length() != 0) {
if (partDiffusionGridFile[i].length() != 0) {
part[i].diffusionGrid->mult(*sigmaT);
} else {
part[i].diffusionGrid = new BaseGrid(*sigmaT);
part[i].diffusionGrid->scale(part[i].diffusion);
// char outFile[256];
// char comment[256]; sprintf(comment,"Diffusion for particle type %d", i);
// sprintf(outFile,"diffusion%d.dx",i);
// part[i].diffusionGrid->write(outFile, comment);
}
}
}
// Load reservoir files if any
for (int i = 0; i < numParts; i++) {
if (partReservoirFile[i].length() != 0) {
printf("\nLoading the reservoirs for %s... \n", part[i].name.val());
part[i].reservoir = new Reservoir(partReservoirFile[i].val());
int nRes = part[i].reservoir->length();
printf("\t -> %d reservoir(s) found in `%s'.\n", nRes, partReservoirFile[i].val());
}
}
// Get the system dimensions
// from the dimensions of supplied 3D potential maps
if (size.length2() > 0) { // use size if it's defined
if (basis1.length2() > 0 || basis2.length2() > 0 || basis3.length2() > 0)
printf("WARNING: both 'size' and 'basis' were specified... using 'size'\n");
basis1 = Vector3(size.x,0,0);
basis2 = Vector3(0,size.y,0);
basis3 = Vector3(0,0,size.z);
}
if (basis1.length2() > 0 && basis2.length2() > 0 && basis3.length2() > 0) {
sys = new BaseGrid( Matrix3(basis1,basis2,basis3), origin, 1, 1, 1 );
} else {
// TODO: use largest system in x,y,z
}
sysDim = sys->getExtent();
// RBTODO: clean this mess up
/* // RigidBodies... */
/* if (numRigidTypes > 0) { */
/* printf("\nCounting rigid bodies specified in the configuration file.\n"); */
/* numRB = 0; */
/* // grow list of rbs */
/* for (int i = 0; i < numRigidTypes; i++) { */
/* numRB += rigidBody[i].num; */
/* std::vector<RigidBody> tmp; */
/* for (int j = 0; j < rigidBody[i].num; j++) { */
/* tmp.push_back( new RigidBody( this, rigidBody[i] ) ); */
/* } */
/* rbs.push_back(tmp); */
/* } */
// // state data
// rbPos = new Vector3[numRB * simNum];
// type = new int[numRB * simNum];
Maxim Belkin
committed
/* } */
/* printf("Initial RigidBodies: %d\n", numRB); */
Maxim Belkin
committed
// Create exclusions from the exclude rule, if it was specified in the config file
if (excludeRule != String("")) {
int oldNumExcludes = numExcludes;
Exclude* newExcludes = makeExcludes(bonds, bondMap, num, numBonds, excludeRule, numExcludes);
if (excludes == NULL) {
Maxim Belkin
committed
excludes = new Exclude[numExcludes];
} else if (numExcludes >= excludeCapacity) {
Exclude* tempExcludes = excludes;
excludes = new Exclude[numExcludes];
for (int i = 0; i < oldNumExcludes; i++)
excludes[i] = tempExcludes[i];
delete tempExcludes;
Maxim Belkin
committed
}
for (int i = oldNumExcludes; i < numExcludes; i++)
excludes[i] = newExcludes[i - oldNumExcludes];
printf("Built %d exclusions.\n",numExcludes);
}
cmaffeo2
committed
buildExcludeMap();
// Count number of particles of each type
numPartsOfType = new int[numParts];
for (int i = 0; i < numParts; ++i) {
numPartsOfType[i] = 0;
}
for (int i = 0; i < num; ++i) {
++numPartsOfType[type[i]];
}
// Some geometric stuff that should be gotten rid of.
Vector3 buffer = (sys->getCenter() + 2.0f*sys->getOrigin())/3.0f;
initialZ = buffer.z;
// Set the initial conditions.
// Do the initial conditions come from restart coordinates?
// inputCoordinates are ignored if restartCoordinates exist.
/*
if (restartCoordinates.length() > 0) {
Maxim Belkin
committed
loadRestart(restartCoordinates.val());
printf("Loaded %d restart coordinates from `%s'.\n", num, restartCoordinates.val());
printf("Particle numbers specified in the configuration file will be ignored.\n");
} else {
// Set the particle types.
// Load coordinates from a file?
if (numPartsFromFile > 0) {
Maxim Belkin
committed
for (int i = 0; i < num; i++) {
Maxim Belkin
committed
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
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
// Break the line down into pieces (tokens) so we can process them individually
String* tokenList = new String[numTokens];
partsFromFile[i].tokenize(tokenList);
int currType = 0;
for (int j = 0; j < numParts; j++)
if (tokenList[2] == part[j].name)
currType = j;
type[i] = currType;
serial[i] = currSerial;
currSerial++;
pos[i] = Vector3(atof(tokenList[3].val()), atof(tokenList[4].val()), atof(tokenList[5].val()));
}
if (partsFromFile != NULL) {
delete[] partsFromFile;
partsFromFile = NULL;
}
} else if (inputCoordinates.length() > 0) {
populate();
printf("Loading coordinates from %s.\n", inputCoordinates.val());
bool loaded = loadCoordinates(inputCoordinates.val());
if (loaded)
printf("Loaded initial coordinates from %s.\n", inputCoordinates.val());
}
}
*/
// Get the maximum particle radius.
minimumSep = 0.0f;
for (int i = 0; i < numParts; ++i)
minimumSep = std::max(minimumSep, part[i].radius);
minimumSep *= 2.5f; // Make it a little bigger.
// Default outputEnergyPeriod
if (outputEnergyPeriod < 0)
outputEnergyPeriod = 10 * outputPeriod;
// If we are running with debug ON, ask the user which force computation to use
if (debug)
getDebugForce();
printf("\n");
switchStart = cutoff - switchLen;
if (fullLongRange == 0)
printf("Cutting off the potential from %.10g to %.10g.\n", switchStart, switchStart+switchLen);
if (fullLongRange != 0)
printf("No cell decomposition created.\n");
cmaffeo2
committed
}
Configuration::~Configuration() {
// System state
delete[] pos;
cmaffeo2
committed
//Han-Yi Chou
if (ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
delete[] momentum;
cmaffeo2
committed
//Han-Yi Chou
if (ParticleDynamicType == String("Langevin") || ParticleDynamicType == String("NoseHooverLangevin"))
delete[] momLast;
delete[] type;
delete[] name;
// Particle parameters
delete[] part;
//delete[] partGridFile;
//delete[] partGridFileScale;
for(int i = 0; i < numParts; ++i)
{
if(partGridFile[i] != NULL)
{
delete[] partGridFile[i];
partGridFile[i] = NULL;
}
if(partGridFileScale[i] != NULL)
{
delete[] partGridFileScale[i];
partGridFileScale[i] = NULL;
}
}
delete partGridFile;
delete partGridFileScale;
//delete numPartGridFiles;
delete[] partForceXGridFile;
delete[] partForceYGridFile;
delete[] partForceZGridFile;
delete[] partDiffusionGridFile;
delete[] partReservoirFile;
partRigidBodyGrid.clear();
// TODO: plug memory leaks
if (partsFromFile != NULL) delete[] partsFromFile;
if (bonds != NULL) delete[] bonds;
if (bondMap != NULL) delete[] bondMap;
if (excludes != NULL) delete[] excludes;
if (excludeMap != NULL) delete[] excludeMap;
if (angles != NULL) delete[] angles;
if (dihedrals != NULL) delete[] dihedrals;
delete[] numPartsOfType;
// Table parameters
delete[] partTableFile;
delete[] partTableIndex0;
delete[] partTableIndex1;
delete[] bondTableFile;
delete[] angleTableFile;
delete[] dihedralTableFile;
eguzman6
committed
//if (type_d != NULL) {
//gpuErrchk(cudaFree(type_d));
gpuErrchk(cudaFree(sys_d));
gpuErrchk(cudaFree(kTGrid_d));
gpuErrchk(cudaFree(part_d));
eguzman6
committed
//gpuErrchk(cudaFree(bonds_d));
//gpuErrchk(cudaFree(bondMap_d));
//gpuErrchk(cudaFree(excludes_d));
//gpuErrchk(cudaFree(excludeMap_d));
//gpuErrchk(cudaFree(angles_d));
//gpuErrchk(cudaFree(dihedrals_d));
//}
printf("Copying particle data to GPU %d\n", GPUManager::current());
BrownianParticleType **part_addr = new BrownianParticleType*[numParts];
// Copy the BaseGrid objects and their member variables/objects
gpuErrchk(cudaMalloc(&part_d, sizeof(BrownianParticleType*) * numParts));
// TODO: The above line fails when there is not enough memory. If it fails, stop.
for (int i = 0; i < numParts; i++)
{
BaseGrid *pmf = NULL, *diffusionGrid = NULL;
BrownianParticleType *b = new BrownianParticleType(part[i]);
// Copy PMF
if (part[i].pmf != NULL)
{
float *tmp;
gpuErrchk(cudaMalloc(&pmf, sizeof(BaseGrid)*part[i].numPartGridFiles));
gpuErrchk(cudaMalloc(&tmp, sizeof(float)*part[i].numPartGridFiles));
gpuErrchk(cudaMemcpy(tmp, part[i].meanPmf, sizeof(float)*part[i].numPartGridFiles,
cudaMemcpyHostToDevice));
b->meanPmf = tmp;
for(int j = 0; j < part[i].numPartGridFiles; ++j)
{
float *val = NULL;
size_t sz = sizeof(float) * part[i].pmf[j].getSize();
//gpuErrchk(cudaMalloc(pmf, sizeof(BaseGrid)));
gpuErrchk(cudaMalloc(&val, sz));
gpuErrchk(cudaMemcpyAsync(val, part[i].pmf[j].val, sz, cudaMemcpyHostToDevice));
BaseGrid *pmf_h = new BaseGrid(part[i].pmf[j]);
pmf_h->val = val;
gpuErrchk(cudaMemcpy(pmf+j, pmf_h, sizeof(BaseGrid), cudaMemcpyHostToDevice));
pmf_h->val = NULL;
}
b->pmf = pmf;
// Copy the diffusion grid
if (part[i].diffusionGrid != NULL) {
float *val = NULL;
size_t sz = sizeof(float) * part[i].diffusionGrid->getSize();
BaseGrid *diffusionGrid_h = new BaseGrid(*part[i].diffusionGrid);
gpuErrchk(cudaMalloc(&diffusionGrid, sizeof(BaseGrid)));
gpuErrchk(cudaMalloc(&val, sz));
diffusionGrid_h->val = val;
gpuErrchk(cudaMemcpyAsync(diffusionGrid, diffusionGrid_h, sizeof(BaseGrid),
cudaMemcpyHostToDevice));
gpuErrchk(cudaMemcpy(val, part[i].diffusionGrid->val, sz, cudaMemcpyHostToDevice));
diffusionGrid_h->val = NULL;
}
//b->pmf = pmf;
b->diffusionGrid = diffusionGrid;
gpuErrchk(cudaMalloc(&part_addr[i], sizeof(BrownianParticleType)));
gpuErrchk(cudaMemcpyAsync(part_addr[i], b, sizeof(BrownianParticleType),
cudaMemcpyHostToDevice));
}
// RBTODO: moved this out of preceding loop; was that correct?
gpuErrchk(cudaMemcpyAsync(part_d, part_addr, sizeof(BrownianParticleType*) * numParts,
cudaMemcpyHostToDevice));
if (temperatureGridFile.length() > 0) {
gpuErrchk(cudaMalloc(&kTGrid_d, sizeof(BaseGrid)));
gpuErrchk(cudaMemcpyAsync(kTGrid_d, kTGrid, sizeof(BaseGrid), cudaMemcpyHostToDevice));
}
// type_d and sys_d
gpuErrchk(cudaMalloc(&sys_d, sizeof(BaseGrid)));
gpuErrchk(cudaMemcpyAsync(sys_d, sys, sizeof(BaseGrid), cudaMemcpyHostToDevice));
eguzman6
committed
/*gpuErrchk(cudaMalloc(&type_d, sizeof(int) * num * simNum));
gpuErrchk(cudaMemcpyAsync(type_d, type, sizeof(int) * num * simNum, cudaMemcpyHostToDevice));
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
if (numBonds > 0) {
// bonds_d
gpuErrchk(cudaMalloc(&bonds_d, sizeof(Bond) * numBonds));
gpuErrchk(cudaMemcpyAsync(bonds_d, bonds, sizeof(Bond) * numBonds, cudaMemcpyHostToDevice));
// bondMap_d
gpuErrchk(cudaMalloc(&bondMap_d, sizeof(int2) * num));
gpuErrchk(cudaMemcpyAsync(bondMap_d, bondMap, sizeof(int2) * num, cudaMemcpyHostToDevice));
}
if (numExcludes > 0) {
// excludes_d
gpuErrchk(cudaMalloc(&excludes_d, sizeof(Exclude) * numExcludes));
gpuErrchk(cudaMemcpyAsync(excludes_d, excludes, sizeof(Exclude) * numExcludes,
cudaMemcpyHostToDevice));
// excludeMap_d
gpuErrchk(cudaMalloc(&excludeMap_d, sizeof(int2) * num));
gpuErrchk(cudaMemcpyAsync(excludeMap_d, excludeMap, sizeof(int2) * num,
cudaMemcpyHostToDevice));
}
if (numAngles > 0) {
// angles_d
gpuErrchk(cudaMalloc(&angles_d, sizeof(Angle) * numAngles));
gpuErrchk(cudaMemcpyAsync(angles_d, angles, sizeof(Angle) * numAngles,
cudaMemcpyHostToDevice));
}
if (numDihedrals > 0) {
// dihedrals_d
gpuErrchk(cudaMalloc(&dihedrals_d, sizeof(Dihedral) * numDihedrals));
gpuErrchk(cudaMemcpyAsync(dihedrals_d, dihedrals,
sizeof(Dihedral) * numDihedrals,
cudaMemcpyHostToDevice));
eguzman6
committed
}*/
gpuErrchk(cudaDeviceSynchronize());
}
void Configuration::setDefaults() {
// System parameters
cmaffeo2
committed
rigidBodyGridGridPeriod = 1;
unsigned long int r0 = clock();
for (int i = 0; i < 4; i++)
r0 *= r0 + 1;
seed = time(NULL) + r0;
origin = Vector3(0,0,0);
size = Vector3(0,0,0);
basis1 = Vector3(0,0,0);
basis2 = Vector3(0,0,0);
basis3 = Vector3(0,0,0);
inputCoordinates = "";
restartCoordinates = "";
cmaffeo2
committed
//Han-Yi Chou
inputMomentum = "";
restartMomentum = "";
copyReplicaCoordinates = 1;
numberFluct = 0;
numberFluctPeriod = 200;
interparticleForce = 1;
tabulatedPotential = 0;
fullLongRange = 0;
// kTGridFile = ""; // Commented out for an unknown reason
temperature = 295.0f;
temperatureGridFile = "";
coulombConst = 566.440698f/92.0f;
electricField = 0.0f;
cutoff = 10.0f;
switchLen = 2.0f;
pairlistDistance = 2.0f;
imdForceScale = 1.0f;
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
outputPeriod = 200;
outputEnergyPeriod = -1;
outputFormat = TrajectoryWriter::formatDcd;
currentSegmentZ = -1.0f;
numCap = 0;
decompPeriod = 10;
readPartsFromFile = 0;
numPartsFromFile = 0;
partsFromFile = NULL;
readBondsFromFile = false;
numBonds = 0;
bonds = NULL;
bondMap = NULL;
numTabBondFiles = 0;
readExcludesFromFile = false;
numExcludes = 0;
excludeCapacity = 256;
excludes = NULL;
excludeMap = NULL;
excludeRule = "";
readAnglesFromFile = false;
numAngles = 0;
angles = NULL;
numTabAngleFiles = 0;
readDihedralsFromFile = false;
numDihedrals = 0;
dihedrals = NULL;
numTabDihedralFiles = 0;
readRestraintsFromFile = false;
numRestraints = 0;
restraints = NULL;
cmaffeo2
committed
//Han-Yi Chou default values
ParticleDynamicType = String("Brown");
RigidBodyDynamicType = String("Brown");
COM_Velocity = Vector3(0.f,0.f,0.f);
ParticleLangevinIntegrator = String("BAOAB"); //The default is BAOAB
// Hidden parameters
// Might be parameters later
numCapFactor = 5;
ParticleInterpolationType = 0;
RigidBodyInterpolationType = 0;
}
int Configuration::readParameters(const char * config_file) {
Reader config(config_file);
printf("Read config file %s\n", config_file);
// Get the number of particles.
const int numParams = config.length();
numParts = config.countParameter("particle");
numRigidTypes = config.countParameter("rigidBody");
// Allocate the particle variables.
part = new BrownianParticleType[numParts];
//partGridFile = new String[numParts];
//partGridFileScale = new float[numParts];
partGridFile = new String*[numParts];
//partGridFileScale = new float[numParts];
partGridFileScale = new float*[numParts];
//int numPartGridFiles = new int[numParts];
partForceXGridFile = new String[numParts];
partForceYGridFile = new String[numParts];
partForceZGridFile = new String[numParts];
partDiffusionGridFile = new String[numParts];
partReservoirFile = new String[numParts];
partRigidBodyGrid.resize(numParts);
// Allocate the table variables.
partTableFile = new String[numParts*numParts];
partTableIndex0 = new int[numParts*numParts];
partTableIndex1 = new int[numParts*numParts];
// Allocate rigid body types
rigidBody = new RigidBodyType[numRigidTypes];
// Set a default
/*
for (int i = 0; i < numParts; ++i) {
partGridFileScale[i] = 1.0f;
}*/
for(int i = 0; i < numParts; ++i)
{
partGridFile[i] = NULL;
partGridFileScale[i] = NULL;
//part[i].numPartGridFiles = -1;
}
//for(int i = 0; i < numParts; ++i)
// cout << part[i].numPartGridFiles << endl;
int btfcap = 10;
bondTableFile = new String[btfcap];
int atfcap = 10;
angleTableFile = new String[atfcap];
int dtfcap = 10;
dihedralTableFile = new String[dtfcap];
int currPart = -1;
int currTab = -1;
int currBond = -1;
int currAngle = -1;
int currDihedral = -1;
int currRB = -1;
int partClassPart = 0;
int partClassRB = 1;
int currPartClass = -1; // 0 => particle, 1 => rigidBody
for (int i = 0; i < numParams; i++) {
String param = config.getParameter(i);
String value = config.getValue(i);
if (param == String("outputName"))
outputName = value;
else if (param == String("timestep"))
timestep = (float) strtod(value.val(), NULL);
cmaffeo2
committed
else if (param == String("rigidBodyGridGridPeriod"))
rigidBodyGridGridPeriod = atoi(value.val());
else if (param == String("steps"))
steps = atol(value.val());
else if (param == String("seed"))
seed = atoi(value.val());
else if (param == String("origin"))
origin = stringToVector3( value );
else if (param == String("systemSize"))
size = stringToVector3( value );
else if (param == String("basis1"))
basis1 = stringToVector3( value );
else if (param == String("basis2"))
basis2 = stringToVector3( value );
else if (param == String("basis3"))
basis3 = stringToVector3( value );
else if (param == String("inputCoordinates"))
inputCoordinates = value;
else if (param == String("restartCoordinates"))
restartCoordinates = value;
cmaffeo2
committed
//Han-Yi Chou
else if (param == String("inputMomentum"))
inputMomentum = value;
else if (param == String("restartMomentum"))
restartMomentum = value;
else if (param == String("copyReplicaCoordinates"))
copyReplicaCoordinates = atoi(value.val());
else if (param == String("temperature"))
temperature = (float) strtod(value.val(),NULL);
temperatureGridFile = value;
else if (param == String("numberFluct"))
numberFluct = atoi(value.val());
else if (param == String("numberFluctPeriod"))
numberFluctPeriod = atoi(value.val());
else if (param == String("interparticleForce"))
interparticleForce = atoi(value.val());
else if (param == String("fullLongRange") || param == String("fullElect") )
fullLongRange = atoi(value.val());
else if (param == String("coulombConst"))
coulombConst = (float) strtod(value.val(), NULL);
else if (param == String("electricField"))
electricField = (float) strtod(value.val(), NULL);
else if (param == String("cutoff"))
cutoff = (float) strtod(value.val(), NULL);
else if (param == String("switchLen"))
switchLen = (float) strtod(value.val(), NULL);
else if (param == String("pairlistDistance"))
pairlistDistance = (float) strtod(value.val(), NULL);
else if (param == String("scaleIMDForce"))
imdForceScale = (float) strtod(value.val(), NULL);
else if (param == String("outputPeriod"))
outputPeriod = atoi(value.val());
else if (param == String("outputEnergyPeriod"))
outputEnergyPeriod = atoi(value.val());
else if (param == String("outputFormat"))
outputFormat = TrajectoryWriter::getFormatCode(value);
else if (param == String("currentSegmentZ"))
currentSegmentZ = (float) strtod(value.val(), NULL);
else if (param == String("numCap"))
numCap = atoi(value.val());
else if (param == String("decompPeriod"))
decompPeriod = atoi(value.val());
cmaffeo2
committed
//Han-Yi Chou
else if (param == String("ParticleDynamicType"))
ParticleDynamicType = value;
else if (param == String("RigidBodyDynamicType"))
RigidBodyDynamicType = value;
else if (param == String("ParticleLangevinIntegrator"))
ParticleLangevinIntegrator = value;
else if (param == String("ParticleInterpolationType"))
ParticleInterpolationType = atoi(value.val());
else if (param == String("RigidBodyInterpolationType"))
RigidBodyInterpolationType = atoi(value.val());
// PARTICLES
else if (param == String("particle")) {
cmaffeo2
committed
else if (param == String("mu")) // for Nose-Hoover Langevin
part[currPart].mu = (float) strtod(value.val(), NULL);
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
else if (param == String("forceXGridFile"))
partForceXGridFile[currPart] = value;
else if (param == String("forceYGridFile"))
partForceYGridFile[currPart] = value;
else if (param == String("forceZGridFile"))
partForceZGridFile[currPart] = value;
else if (param == String("diffusionGridFile"))
partDiffusionGridFile[currPart] = value;
else if (param == String("diffusion"))
part[currPart].diffusion = (float) strtod(value.val(), NULL);
else if (param == String("charge"))
part[currPart].charge = (float) strtod(value.val(), NULL);
else if (param == String("radius"))
part[currPart].radius = (float) strtod(value.val(), NULL);
else if (param == String("eps"))
part[currPart].eps = (float) strtod(value.val(), NULL);
else if (param == String("reservoirFile"))
partReservoirFile[currPart] = value;
else if (param == String("tabulatedPotential"))
tabulatedPotential = atoi(value.val());
else if (param == String("tabulatedFile"))
readTableFile(value, ++currTab);
else if (param == String("tabulatedBondFile")) {
if (numTabBondFiles >= btfcap) {
String* temp = bondTableFile;
btfcap *= 2;
bondTableFile = new String[btfcap];
for (int j = 0; j < numTabBondFiles; j++)
bondTableFile[j] = temp[j];
delete[] temp;
}
if (readBondFile(value, ++currBond))
numTabBondFiles++;
} else if (param == String("inputParticles")) {
if (readPartsFromFile) {
printf("WARNING: More than one particle file specified. Ignoring new file.\n");
} else {
partFile = value;
readPartsFromFile = true;
loadedCoordinates = true;
}
} else if (param == String("inputBonds")) {
if (readBondsFromFile) {
printf("WARNING: More than one bond file specified. Ignoring new bond file.\n");
} else {
bondFile = value;
readBondsFromFile = true;
}
} else if (param == String("inputExcludes")) {
if (readExcludesFromFile) {
printf("WARNING: More than one exclude file specified. Ignoring new exclude file.\n");