Skip to content
Snippets Groups Projects
Configuration.cpp 46.7 KiB
Newer Older
  • Learn to ignore specific revisions
  • cmaffeo2's avatar
    cmaffeo2 committed
    1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 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 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 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 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 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 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 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 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 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 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000
    #include "Configuration.h"
    
    #define gpuErrchk(ans) { gpuAssert((ans), __FILE__, __LINE__); }
    inline void gpuAssert(cudaError_t code, 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);
       }
    }
    
    Configuration::Configuration(const char* config_file, int simNum, bool debug) :
    		simNum(simNum) {
    	// Read the parameters.
    	type_d = NULL;
    	kTGrid_d = NULL;
    	bonds_d = NULL;
    	bondMap_d = NULL;
    	excludes_d = NULL;
    	excludeMap_d = NULL;
    	angles_d = NULL;
    	dihedrals_d = NULL;
    	setDefaults();
    	readParameters(config_file);
    	if (readPartsFromFile) readAtoms();
    	if (readBondsFromFile) readBonds();
    	if (readExcludesFromFile) readExcludes();
    	if (readAnglesFromFile) readAngles();
    	if (readDihedralsFromFile) readDihedrals();
    	/*
    	if (kTGridFile.length() != 0) {
    		printf("\nFound kT grid file: %s\n", kTGridFile.val());
    		kTGrid = new BaseGrid(kTGridFile.val());
    		printf("Loaded `%s'.\n", kTGridFile.val());
    		printf("System size %s.\n", kTGrid->getExtent().toString().val());
    	}
    	*/
    	if (temperatureGrid.length() != 0) {
    		printf("\nFound temperature grid file: %s\n", temperatureGrid.val());
    		tGrid = new BaseGrid(temperatureGrid.val());
    		printf("Loaded `%s'.\n", temperatureGrid.val());
    		printf("System size %s.\n", tGrid->getExtent().toString().val());
    
    		float ToSo = 1.0f / (295.0f * 4.634248239f); // 1 / (To * sigma(To)) 
    		sigmaT = new BaseGrid(temperatureGrid.val());
    		sigmaT->shift(-122.8305f);
    		sigmaT->scale(0.0269167f);
    		sigmaT->mult(*tGrid);
    		sigmaT->scale(ToSo);
    		kTGrid = new BaseGrid(temperatureGrid.val());
    		float factor = kT/295.0f;
    		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");
    	for (int i = 0; i < numParts; i++) {
    		// Decide which type of grid is given.
    		String map = partGridFile[i];
    		int len = map.length();
    		if (len >= 3 && map[len-3]=='.' && map[len-2]=='d' && map[len-1]=='x') {
    			// A dx file. Load the old-fashioned way.
    			part[i].pmf = new BaseGrid(map.val());
    			part[i].meanPmf = part[i].pmf->mean();
    			printf("Loaded dx grid `%s'.\n", map.val());
    			printf("System size %s.\n", part[i].pmf->getExtent().toString().val());
    		} else if  (len >= 4 && map[len-4]=='.' && map[len-3]=='d' && map[len-2]=='e' && map[len-1]=='f') {
    			// A system definition file.
    			String rootGrid = OverlordGrid::readDefFirst(map);
    			OverlordGrid* over = new 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].pmf = static_cast<BaseGrid*>(over);
    			part[i].meanPmf = part[i].pmf->mean();
    		} 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("System 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("System size %s.\n", part[i].forceYGrid->getExtent().toString().val());
    		} 
    
    		if (partForceZGridFile[i].length() != 0) {
    			part[i].forceZGrid = new BaseGrid(partForceZGridFile[i].val());
    			printf("Loaded `%s'.\n", partForceZGridFile[i].val());
    			printf("System 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("System size %s.\n", part[i].diffusionGrid->getExtent().toString().val());
    		}
    
    		if (temperatureGrid.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
    	sys = part[0].pmf;
    	sysDim = part[0].pmf->getExtent();
    
    
    	// Get the initial number of particles.
        //
        printf("\nCounting particles specified in the ");
    	if (restartCoordinates.length() > 0) {
            // Read them from the restart file.
    		printf("restart file.\n");
    		num = countRestart(restartCoordinates.val());
    	} else if (numPartsFromFile == 0) {
            // Sum up all particles in config file
    		printf("configuration file.\n");
    		//int num0 = 0;
            num = 0;
    		for (int i = 0; i < numParts; i++) num += part[i].num;
    		//num = num0;	
    	} else {
            // Determine number of particles from input file (PDB-style) 
            printf("input file.\n");
    		num = numPartsFromFile;
    	}
    
    	// Set the number capacity.    
    	printf("\n");
    	printf("Initial particles: %d\n", num);
    	if (numCap <= 0) numCap = numCapFactor*num; // max number of particles
    	if (numCap <= 0) numCap = 20;
    
    	// Allocate particle variables.
    	pos = new Vector3[num * simNum];
    	type = new int[num * simNum];
    	serial = new int[num * simNum];
    	posLast = new Vector3[num * simNum];
    	name = new String[num * simNum];
    	currSerial = 0;
    
    	// 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) {
    			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;
    		}	
    		for (int i = oldNumExcludes; i < numExcludes; i++)
    			excludes[i] = newExcludes[i - oldNumExcludes];
    		printf("Built %d exclusions.\n",numExcludes);
    
    		// Call compareExcludeIndex with qsort to sort the excludes by BOTH ind1 AND ind2
    		std::sort(excludes, excludes + numExcludes, compare());
    
    		/* Each particle may have a varying number of excludes
    		 * excludeMap is an array with one element for each particle
    		 * which keeps track of where a particle's excludes are stored
    		 * in the excludes array.
    		 * excludeMap[i].x is the index in the excludes array where the ith particle's excludes begin
    		 * excludeMap[i].y is the index in the excludes array where the ith particle's excludes end
    		 */
    		excludeMap = new int2[numPartsFromFile];
    		for (int i = 0; i < numPartsFromFile; i++) {	
    			excludeMap[i].x = -1;
    			excludeMap[i].y = -1;
    		}
    		int currPart = -1;
    		int lastPart = -1;
    		for (int i = 0; i < numExcludes; i++) {
    			if (excludes[i].ind1 != currPart) {
    				currPart = excludes[i].ind1;
    				excludeMap[currPart].x = i;
    				if (lastPart >= 0)
    					excludeMap[lastPart].y = i;
    				lastPart = currPart;
    			}
    		}
    	}
    
    	// 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) {
    		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) {
    			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;
    				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());
    		}
    	}
    	*/
    	
    	loadedCoordinates = false;
    
      // If we have a restart file - use it
    	if (restartCoordinates.length() > 0) {
    		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");
    		loadedCoordinates = true;
    	} else {
    		// 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)
    					type[i + s*num] = currType;
    
    				serial[i] = currSerial++;
    
    				pos[i] = Vector3(atof(tokenList[3].val()),
    												 atof(tokenList[4].val()),
    												 atof(tokenList[5].val()));
    			}
    			delete[] partsFromFile;
    			partsFromFile = NULL;
    		} 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");
    			}
    		}
    	}
    	// 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");
    }
    
    Configuration::~Configuration() {
    	// System state
    	delete[] pos;
    	delete[] posLast;
    	delete[] type;
    	delete[] name;
    	
    	// Particle parameters
    	delete[] part;
    	delete[] partGridFile;
    	delete[] partForceXGridFile;
    	delete[] partForceYGridFile;
    	delete[] partForceZGridFile;
    	delete[] partDiffusionGridFile;
    	delete[] partReservoirFile;
    	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;
    
    	// Table parameters
    	delete[] partTableFile;
    	delete[] partTableIndex0;
    	delete[] partTableIndex1;
    
    	delete[] bondTableFile;
    
    	delete[] angleTableFile;
    
    	delete[] dihedralTableFile;
    
    	if (type_d != NULL) {
    		gpuErrchk(cudaFree(type_d));
    		gpuErrchk(cudaFree(sys_d));
    		gpuErrchk(cudaFree(kTGrid_d));
    		gpuErrchk(cudaFree(part_d));
    		gpuErrchk(cudaFree(bonds_d));
    		gpuErrchk(cudaFree(bondMap_d));
    		gpuErrchk(cudaFree(excludes_d));
    		gpuErrchk(cudaFree(excludeMap_d));
    		gpuErrchk(cudaFree(angles_d));
    		gpuErrchk(cudaFree(dihedrals_d));
    	}
    }
    
    void Configuration::copyToCUDA() {
    	printf("Copying 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 *val = NULL;
    			size_t sz = sizeof(float) * part[i].pmf->getSize();
    		  gpuErrchk(cudaMalloc(&pmf, sizeof(BaseGrid)));
    		  gpuErrchk(cudaMalloc(&val, sz));
    		  gpuErrchk(cudaMemcpyAsync(val, part[i].pmf->val, sz, cudaMemcpyHostToDevice));
    		  BaseGrid *pmf_h = new BaseGrid(*part[i].pmf);
    			pmf_h->val = val;
    			gpuErrchk(cudaMemcpy(pmf, pmf_h, sizeof(BaseGrid), cudaMemcpyHostToDevice));
    			pmf_h->val = NULL;
    		}
    		
    		// 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));
    		
    		gpuErrchk(cudaMemcpyAsync(part_d, part_addr, sizeof(BrownianParticleType*) * numParts,
    				cudaMemcpyHostToDevice));
    	}
    
    	// kTGrid_d
    	kTGrid_d = NULL;
    	if (temperatureGrid.length() > 0) {
    		gpuErrchk(cudaMalloc(&kTGrid_d, sizeof(BaseGrid)));
    		gpuErrchk(cudaMemcpyAsync(kTGrid_d, kTGrid, sizeof(BaseGrid), cudaMemcpyHostToDevice));
    	}
    
    	// type_d and sys_d
    	gpuErrchk(cudaMalloc(&type_d, sizeof(int) * num * simNum));
    	gpuErrchk(cudaMemcpyAsync(type_d, type, sizeof(int) * num * simNum, cudaMemcpyHostToDevice));
    	gpuErrchk(cudaMalloc(&sys_d, sizeof(BaseGrid)));
    	gpuErrchk(cudaMemcpyAsync(sys_d, sys, sizeof(BaseGrid), cudaMemcpyHostToDevice));
    
    	
    	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));
    	}
    	gpuErrchk(cudaDeviceSynchronize());
    }
    
    void Configuration::setDefaults() {
    	// System parameters
    	outputName = "out";
    	timestep = 1e-5f;
    	steps = 100;
    	seed = 0;
    	inputCoordinates = "";
    	restartCoordinates = "";
    	numberFluct = 0;
    	numberFluctPeriod = 200;
    	interparticleForce = 1;
    	tabulatedPotential = 0;
    	fullLongRange = 1;
    	kT = 1.0f;
    	//	kTGridFile = ""; // Commented out for an unknown reason
    	temperatureGrid = "";
    	coulombConst = 566.440698f/92.0f;
    	electricField = 0.0f;
    	cutoff = 10.0f;
    	switchLen = 2.0f;
    	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;
    
    	// Hidden parameters
    	// Might be parameters later
    	numCapFactor = 5;
    }
    
    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");
    
    	// Allocate the particle variables.
    	part = new BrownianParticleType[numParts];
    	partGridFile = new String[numParts];
    	partForceXGridFile = new String[numParts];
    	partForceYGridFile = new String[numParts];
    	partForceZGridFile = new String[numParts];
    	partDiffusionGridFile = new String[numParts];
    	partReservoirFile = new String[numParts];
    
    	// Allocate the table variables.
    	partTableFile = new String[numParts*numParts];
    	partTableIndex0 = new int[numParts*numParts];
    	partTableIndex1 = new int[numParts*numParts];
    	
    	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;
    	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);
    		else if (param == String("steps"))
    			steps = atol(value.val());
    		else if (param == String("seed"))
    			seed = atoi(value.val());
    		else if (param == String("inputCoordinates"))
    			inputCoordinates = value;
    		else if (param == String("restartCoordinates"))
    			restartCoordinates = value;
    		else if (param == String("kT"))
    			kT = (float) strtod(value.val(), NULL);
    		// else if (param == String("kTGridFile")) kTGridFile = value;
    		else if (param == String("temperatureGrid"))
    			temperatureGrid = 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("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());
    		else if (param == String("particle"))
    			part[++currPart] = BrownianParticleType(value);
    		else if (param == String("num"))
    			part[currPart].num = atoi(value.val());
    		else if (param == String("gridFile"))
    			partGridFile[currPart] = value;
    		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[i] = temp[i];
    				delete[] temp;
    			}
    			if (readBondFile(value, ++currBond))
    				numTabBondFiles++;
    		} else if (param == String("inputParticles")) {
    			if (readPartsFromFile) {
    				printf("WARNING: More than one particle file specified. Discarding 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. Discarding new bond file.\n");
    			} else {
    				bondFile = value;				
    				readBondsFromFile = true;
    			}
    		} else if (param == String("inputExcludes")) {
    			if (readExcludesFromFile) {
    				printf("WARNING: More than one exclude file specified. Discarding new exclude file.\n");
    			} else {
    				excludeFile = value;				
    				readExcludesFromFile = true;
    			}
    		} else if (param == String("exclude") or param == String("exclusion")) {
    			excludeRule = value; 
    		} else if (param == String("inputAngles")) {
    			if (readAnglesFromFile) {
    				printf("WARNING: More than one angle file specified. Discarding new angle file.\n");
    			} else {
    				angleFile = value;
    				readAnglesFromFile = true;
    			}
    		} else if (param == String("tabulatedAngleFile")) {
    			if (numTabAngleFiles >= atfcap) {
    				String* temp = angleTableFile;
    				atfcap *= 2;	
    				angleTableFile = new String[atfcap];
    				for (int j = 0; j < numTabAngleFiles; j++)
    					angleTableFile[i] = temp[i];
    				delete[] temp;
    			}
    			if (readAngleFile(value, ++currAngle))
    				numTabAngleFiles++;
    		} else if (param == String("inputDihedrals")) {
    			if (readDihedralsFromFile) {
    				printf("WARNING: More than one dihedral file specified. Discarding new dihedral file.\n");
    			} else {
    				dihedralFile = value;
    				readDihedralsFromFile = true;
    			}
    		} else if (param == String("tabulatedDihedralFile")) {
    			if (numTabDihedralFiles >= dtfcap) {
    				String * temp = dihedralTableFile;
    				dtfcap *= 2;
    				dihedralTableFile = new String[dtfcap];
    				for (int j = 0; j < numTabDihedralFiles; j++)
    					dihedralTableFile[i] = temp[i];
    				delete[] temp;
    			}
    			if (readDihedralFile(value, ++currDihedral))
    				numTabDihedralFiles++;
    		} else {
    			printf("WARNING: Unrecognized keyword `%s'.\n", param.val());
    		}
    	}
    	return numParams;
    }
    
    void Configuration::readAtoms() {
    	// Open the file
    	FILE* inp = fopen(partFile.val(), "r");
    	char line[256];
    
    	// If the particle file cannot be found, exit the program
    	if (inp == NULL) {
    		printf("ERROR: Could not open `%s'.\n", partFile.val());
    		bool found = true;
    		for (int i = 0; i < numParts; i++)
    			if (part[i].num == 0)
    				found = false;
    		if (!found) {
    			printf("ERROR: Number of particles not specified in config file.\n");
    			exit(1);
    		}
    		printf("Using default coordinates file\n");
    		return;
    	}
    
    	// Our particle array has a starting capacity of 256
    	// We will expand this later if we need to.
    	int capacity = 256;
    	numPartsFromFile = 0;
    	partsFromFile = new String[capacity];
    	indices = new int[capacity];
    	indices[0] = 0;
    
    	// Get and process all lines of input
    	while (fgets(line, 256, inp) != NULL) {
    		// Lines in the particle file that begin with # are comments
    		if (line[0] == '#') continue;
    		      
    		String s(line);
    		int numTokens = s.tokenCount();
    		      
    		// Break the line down into pieces (tokens) so we can process them individually
    		String* tokenList = new String[numTokens];
    		s.tokenize(tokenList);
    
    		// Legitimate ATOM input lines have 6 tokens: 
    		// ATOM | Index | Name | X-coord | Y-coord | Z-coord
    		// A line without exactly six tokens should be discarded.
    		if (numTokens != 6) {
    			printf("Warning: Invalid particle file line: %s\n", line);
    			return;
    		}
    
    		// Ensure that this particle's type was defined in the config file.
    		// If not, discard this line.
    		bool found;
    		for (int j = 0; j < numParts; j++) {
    			// If this particle type exists, add a new one to the list
    			if (part[j].name == tokenList[2]) {
    				found = true;
    				part[j].num++;
    			}
    		}
    
    		// If the particle's type does not exist according to the config file, discard it.
    		if (!found) {
    			printf("WARNING Unknown particle type %s found and discarded.\n", tokenList[2].val());
    			continue;
    		}
    
    		// If we don't have enough room in our particle array, we need to expand it.
    		if (numPartsFromFile >= capacity) {
    			// Temporary pointers to the old arrays
    			String* temp = partsFromFile;	
    			int* temp2 = indices;
    
    			// Double the capacity
    			capacity *= 2;
    
    			// Create pointers to new arrays which are twice the size of the old ones
    			partsFromFile = new String[capacity];
    			indices = new int[capacity];
    		
    			// Copy the old values into the new arrays
    			for (int j = 0; j < numPartsFromFile; j++) {
    				partsFromFile[j] = temp[j];
    				indices[j] = temp2[j];
    			}
    
    			// delete the old arrays
    			delete[] temp;
    			delete[] temp2;
    		}
    		// Make sure the index of this particle is unique.
    		// NOTE: The particle list is sorted by index. 
    		bool uniqueID = true;		
    		int key = atoi(tokenList[1].val());
    		int mid = 0;
    
    		// If the index is greater than the last index in the list, 
    		// this particle belongs at the end of the list. Since the 
    		// list is kept sorted, we know this is okay.
    		if (numPartsFromFile == 0 || key > indices[numPartsFromFile - 1]) {
    			indices[numPartsFromFile] = key;
    			partsFromFile[numPartsFromFile++] = line;
    		}
    		// We need to do a binary search to figure out if
    		// the index already exists in the list. 
    		// The assumption is that input files SHOULD have their indices sorted in 
    		// ascending order, so we shouldn't actually use the binary search 
    		// or the sort (which is pretty time consuming) very often.
    		else {
    			int low = 0, high = numPartsFromFile - 1;
    			
    			while (low <= high) {
    				mid = (int)((high - low) / 2 + low);
    				int curr = indices[mid];
    				if (curr < key) {
    					low = mid + 1;
    				} else if (curr > key) {
    					high = mid - 1;
    				} else {
    					// For now, particles with non-unique IDs are simply not added to the array
    					// Other possible approaches which are not yet implemented:
    					// 1: Keep track of these particles and assign them new IDs after you have
    					//    already added all of the other particles. 	
    					// 2: Get rid of ALL particles with that ID, even the ones that have already 
    					//    been added.
    					printf("WARNING: Non-unique ID found: %s\n", line);
    					uniqueID = false;
    					break;
    				}
    			}
    			if (uniqueID) {
    				// Add the particle to the end of the array, then sort it. 
    				indices[numPartsFromFile] = key;
    				partsFromFile[numPartsFromFile++] = line;
    				std::sort(indices, indices + numPartsFromFile);
    				std::sort(partsFromFile, partsFromFile + numPartsFromFile, compare());		
    			}
    		}
    	}
    }
    
    void Configuration::readBonds() {
    	// Open the file
    	FILE* inp = fopen(bondFile.val(), "r");
    	char line[256];
    
    	// If the particle file cannot be found, exit the program
    	if (inp == NULL) {
    		printf("WARNING: Could not open `%s'.\n", bondFile.val());
    		printf("         This simulation will not use particle bonds.\n");
    		return;
    	}
    
    	// Our particle array has a starting capacity of 256
    	// We will expand this later if we need to.
    	int capacity = 256;
    	numBonds = 0;
    	bonds = new Bond[capacity];
    
    	// Get and process all lines of input
    	while (fgets(line, 256, inp) != NULL) {
    		
    		// Lines in the particle file that begin with # are comments
    		if (line[0] == '#') continue;
    		      
    		String s(line);
    		int numTokens = s.tokenCount();
    		      
    		// Break the line down into pieces (tokens) so we can process them individually
    		String* tokenList = new String[numTokens];
    		s.tokenize(tokenList);
    
    		// Legitimate BOND input lines have 4 tokens: 
    		// BOND | OPERATION_FLAG | INDEX1 | INDEX2 | FILENAME 
    		// A line without exactly five tokens should be discarded.
    		if (numTokens != 5) {
    			printf("WARNING: Invalid bond file line: %s\n", line);
    			continue;
    		}
    
    		String op = tokenList[1];
    		int ind1 = atoi(tokenList[2].val());
    		int ind2 = atoi(tokenList[3].val());
    		String file_name = tokenList[4];
    		
    		// If we don't have enough room in our bond array, we need to expand it.
    		if (numBonds >= capacity) {
    			// Temporary pointer to the old array
    			Bond* temp = bonds;	
    
    			// Double the capacity
    			capacity *= 2;
    
    			// Create pointer to new array which is twice the size of the old one
    			bonds = new Bond[capacity];
    		
    			// Copy the old values into the new array
    			for (int j = 0; j < numBonds; j++)
    				bonds[j] = temp[j];
    
    			// delete the old array
    			delete[] temp;
    		}
    		// Add the bond to the bond array
    		// We must add it twice: Once for (ind1, ind2) and once for (ind2, ind1)
    		Bond* b = new Bond(op, ind1, ind2, file_name);
    		bonds[numBonds++] = *b;
    		b = new Bond(op, ind2, ind1, file_name);
    		bonds[numBonds++] = *b;
    		delete[] tokenList;
    	}	
    	// Call compareBondIndex with qsort to sort the bonds by BOTH ind1 AND ind2
    	std::sort(bonds, bonds + numBonds, compare());
    
    	/* Each particle may have a varying number of bonds
    	 * bondMap is an array with one element for each particle
    	 * which keeps track of where a particle's bonds are stored
    	 * in the bonds array.
    	 * bondMap[i].x is the index in the bonds array where the ith particle's bonds begin
    	 * bondMap[i].y is the index in the bonds array where the ith particle's bonds end
    	 */
    	bondMap = new int2[numPartsFromFile];
    	for (int i = 0; i < numPartsFromFile; i++) {	
    		bondMap[i].x = -1;
    		bondMap[i].y = -1;
    	}
    	int currPart = -1;
    	int lastPart = -1;
    	for (int i = 0; i < numBonds; i++) {
    		if (bonds[i].ind1 != currPart) {
    			currPart = bonds[i].ind1;
    			bondMap[currPart].x = i;
    			if (lastPart >= 0) bondMap[lastPart].y = i;
    			lastPart = currPart;
    		}
    	}
    	if (bondMap[lastPart].x > 0)
    		bondMap[lastPart].y = numBonds;
    }
    
    void Configuration::readExcludes()
    {
    	// Open the file
    	FILE* inp = fopen(excludeFile.val(), "r");
    	char line[256];
    
    	// If the exclusion file cannot be found, exit the program
    	if (inp == NULL) {
    		printf("WARNING: Could not open `%s'.\n", excludeFile.val());
    		printf("This simulation will not use exclusions.\n");
    		return;
    	}
    
    	// Our particle array has a starting capacity of 256
    	// We will expand this later if we need to.
    	excludeCapacity = 256;
    	numExcludes = 0;
    	excludes = new Exclude[excludeCapacity];
    
    	// Get and process all lines of input
    	while (fgets(line, 256, inp) != NULL) {
    		// Lines in the particle file that begin with # are comments
    		if (line[0] == '#') continue;
    		      
    		String s(line);
    		int numTokens = s.tokenCount();
    		      
    		// Break the line down into pieces (tokens) so we can process them individually
    		String* tokenList = new String[numTokens];
    		s.tokenize(tokenList);
    
    		// Legitimate EXCLUDE input lines have 3 tokens: