bil.h 57.6 KB
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 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542
#ifndef STIM_BIL_H
#define STIM_BIL_H

#include "../envi/envi_header.h"
#include "../envi/hsi.h"
#include "../math/fd_coefficients.h"
#include <stim/cuda/cudatools/error.h>
#include <cstring>
#include <utility>
#include <deque>

namespace stim{

/**
	The BIL class represents a 3-dimensional binary file stored using band interleaved by line (BIL) image encoding. The binary file is stored
	such that X-Z "frames" are stored sequentially to form an image stack along the y-axis. When accessing the data sequentially on disk,
	the dimensions read, from fastest to slowest, are X, Z, Y.

	This class is optimized for data streaming, and therefore supports extremely large (terabyte-scale) files. Data is loaded from disk
	on request. Functions used to access data are written to support efficient reading.
*/

template <typename T>

class bil: public hsi<T> {

protected:


	using hsi<T>::w;				//use the wavelength array in stim::hsi
	using hsi<T>::nnz;
	using binary<T>::progress;
	using hsi<T>::X;
	using hsi<T>::Y;
	using hsi<T>::Z;

public:

	using binary<T>::open;
	using binary<T>::file;
	using binary<T>::R;

	bil(){ hsi<T>::init_bil(); }

	/// Open a data file for reading using the class interface.

	/// @param filename is the name of the binary file on disk
	/// @param X is the number of samples along dimension 1
	/// @param Y is the number of samples (lines) along dimension 2
	/// @param B is the number of samples (bands) along dimension 3
	/// @param header_offset is the number of bytes (if any) in the binary header
	/// @param wavelengths is an optional STL vector of size B specifying a numerical label for each band
	bool open(std::string filename,
			  unsigned long long X,
			  unsigned long long Y,
			  unsigned long long B,
			  unsigned long long header_offset,
			  std::vector<double> wavelengths,
			  stim::iotype io = stim::io_in){

		w = wavelengths;

		return open(filename, vec<unsigned long long>(X, B, Y), header_offset, io);

	}

	/// Retrieve a single band (based on index) and stores it in pre-allocated memory.

	/// @param p is a pointer to an allocated region of memory at least X * Y * sizeof(T) in size.
	/// @param page <= B is the integer number of the band to be copied.
	bool band_index( T * p, unsigned long long page, bool PROGRESS = false){
		//return binary<T>::read_plane_1(p, page);

		if(PROGRESS) progress = 0;
		unsigned long long L = X() * sizeof(T);		//caculate the number of bytes in a sample line
		unsigned long long jump = X() * (Z() - 1) * sizeof(T);

		if (page >= Z()){										//make sure the bank number is right
			std::cout<<"ERROR: page out of range"<<std::endl;
			return false;
		}

		file.seekg(X() * page * sizeof(T), std::ios::beg);
		for (unsigned long long i = 0; i < Y(); i++)
		{
			file.read((char *)(p + i * X()), L);
			file.seekg( jump, std::ios::cur);
			if(PROGRESS) progress = (double)(i+1) / Y() * 100;
		}

		return true;
	}

	/// Retrieve a single band (by numerical label) and stores it in pre-allocated memory.

	/// @param p is a pointer to an allocated region of memory at least X * Y * sizeof(T) in size.
	/// @param wavelength is a floating point value (usually a wavelength in spectral data) used as a label for the band to be copied.
	bool band( T * p, double wavelength, bool PROGRESS = false){

		//if there are no wavelengths in the BSQ file
		if(w.size() == 0)
			return band_index(p, (unsigned long long)wavelength);

		unsigned long long XY = X() * Y();	//calculate the number of pixels in a band
		unsigned long long S = XY * sizeof(T);		//calculate the number of bytes of a band

		unsigned long long page=0;                      //bands around the wavelength


		//get the bands numbers around the wavelength

		//if wavelength is smaller than the first one in header file
		if ( w[page] > wavelength ){
			band_index(p, page, PROGRESS);
			return true;
		}

		while( w[page] < wavelength )
		{
			page++;
			//if wavelength is larger than the last wavelength in header file
			if (page == Z()) {
				band_index(p, Z()-1, PROGRESS);
				return true;
			}
		}
		if ( wavelength < w[page] )
		{
			T * p1;
			T * p2;
			p1=(T*)malloc(S);                     //memory allocation
			p2=(T*)malloc(S);
			band_index(p1, page - 1);
			band_index(p2, page, PROGRESS);
			for(unsigned long long i=0; i < XY; i++){
				double r = (wavelength - w[page-1]) / (w[page] - w[page-1]);
				p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]);
			}
			free(p1);
			free(p2);
		}
		else                           //if the wavelength is equal to a wavelength in header file
		{
			band_index(p, page, PROGRESS);
		}

		return true;
	}

	/// Retrieves a band of x values from a given xz plane.

	/// @param p is a pointer to pre-allocated memory at least Z * sizeof(T) in size
	/// @param c is a pointer to an existing XZ plane (size X*Z*sizeof(T))
	/// @param wavelength is the wavelength to retrieve
	bool read_x_from_xz(T* p, T* c, double wavelength)
	{
		unsigned long long B = Z();
		unsigned long long L = X() * sizeof(T);

		unsigned long long page=0;                      //samples around the wavelength
		T * p1;
		T * p2;

		//get the bands numbers around the wavelength

		//if wavelength is smaller than the first one in header file
		if ( w[page] > wavelength ){
			memcpy(p, c, L);
			return true;
		}

		while( w[page] < wavelength )
		{
			page++;
			//if wavelength is larger than the last wavelength in header file
			if (page == B) {
				memcpy(p, c + (B - 1) * X(), L);
				return true;
			}
		}
		if ( wavelength < w[page] )
		{
			p1=(T*)malloc( L );                     //memory allocation
			p2=(T*)malloc( L );

			memcpy(p1, c + (page - 1) * X(), L);
			memcpy(p2, c + page * X(), L);

			for(unsigned long long i=0; i < X(); i++){
				double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
				p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]);
			}
		}
		else                           //if the wavelength is equal to a wavelength in header file
			memcpy(p, c + page * X(), L);

		return true;
	}

	/// Retrieve a single spectrum (B-axis line) at a given (x, y) location and stores it in pre-allocated memory.

	/// @param p is a pointer to pre-allocated memory at least B * sizeof(T) in size.
	/// @param x is the x-coordinate (dimension 1) of the spectrum.
	/// @param y is the y-coordinate (dimension 2) of the spectrum.
	bool spectrum(T* p, size_t n, bool PROGRESS = false){
		size_t y = n / X();
		size_t x = n - y * X();
		return binary<T>::read_line_1(p, x, y, PROGRESS);
	}
	bool spectrum(T * p, unsigned long long x, unsigned long long y, bool PROGRESS = false){
		return binary<T>::read_line_1(p, x, y, PROGRESS);
	}

	/// Retrieve a single pixel and stores it in pre-allocated memory.

	/// @param p is a pointer to pre-allocated memory at least sizeof(T) in size.
	/// @param n is an integer index to the pixel using linear array indexing.
	bool pixel(T * p, unsigned long long n){

		//calculate the corresponding x, y
		unsigned long long x = n % X();
		unsigned long long y = n / X();

		//get the pixel
		return spectrum(p, x, y);
	}

	//given a Y ,return a XZ slice
	bool read_plane_xz(T * p, size_t y){
		return binary<T>::read_plane_2(p, y);
	}

	//given a Y, return ZX slice (transposed such that the spectrum is the leading dimension)
	bool read_plane_zx(T* p, size_t y){
		T* temp = (T*) malloc(X() * Z() * sizeof(T));		//allocate space to store the temporary xz plane
		if(!binary<T>::read_plane_2(temp, y))				//load the plane from disk
			return false;

		size_t z, x;
		for(z = 0; z < Z(); z++){
			for(x = 0; x <= z; x++){
				p[x * Z() + z] = temp[z * X() + x];		//copy to the destination frame
			}
		}
		return true;
	}

	//load a frame y into a pre-allocated double-precision array
	int read_plane_xzd(double* f, size_t y){		
		size_t XB = X() * Z();
		T* temp = (T*) malloc(XB * sizeof(T));			//create a temporary location to store the plane at current precision
		if(!read_plane_y(temp, y)) return 1;			//read the plane in its native format, if it fails return a 1
		for(size_t i = 0; i < XB; i++) f[i] = temp[i];	//convert the plane to a double
		return 0;
	}

	//given a Y, return ZX slice (transposed such that the spectrum is the leading dimension)
	int read_plane_zxd(double* p, size_t y){
		T* temp = (T*) malloc(X() * Z() * sizeof(T));		//allocate space to store the temporary xz plane
		binary<T>::read_plane_2(temp, y);					//load the plane from disk
		size_t z, x;
		for(z = 0; z < Z(); z++){
			for(x = 0; x < X(); x++){
				p[x * Z() + z] = (double)temp[z * X() + x];	//copy to the destination frame
			}
		}
		return 0;
	}


	/// Perform baseline correction given a list of baseline points and stores the result in a new BSQ file.

	/// @param outname is the name of the output file used to store the resulting baseline-corrected data.
	/// @param wls is the list of baseline points based on band labels.
	bool baseline(std::string outname, std::vector<double> wls, unsigned char* mask = NULL, bool PROGRESS = false){

		unsigned long long N = wls.size();			//get the number of baseline points

		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file
		std::string headername = outname + ".hdr";              //the header file name

		//simplify image resolution
		unsigned long long ZX = Z() * X();		//calculate the number of points in a Y slice
		unsigned long long L = ZX * sizeof(T);			//calculate the number of bytes of a Y slice
		unsigned long long B = Z();

		T* c;			//pointer to the current Y slice
		c = (T*)malloc(L);  //memory allocation

		T* a;			//pointer to the two YZ lines surrounding the current YZ line
		T* b;

		a = (T*)malloc(X() * sizeof(T));
		b = (T*)malloc(X() * sizeof(T));


		double ai, bi;	//stores the two baseline points wavelength surrounding the current band
		double ci;		//stores the current band's wavelength
		unsigned long long control;

		if (a == NULL || b == NULL || c == NULL){
			std::cout<<"ERROR: error allocating memory";
			exit(1);
		}
	//	loop start	correct every y slice

		for (unsigned long long k =0; k < Y(); k++)
		{
			//get the current y slice
			read_plane_xz(c, k);

			//initialize lownum, highnum, low, high
			ai = w[0];
			control=0;
			//if no baseline point is specified at band 0,
			//set the baseline point at band 0 to 0
			if(wls[0] != w[0]){
				bi = wls[control];
				memset(a, (char)0, X() * sizeof(T) );
			}
			//else get the low band
			else{
				control++;
				read_x_from_xz(a, c, ai);
				bi = wls[control];
			}
			//get the high band
			read_x_from_xz(b, c, bi);

			//correct every YZ line

			for(unsigned long long cii = 0; cii < B; cii++){

				//update baseline points, if necessary
				if( w[cii] >= bi && cii != B - 1) {
					//if the high band is now on the last BL point
					if (control != N-1) {

						control++;		//increment the index

						std::swap(a, b);	//swap the baseline band pointers

						ai = bi;
						bi = wls[control];
						read_x_from_xz(b, c, bi);

					}
					//if the last BL point on the last band of the file?
					else if ( wls[control] < w[B - 1]) {

						std::swap(a, b);	//swap the baseline band pointers

						memset(b, (char)0, X() * sizeof(T) );	//clear the high band

						ai = bi;
						bi = w[B - 1];
					}
				}

				ci = w[cii];

				unsigned long long jump = cii * X();
				//perform the baseline correction
				for(unsigned i=0; i < X(); i++)
				{
					double r = (double) (ci - ai) / (double) (bi - ai);
					c[i + jump] =(T) ( c[i + jump] - (b[i] - a[i]) * r - a[i] );
				}

			}//loop for YZ line end
			target.write(reinterpret_cast<const char*>(c), L);   //write the corrected data into destination

			if(PROGRESS) progress = (double)(k+1) / Y() * 100;
		}//loop for Y slice end

		free(a);
		free(b);
		free(c);
		target.close();

		return true;

	}

	/// Normalize all spectra based on the value of a single band, storing the result in a new BSQ file.

	/// @param outname is the name of the output file used to store the resulting baseline-corrected data.
	///	@param w is the label specifying the band that the hyperspectral image will be normalized to.
	///	@param t is a threshold specified such that a spectrum with a value at w less than t is set to zero. Setting this threshold allows the user to limit division by extremely small numbers.
	bool ratio(std::string outname, double w, unsigned char* mask = NULL, bool PROGRESS = false)
	{
		unsigned long long B = Z();		//calculate the number of bands
		unsigned long long ZX = Z() * X();
		unsigned long long XY = X() * Y();	//calculate the number of pixels in a band
		unsigned long long S = XY * sizeof(T);		//calculate the number of bytes in a band
		unsigned long long L = ZX * sizeof(T);

		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file
		std::string headername = outname + ".hdr";              //the header file name

		T * c;				//pointer to the current ZX slice
		T * b;				//pointer to the standard band

		b = (T*)malloc( S );     //memory allocation
		c = (T*)malloc( L );

		band(b, w);             //get the certain band into memory

		for(unsigned long long j = 0; j < Y(); j++)
		{
			read_plane_xz(c, j);
			for(unsigned long long i = 0; i < B; i++)
			{
				for(unsigned long long m = 0; m < X(); m++)
				{
					if( mask != NULL && !mask[m + j * X()] )
						c[m + i * X()] = (T)0.0;
					else
						c[m + i * X()] = c[m + i * X()] / b[m + j * X()];
				}
			}
			target.write(reinterpret_cast<const char*>(c), L);   //write normalized data into destination

			if(PROGRESS) progress = (double)(j+1) / Y() * 100;
		}

		free(b);
		free(c);
		target.close();
		return true;
	}

	void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){
		std::cout<<"ERROR: algorithm not implemented"<<std::endl;
		exit(1);
		//This code is almost done but has to be debugged
		/*std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file		

		size_t B = Z();									//calculate the number of bands
		size_t L = X();									//get the number of items in a line
		size_t Lb = L * sizeof(T);						//calculate the number of bytes in a line
		size_t XY = X() * Y();							//calculate the number of pixels in an image

		T* line = (T*) malloc(Lb);						//allocate space for a line
		T* len = (T*) malloc(Lb);						//allocate space for the lengths in a line

		for(size_t y = 0; y < Y(); y++){						//for each line in the image
			file.seekg(y * Lb, std::ios::beg);			//move the pointer to the current file to the beginning
			for(size_t b = 0; b < B; b++){				//for each band in this line
				file.read((char*)line, Lb);				//read a line of pixels
				for(size_t l = 0; l < L; l++)			//for each element in the line for this band
					len[l] += line[l] * line[l];		//add the square of the spectral component
			}
			for(size_t l = 0; l < L; l++)				//for each length element in the line
				len[l] = sqrt(len[l]);					//calculate the square root of the sum of squares
			
			file.seekg(y * Lb, std::ios::beg);			//move the pointer to the current file to the beginning
			for(size_t b = 0; b < B; b++){				//for each band in this line
				file.read((char*)line, Lb);				//read a line of pixels
				for(size_t l = 0; l < L; l++)			//for each element in the line for this band
					line[l] /= len[l];					//divide each element by the length
			}
			target.write((char*)line, Lb);				//write the normalized line to the target file
			if(PROGRESS) progress = (double)(y+1) / (double)Y();
		}*/
	}

	bool select(std::string outfile, std::vector<double> bandlist, unsigned char* mask = NULL, bool PROGRESS = NULL) {
		std::cout << "ERROR: select() not implemented for BIL" << std::endl;
		exit(1);
	}

	/// Convert the current BIL file to a BSQ file with the specified file name.

	/// @param outname is the name of the output BSQ file to be saved to disk.
	bool bsq(std::string outname, bool PROGRESS = false)
	{
		unsigned long long S = X() * Y() * sizeof(T);		//calculate the number of bytes in a band

		std::ofstream target(outname.c_str(), std::ios::binary);
		std::string headername = outname + ".hdr";

		T * p;			//pointer to the current band
		p = (T*)malloc(S);

		for ( unsigned long long i = 0; i < Z(); i++)
		{
				band_index(p, i);
				target.write(reinterpret_cast<const char*>(p), S);   //write a band data into target file

				if(PROGRESS) progress = (double)(i+1) / Z() * 100;	//store the progress for the current operation
		}

		free(p);
		target.close();
		return true;
	}

	/// Convert the current BIL file to a BIP file with the specified file name.

	/// @param outname is the name of the output BIP file to be saved to disk.
	bool bip(std::string outname, bool PROGRESS = false)
	{
		unsigned long long S = X() * Z() * sizeof(T);		//calculate the number of bytes in a ZX slice

		std::ofstream target(outname.c_str(), std::ios::binary);
		std::string headername = outname + ".hdr";

		T * p;			//pointer to the current XZ slice for bil file
		p = (T*)malloc(S);
		T * q;			//pointer to the current ZX slice for bip file
		q = (T*)malloc(S);

		for ( unsigned long long i = 0; i < Y(); i++)
		{
			read_plane_xz(p, i);
			for ( unsigned long long k = 0; k < Z(); k++)
			{
				unsigned long long ks = k * X();
				for ( unsigned long long j = 0; j < X(); j++)
					q[k + j * Z()] = p[ks + j];

				if(PROGRESS) progress = (double)((i+1) * Z() + k+1) / (Z() * Y()) * 100;	//store the progress for the current operation
			}

			target.write(reinterpret_cast<const char*>(q), S);   //write a band data into target file
		}

		free(p);
		free(q);
		target.close();
		return true;
	}


	/// Return a baseline corrected band given two adjacent baseline points and their bands. The result is stored in a pre-allocated array.

	/// @param lb is the label value for the left baseline point
	/// @param rb is the label value for the right baseline point
	/// @param lp is a pointer to an array holding the band image for the left baseline point
	/// @param rp is a pointer to an array holding the band image for the right baseline point
	/// @param wavelength is the label value for the requested baseline-corrected band
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size.
	bool baseline_band(double lb, double rb, T* lp, T* rp, double wavelength, T* result){

		unsigned long long XY = X() * Y();
		band(result, wavelength);		//get band

		//perform the baseline correction
		double r = (double) (wavelength - lb) / (double) (rb - lb);
		for(unsigned long long i=0; i < XY; i++){
			result[i] =(T) (result[i] - (rp[i] - lp[i]) * r - lp[i] );
		}
		return true;
	}

	/// Return a baseline corrected band given two adjacent baseline points. The result is stored in a pre-allocated array.

	/// @param lb is the label value for the left baseline point
	/// @param rb is the label value for the right baseline point
	/// @param bandwavelength is the label value for the desired baseline-corrected band
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size.
	bool height(double lb, double rb, double bandwavelength, T* result){

		T* lp;
		T* rp;
		unsigned long long XY = X() * Y();
		unsigned long long S = XY * sizeof(T);
		lp = (T*) malloc(S);			//memory allocation
		rp = (T*) malloc(S);

		band(lp, lb);
		band(rp, rb);

		baseline_band(lb, rb, lp, rp, bandwavelength, result);

		free(lp);
		free(rp);
		return true;
	}



	/// Calculate the area under the spectrum between two specified points and stores the result in a pre-allocated array.

	/// @param lb is the label value for the left baseline point
	/// @param rb is the label value for the right baseline point
	/// @param lab is the label value for the left bound (start of the integration)
	/// @param rab is the label value for the right bound (end of the integration)
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
	bool area(double lb, double rb, double lab, double rab, T* result){

		T* lp;	//left band pointer
		T* rp;	//right band pointer
		T* cur;		//current band 1
		T* cur2;	//current band 2

		unsigned long long XY = X() * Y();
		unsigned long long S = XY * sizeof(T);

		lp = (T*) malloc(S);			//memory allocation
		rp = (T*) malloc(S);
		cur = (T*) malloc(S);
		cur2 = (T*) malloc(S);

		memset(result, (char)0, S);

		//find the wavelenght position in the whole band
		unsigned long long n = w.size();
		unsigned long long ai = 0;		//left bound position
		unsigned long long bi = n - 1;		//right bound position



		//to make sure the left and the right bound are in the bandwidth
		if (lb < w[0] || rb < w[0] || lb > w[n-1] || rb >w[n-1]){
			std::cout<<"ERROR: left bound or right bound out of bandwidth"<<std::endl;
			exit(1);
		}
		//to make sure rigth bound is bigger than left bound
		else if(lb > rb){
			std::cout<<"ERROR: right bound should be bigger than left bound"<<std::endl;
			exit(1);
		}

		//get the position of lb and rb
		while (lab >= w[ai]){
			ai++;
		}
		while (rab <= w[bi]){
			bi--;
		}

		band(lp, lb);
		band(rp, rb);

		//calculate the beginning and the ending part
		baseline_band(lb, rb, lp, rp, rab, cur2);		//ending part
		baseline_band(lb, rb, lp, rp, w[bi], cur);
		for(unsigned long long j = 0; j < XY; j++){
			result[j] += (T)((rab - w[bi]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
		}
		baseline_band(lb, rb, lp, rp, lab, cur2);		//beginnning part
		baseline_band(lb, rb, lp, rp, w[ai], cur);
		for(unsigned long long j = 0; j < XY; j++){
			result[j] += (T)((w[ai] - lab) * ((double)cur[j] + (double)cur2[j]) / 2.0);
		}

		//calculate the area
		ai++;
		for(unsigned long long i = ai; i <= bi ;i++)
		{
			baseline_band(lb, rb, lp, rp, w[ai], cur2);
			for(unsigned long long j = 0; j < XY; j++)
			{
				result[j] += (T)((w[ai] - w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
			}
			std::swap(cur,cur2);		//swap the band pointers
		}

		free(lp);
		free(rp);
		free(cur);
		free(cur2);
		return true;
	}

	/// Compute the ratio of two baseline-corrected peaks. The result is stored in a pre-allocated array.

	/// @param lb1 is the label value for the left baseline point for the first peak (numerator)
	/// @param rb1 is the label value for the right baseline point for the first peak (numerator)
	/// @param pos1 is the label value for the first peak (numerator) position
	/// @param lb2 is the label value for the left baseline point for the second peak (denominator)
	/// @param rb2 is the label value for the right baseline point for the second peak (denominator)
	/// @param pos2 is the label value for the second peak (denominator) position
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
	bool ph_to_ph(T* result, double lb1, double rb1, double pos1, double lb2, double rb2, double pos2, unsigned char* mask = NULL){

		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
		T* p2 = (T*)malloc(X() * Y() * sizeof(T));

		//get the two peak band
		height(lb1, rb1, pos1, p1);
		height(lb2, rb2, pos2, p2);
		//calculate the ratio in result
		for(unsigned long long i = 0; i < X() * Y(); i++){
			if(p1[i] == 0 && p2[i] ==0)
				result[i] = 1;
			else
				result[i] = p1[i] / p2[i];
		}

		free(p1);
		free(p2);
		return true;
	}

	/// Compute the ratio between a peak area and peak height.

	/// @param lb1 is the label value for the left baseline point for the first peak (numerator)
	/// @param rb1 is the label value for the right baseline point for the first peak (numerator)
	/// @param pos1 is the label value for the first peak (numerator) position
	/// @param lb2 is the label value for the left baseline point for the second peak (denominator)
	/// @param rb2 is the label value for the right baseline point for the second peak (denominator)
	/// @param pos2 is the label value for the second peak (denominator) position
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
	bool pa_to_ph(T* result, double lb1, double rb1, double lab1, double rab1,
					double lb2, double rb2, double pos, unsigned char* mask = NULL){

		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
		T* p2 = (T*)malloc(X() * Y() * sizeof(T));

		//get the area and the peak band
		area(lb1, rb1, lab1, rab1, p1);
		height(lb2, rb2, pos, p2);
		//calculate the ratio in result
		for(unsigned long long i = 0; i < X() * Y(); i++){
			if(p1[i] == 0 && p2[i] ==0)
				result[i] = 1;
			else
				result[i] = p1[i] / p2[i];
		}

		free(p1);
		free(p2);
		return true;
	}

	/// Compute the ratio between two peak areas.

	/// @param lb1 is the label value for the left baseline point for the first peak (numerator)
	/// @param rb1 is the label value for the right baseline point for the first peak (numerator)
	/// @param lab1 is the label value for the left bound (start of the integration) of the first peak (numerator)
	/// @param rab1 is the label value for the right bound (end of the integration) of the first peak (numerator)
	/// @param lb2 is the label value for the left baseline point for the second peak (denominator)
	/// @param rb2 is the label value for the right baseline point for the second peak (denominator)
	/// @param lab2 is the label value for the left bound (start of the integration) of the second peak (denominator)
	/// @param rab2 is the label value for the right bound (end of the integration) of the second peak (denominator)
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
	bool pa_to_pa(T* result, double lb1, double rb1, double lab1, double rab1,
					double lb2, double rb2, double lab2, double rab2, unsigned char* mask = NULL){

		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
		T* p2 = (T*)malloc(X() * Y() * sizeof(T));

		//get the area and the peak band
		area(lb1, rb1, lab1, rab1, p1);
		area(lb2, rb2, lab2, rab2, p2);
		//calculate the ratio in result
		for(unsigned long long i = 0; i < X() * Y(); i++){
			if(p1[i] == 0 && p2[i] ==0)
				result[i] = 1;
			else
				result[i] = p1[i] / p2[i];
		}

		free(p1);
		free(p2);
		return true;
	}

	/// Compute the definite integral of a baseline corrected peak.

	/// @param lb is the label value for the left baseline point
	/// @param rb is the label value for the right baseline point
	/// @param lab is the label for the start of the definite integral
	/// @param rab is the label for the end of the definite integral
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
	bool x_area(double lb, double rb, double lab, double rab, T* result){
		T* lp;	//left band pointer
		T* rp;	//right band pointer
		T* cur;		//current band 1
		T* cur2;	//current band 2

		unsigned long long XY = X() * Y();
		unsigned long long S = XY * sizeof(T);

		lp = (T*) malloc(S);			//memory allocation
		rp = (T*) malloc(S);
		cur = (T*) malloc(S);
		cur2 = (T*) malloc(S);

		memset(result, (char)0, S);

		//find the wavelenght position in the whole band
		unsigned long long n = w.size();
		unsigned long long ai = 0;		//left bound position
		unsigned long long bi = n - 1;		//right bound position

		//to make sure the left and the right bound are in the bandwidth
		if (lb < w[0] || rb < w[0] || lb > w[n-1] || rb >w[n-1]){
			std::cout<<"ERROR: left bound or right bound out of bandwidth"<<std::endl;
			exit(1);
		}
		//to make sure rigth bound is bigger than left bound
		else if(lb > rb){
			std::cout<<"ERROR: right bound should be bigger than left bound"<<std::endl;
			exit(1);
		}

		//get the position of lb and rb
		while (lab >= w[ai]){
			ai++;
		}
		while (rab <= w[bi]){
			bi--;
		}

		band(lp, lb);
		band(rp, rb);

		//calculate the beginning and the ending part
		baseline_band(lb, rb, lp, rp, rab, cur2);		//ending part
		baseline_band(lb, rb, lp, rp, w[bi], cur);
		for(unsigned long long j = 0; j < XY; j++){
			result[j] += (T)((rab - w[bi]) * (rab + w[bi]) * ((double)cur[j] + (double)cur2[j]) / 4.0);
		}
		baseline_band(lb, rb, lp, rp, lab, cur2);		//beginnning part
		baseline_band(lb, rb, lp, rp, w[ai], cur);
		for(unsigned long long j = 0; j < XY; j++){
			result[j] += (T)((w[ai] - lab) * (w[ai] + lab) * ((double)cur[j] + (double)cur2[j]) / 4.0);
		}

		//calculate f(x) times x
		ai++;
		for(unsigned long long i = ai; i <= bi ;i++)
		{
			baseline_band(lb, rb, lp, rp, w[ai], cur2);
			for(unsigned long long j = 0; j < XY; j++)
			{
				result[j] += (T)((w[ai] - w[ai-1]) * (w[ai] + w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 4.0);
			}
			std::swap(cur,cur2);		//swap the band pointers
		}

		free(lp);
		free(rp);
		free(cur);
		free(cur2);
		return true;
	}

	/// Compute the centroid of a baseline corrected peak.

	/// @param lb is the label value for the left baseline point
	/// @param rb is the label value for the right baseline point
	/// @param lab is the label for the start of the peak
	/// @param rab is the label for the end of the peak
	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
	bool centroid(T* result, double lb, double rb, double lab, double rab, unsigned char* mask = NULL){
		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
		T* p2 = (T*)malloc(X() * Y() * sizeof(T));

		//get the area and the peak band
		x_area(lb, rb, lab, rab, p1);
		area(lb, rb, lab, rab, p2);
		//calculate the ratio in result
		for(unsigned long long i = 0; i < X() * Y(); i++){
			if(mask == NULL || mask[i])
				result[i] = p1[i] / p2[i];
		}

		free(p1);
		free(p2);
		return true;
	}

	/// Create a mask based on a given band and threshold value.

	/// All pixels in the
	/// specified band greater than the threshold are true and all pixels less than the threshold are false.
	/// @param mask_band is the band used to specify the mask
	/// @param threshold is the threshold used to determine if the mask value is true or false
	/// @param p is a pointer to a pre-allocated array at least X * Y in size
	bool build_mask(unsigned char* out_mask, double mask_band, double lower, double upper, unsigned char* mask = NULL, bool PROGRESS = false){

		T* temp = (T*)malloc(X() * Y() * sizeof(T));		//allocate memory for the certain band
		band(temp, mask_band, PROGRESS);

		for (unsigned long long i = 0; i < X() * Y(); i++) {
			if(mask == NULL || mask[i] != 0){
				if(temp[i] > lower && temp[i] < upper){
					out_mask[i] = 255;
				}
				else
					out_mask[i] = 0;
			}
		}

		free(temp);
		return true;
	}

	/// Apply a mask file to the BSQ image, setting all values outside the mask to zero.

	/// @param outfile is the name of the masked output file
	/// @param p is a pointer to memory of size X * Y, where p(i) = 0 for pixels that will be set to zero.
	bool apply_mask(std::string outfile, unsigned char* p, bool PROGRESS = false){

		std::ofstream target(outfile.c_str(), std::ios::binary);

		//I THINK THIS IS WRONG
		unsigned long long XZ = X() * Z();		//calculate the number of values in a page on disk
		unsigned long long L = XZ * sizeof(T);	//calculate the size of the page (in bytes)

		T * temp = (T*)malloc(L);		//allocate memory for a temporary page

		for (unsigned long long i = 0; i < Y(); i++)			//for each value in Y() (BIP should be X)
		{
			read_plane_xz(temp, i);							//retrieve an ZX slice, stored in temp
			for ( unsigned long long j = 0; j < Z(); j++)	//for each Z() (Y)
			{
				for (unsigned long long k = 0; k < X(); k++) //for each band
				{
				if(p[i * X() + k] == 0)
					temp[j * X() + k] = 0;
				else
					continue;
				}
			}
			target.write(reinterpret_cast<const char*>(temp), L);   //write a band data into target file
			if(PROGRESS) progress = (double)(i+1) / (double)Y() * 100;
		}
		target.close();
		free(temp);
		return true;
	}

	/// Copies all spectra corresponding to nonzero values of a mask into a pre-allocated matrix of size (B x P)
	///		where P is the number of masked pixels and B is the number of bands. The allocated memory can be accessed
	///		using the following indexing: i = p*B + b
	/// @param matrix is the destination for the pixel data
	/// @param mask is the mask
	bool sift(T* matrix, unsigned char* mask = NULL, bool PROGRESS = false){
		size_t Lbytes = sizeof(T) * X();
		T* line = (T*) malloc( Lbytes );					//allocate space for a line

		file.seekg(0, std::ios::beg);	//seek to the beginning of the file

		size_t pl;
		size_t p = 0;										//create counter variables
		for(size_t y = 0; y < Y(); y++){					//for each line in the data set
			for(size_t b = 0; b < Z(); b++){				//for each band in the data set
				pl = 0;										//initialize the pixel offset for the current line to zero (0)
				file.read( (char*)line, Lbytes );			//read the current line
				for(size_t x = 0; x < X(); x++){
					if(mask == NULL || mask[y * X() + x]){					//if the current pixel is masked
						size_t i = (p + pl) * Z() + b;		//calculate the index in the sifted matrix
						matrix[i] = line[x];				//store the current value in the line at the correct matrix location
						pl++;								//increment the pixel pointer
					}
				}
				if(PROGRESS) progress = (double)( (y+1)*Z() + 1) / (double)(Y() * Z()) * 100;
			}
			p += pl;										//add the line increment to the running pixel index
		}
		return true;
	}

	/// Saves to disk only those spectra corresponding to mask values != 0
	/// Unlike the BIP and BSQ versions of this function, this version saves a different format: the BIL file is saved as a BIP
	bool sift(std::string outfile, unsigned char* p, bool PROGRESS = false){
		// Assume X() = X, Y() = Y, Z() = Z.
		std::ofstream target(outfile.c_str(), std::ios::binary);

		//for loading pages:
		unsigned long long XZ = X() * Z();		//calculate the number of values in an XZ page on disk
		unsigned long long B = Z();			//calculate the number of bands
		unsigned long long L = XZ * sizeof(T);	//calculate the size of the page (in bytes)

		//allocate temporary memory for a XZ slice
		T* slice = (T*) malloc(L);

		//allocate temporary memory for a spectrum
		T* spec = (T*) malloc(B * sizeof(T));

		//for each slice along the y axis
		for (unsigned long long y = 0; y < Y(); y++)			//Select a page by choosing Y coordinate, Y()
		{
			read_plane_xz(slice, y);							//retrieve an ZX page, store in "slice"

			//for each sample along X
			for (unsigned long long x = 0; x < X(); x++)		//Select a pixel by choosing X coordinate in the page, X()
			{
				//if the mask != 0 at that xy pixel
				if (p[y * X() + x] != 0)					//if the mask != 0 at that XY pixel
				{
					//for each band at that pixel
					for (unsigned long long b = 0; b < B; b++)		//Select a voxel by choosing Z coordinate at the pixel
					{
						spec[b] = slice[b*X() + x];		//Pass the correct spectral value from XZ page into the spectrum to be saved.
					}
					target.write((char*)spec, B * sizeof(T));		//write that spectrum to disk. Size is L2.
				}

				if(PROGRESS) progress = (double) ((y+1) * X() + x+1) / (Y() * X()) * 100;
			}
		}
		target.close();
		free(slice);
		free(spec);

		return true;
	}

	/// Calculate the mean band value (average along B) at each pixel location.

	/// @param p is a pointer to memory of size X * Y * sizeof(T) that will store the band averages.
	bool band_avg(T* p){
		unsigned long long XZ = X() * Z();
		T* temp = (T*)malloc(sizeof(T) * XZ);
		T* line = (T*)malloc(sizeof(T) * X());

		for (unsigned long long i = 0; i < Y(); i++){
			getY(temp, i);
			//initialize x-line
			for (unsigned long long j = 0; j < X(); j++){
				line[j] = 0;
			}
			unsigned long long c = 0;
			for (unsigned long long j = 0; j < Z(); j++){
				for (unsigned long long k = 0; k < X(); k++){
					line[k] += temp[c] / (T)Z();
					c++;
				}
			}
			for (unsigned long long j = 0; j < X(); j++){
				p[j + i * X()] = line[j];
			}
		}
		free(temp);
		return true;
	}

	/// Calculate the mean value for all masked (or valid) pixels in a band and returns the average spectrum

	/// @param p is a pointer to pre-allocated memory of size [B * sizeof(T)] that stores the mean spectrum
	/// @param mask is a pointer to memory of size [X * Y] that stores the mask value at each pixel location
	bool mean_spectrum(double* m, double* std, unsigned char* mask = NULL, bool PROGRESS = false){
		unsigned long long XZ = X() * Z();
		unsigned long long XY = X() * Y();
		T* temp = (T*)malloc(sizeof(T) * XZ);
		memset(m, 0, Z() * sizeof(double));							//initialize the mean to zero
		double* e_x2 = (double*)malloc(Z() * sizeof(double));		//allocate space for E[x^2]
		memset(e_x2, 0, Z() * sizeof(double));						//initialize E[x^2] to zero
		//calculate vaild number in a band
		size_t count = nnz(mask);							//count the number of pixels in the mask

		double x;											//create a register to store the pixel value
		for (unsigned long long k = 0; k < Y(); k++){
			read_plane_xz(temp, k);
			unsigned long long kx = k * X();
			for (unsigned long long i = 0; i < X(); i++){
				if (mask == NULL || mask[kx + i] != 0){
					for (unsigned long long j = 0; j < Z(); j++){
						x = temp[j * X() + i];
						m[j] += x / (double)count;
						e_x2[j] += x*x / (double)count;
					}
				}
			}
			if(PROGRESS) progress = (double)(k+1) / Y() * 100;
		}

		for(size_t i = 0; i < Z(); i++)							//calculate the standard deviation
			std[i] = sqrt(e_x2[i] - m[i] * m[i]);

		free(temp);
		return true;
	}

	int co_matrix_cublas(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
		cublasStatus_t stat;
		cublasHandle_t handle;

		progress = 0;														//initialize the progress to zero (0)
		size_t XY = X() * Y();												//calculate the number of elements in a band image
		size_t XB = X() * Z();
		size_t B = Z();														//calculate the number of spectral elements

		double* F = (double*)malloc(sizeof(double) * B * X());				//allocate space for the frame that will be pulled from the file
		double* F_dev;
		HANDLE_ERROR(cudaMalloc(&F_dev, X() * B * sizeof(double)));			//allocate space for the frame on the GPU
		double* s_dev;														//declare a device pointer that will store the spectrum on the GPU
		double* A_dev;														//declare a device pointer that will store the covariance matrix on the GPU
		double* avg_dev;													//declare a device pointer that will store the average spectrum
		HANDLE_ERROR(cudaMalloc(&s_dev, B * sizeof(double)));				//allocate space on the CUDA device for a spectrum
		HANDLE_ERROR(cudaMalloc(&A_dev, B * B * sizeof(double)));			//allocate space on the CUDA device for the covariance matrix
		HANDLE_ERROR(cudaMemset(A_dev, 0, B * B * sizeof(double)));			//initialize the covariance matrix to zero (0)
		HANDLE_ERROR(cudaMalloc(&avg_dev, XB * sizeof(double)));				//allocate space on the CUDA device for the average spectrum
		for(size_t x = 0; x < X(); x++)											//make multiple copies of the average spectrum in order to build a matrix
			HANDLE_ERROR(cudaMemcpy(&avg_dev[x * B], avg, B * sizeof(double), cudaMemcpyHostToDevice));	
		//stat = cublasSetVector((int)B, sizeof(double), avg, 1, avg_dev, 1);	//copy the average spectrum to the CUDA device

		double ger_alpha = 1.0/(double)XY;										//scale the outer product by the inverse of the number of samples (mean outer product)
		double axpy_alpha = -1;													//multiplication factor for the average spectrum (in order to perform a subtraction)

		CUBLAS_HANDLE_ERROR(stat = cublasCreate(&handle));								//create a cuBLAS instance
		if (stat != CUBLAS_STATUS_SUCCESS) return 1;									//test the cuBLAS instance to make sure it is valid

		else std::cout<<"Using cuBLAS to calculate the mean covariance matrix..."<<std::endl;
		double beta = 1.0;
		size_t x, y;
		for(y = 0; y < Y(); y++){										//for each line
			read_plane_zxd(F, y);												//read a frame from the file
			HANDLE_ERROR(cudaMemcpy(F_dev, F, XB * sizeof(double), cudaMemcpyHostToDevice));	//copy the frame to the GPU
			CUBLAS_HANDLE_ERROR(cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, (int)B, (int)X(), &axpy_alpha, avg_dev, (int)B, &beta, F_dev, (int)B, F_dev, (int)B));//subtract the mean spectrum

			for(x = 0; x < X(); x++)
				CUBLAS_HANDLE_ERROR(cublasDsyr(handle, CUBLAS_FILL_MODE_UPPER, (int)B, &ger_alpha, &F_dev[x*B], 1, A_dev, (int)B));			//perform an outer product
			if(PROGRESS) progress = (double)(y + 1) / Y() * 100;
		}

		cublasGetMatrix((int)B, (int)B, sizeof(double), A_dev, (int)B, co, (int)B);			//copy the result from the GPU to the CPU

		cudaFree(A_dev);																	//clean up allocated device memory
		cudaFree(s_dev);
		cudaFree(avg_dev);

		for(unsigned long long i = 0; i < B; i++){										//copy the upper triangular portion to the lower triangular portion
			for(unsigned long long j = i+1; j < B; j++){
				co[B * i + j] = co[B * j + i];
			}
		}

		return 0;



	}


	/// Calculate the covariance matrix for all masked pixels in the image.

	/// @param co is a pointer to pre-allocated memory of size [B * B] that stores the resulting covariance matrix
	/// @param avg is a pointer to memory of size B that stores the average spectrum
	/// @param mask is a pointer to memory of size [X * Y] that stores the mask value at each pixel location
	bool co_matrix(double* co, double* avg, unsigned char *mask, int cuda_device = 0, bool PROGRESS = false){
		progress = 0;

		if (cuda_device >= 0) {													//if a CUDA device is specified
			int dev_count;
			HANDLE_ERROR(cudaGetDeviceCount(&dev_count));						//get the number of CUDA devices
			std::cout << "Number of CUDA devices: " << dev_count << std::endl;		//output the number of CUDA devices
			cudaDeviceProp prop;
			//int best_device_id = 0;													//stores the best CUDA device
			//float best_device_cc = 0.0f;												//stores the compute capability of the best device
			std::cout << "CUDA devices----" << std::endl;
			for (int d = 0; d < dev_count; d++) {									//for each CUDA device
				cudaGetDeviceProperties(&prop, d);								//get the property of the first device
				//float cc = prop.major + prop.minor / 10.0f;						//calculate the compute capability
				std::cout << d << ":  [" << prop.major << "." << prop.minor << "]      " << prop.name << std::endl;	//display the device information
				//if(cc > best_device_cc){
				//	best_device_cc = cc;										//if this is better than the previous device, use it
				//	best_device_id = d;
				//}
			}
			if (dev_count > 0 && dev_count > cuda_device) {							//if the first device is not an emulator
				cudaGetDeviceProperties(&prop, cuda_device);									//get the property of the requested CUDA device
				if (prop.major != 9999) {
					std::cout << "Using device " << cuda_device << std::endl;
					HANDLE_ERROR(cudaSetDevice(cuda_device));
					int status = co_matrix_cublas(co, avg, mask, PROGRESS);			//use cuBLAS to calculate the covariance matrix
					if (status == 0) return true;									//if the cuBLAS function returned correctly, we're done
				}
			}
		}

		//memory allocation
		unsigned long long xy = X() * Y();
		unsigned long long B = Z();
		T* temp = (T*)malloc(sizeof(T) * B);
		//count vaild pixels in a band
		unsigned long long count = 0;
		for (unsigned long long j = 0; j < xy; j++){
			if (mask == NULL || mask[j] != 0){
				count++;
			}
		}
		//initialize correlation matrix
		for (unsigned long long i = 0; i < B; i++){
			for (unsigned long long k = 0; k < B; k++){
				co[i * B + k] = 0;
			}
		}
		//calculate correlation coefficient matrix
		for (unsigned long long j = 0; j < xy; j++){
			if (mask == NULL || mask[j] != 0){
				pixel(temp, j);
				for (unsigned long long i = 0; i < B; i++){
					for (unsigned long long k = i; k < B; k++){
						co[i * B + k] += ((double)temp[i] - (double)avg[i]) * ((double)temp[k] - (double)avg[k]) / (double)count;
					}
				}
			}
			if(PROGRESS) progress = (double)(j+1) / xy * 100;
		}
		//because correlation matrix is symmetric
		for (unsigned long long i = 0; i < B; i++){
			for (unsigned long long k = i + 1; k < B; k++){
				co[k * B + i] = co[i * B + k];
			}
		}

		free(temp);
		return true;
	}


	/// Crop a region of the image and save it to a new file.

	/// @param outfile is the file name for the new cropped image
	/// @param x0 is the lower-left x pixel coordinate to be included in the cropped image
	/// @param y0 is the lower-left y pixel coordinate to be included in the cropped image
	/// @param x1 is the upper-right x pixel coordinate to be included in the cropped image
	/// @param y1 is the upper-right y pixel coordinate to be included in the cropped image
	bool crop(std::string outfile, unsigned long long x0,
								   unsigned long long y0,
								   unsigned long long x1,
								   unsigned long long y1,
								   unsigned long long b0,
								   unsigned long long b1,
								   bool PROGRESS = false){

		//calculate the new image parameters
		unsigned long long samples = x1 - x0 + 1;
		unsigned long long lines = y1 - y0 + 1;
		unsigned long long bands = b1 - b0 + 1;

		//calculate the size of a line
		unsigned long long L = samples * sizeof(T);


		//allocate space for a line
		T* temp = (T*)malloc(bands * L);

		//create an output stream to store the cropped file
		std::ofstream out(outfile.c_str(), std::ios::binary);

		//calculate the distance between bands
		unsigned long long jumpb = (X() - samples) * sizeof(T);

		//distance needed to jump from the previous line of the last band to the next line of the first band
		unsigned long long longjump = ((Z() - bands) * X()) * sizeof(T);

		//set the start position for the cropped region
		file.seekg((y0 * X() * Z() + b0 * X() + x0) * sizeof(T), std::ios::beg);

		for (unsigned long long x = 0; x < lines; x++)
		{
			for (unsigned long long z = b0; z <= b1; z++)
			{
				file.read((char *)(temp + z * samples), sizeof(T) * samples);
				file.seekg(jumpb, std::ios::cur);    //go to the next band

				if(PROGRESS) progress = (double)(x * (b1 - b0 + 1) + z + 1) / (lines * (b1 - b0 + 1)) * 100;
			}

			//write slice data into target file
			out.write(reinterpret_cast<const char*>(temp), bands * L);

			//seek to the beginning of the next X-Z slice
			file.seekg(longjump, std::ios::cur);
		}

		//free the temporary frame
		free(temp);

		return true;
	}

	/// Remove a list of bands from the ENVI file

	/// @param outfile is the file name for the output hyperspectral image (with trimmed bands)
	/// @param b is an array of bands to be eliminated
	void trim(std::string outfile, std::vector<size_t> band_array, bool PROGRESS = false){
		std::ofstream out(outfile.c_str(), std::ios::binary);	//open the output file for writing
		file.seekg(0, std::ios::beg);							//move to the beginning of the input file

		size_t Xb = X() * sizeof(T);							//calculate the number of bytes in a line
		T* line = (T*)malloc(Xb);								//allocate space for a line

		size_t i;												//create an index into the band array
		for(size_t y = 0; y < Y(); y++){						//for each Y plane
			i = 0;
			for(size_t b = 0; b < Z(); b++){					//for each band
				if(b != band_array[i]){							//if this band isn't trimmed
					file.read((char*)line, Xb);					//read a line
					out.write((char*)line, Xb);					//write the line
				}
				else{
					file.seekg(Xb, std::ios::cur);				//if this band is trimmed, skip it
					i++;
				}
			}
			if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100;
		}
		free(line);
	}

	/// Combine two BSQ images along the Y axis

	/// @param outfile is the combined file to be output
	/// @param infile is the input file stream for the image to combine with this one
	/// @param Sx is the size of the second image along X
	/// @param Sy is the size of the second image along Y
	/// @param offset is a shift (negative or positive) in the combined image to the left or right
	void combine(std::string outfile, bil<T>* C, long long xp, long long yp, bool PROGRESS = false){
		std::ofstream out(outfile.c_str(), std::ios::binary);	//open the output file for writing
		file.seekg(0, std::ios::beg);								//move to the beginning of both files
		C->file.seekg(0, std::ios::beg);

		size_t S[2];				//size of the output band image
		size_t p0[2];				//position of the current image in the output
		size_t p1[2];				//position of the source image in the output

		hsi<T>::calc_combined_size(xp, yp, C->X(), C->Y(), S[0], S[1], p0[0], p0[1], p1[0], p1[1]);	//calculate the image placement parameters

		size_t line_bytes = X() * sizeof(T);
		size_t band_bytes = X() * Y() * sizeof(T);
		T* cur = (T*)malloc(line_bytes);					//allocate space for a band of the current image

		size_t line_src_bytes = C->X() * sizeof(T);
		size_t band_src_bytes = C->X() * C->Y() * sizeof(T);
		T* src = (T*)malloc(line_src_bytes);			//allocate space for a band of the source image

		size_t line_dst_bytes = S[0] * sizeof(T);
		size_t band_dst_bytes = S[0] * S[1] * sizeof(T);
		T* dst = (T*)malloc(line_dst_bytes);						//allocate space for a band of the destination image

		for(size_t y = 0; y < S[1]; y++){							//for each line in the destination file
			memset(dst, 0, line_dst_bytes);							//set all values to zero (0) in the destination image
			for(size_t b = 0; b < Z(); b++){						//for each band in both images
				if(y >= p0[1] && y < p0[1] + Y()){					//if the current image crosses this line
					file.read((char*)cur, line_bytes);				//read a line from the current image
					memcpy(&dst[p0[0]], cur, line_bytes);			//copy the current line to the correct spot in the destination line
				}
				if(y >= p1[1] && y < p1[1] + C->Y()){				//if the source image crosses this line
					C->file.read((char*)src, line_src_bytes);		//read a line from the source image
					memcpy(&dst[p1[0]], src, line_src_bytes);		//copy the source line into the correct spot in the destination line
				}
				out.write((char*)dst, line_dst_bytes);
				if(PROGRESS) progress = (double)((y + 1)*Z() + b + 1) / (double) (Z() * S[1]) * 100;
			}
		}
	}

	///Append two files together along the band dimension
	void append(std::string outfile, bil<T>* C, bool PROGRESS = false) {
		std::ofstream out(outfile.c_str(), std::ios::binary);	//open the output file for writing
		file.seekg(0, std::ios::beg);							//move to the beginning of both files
		C->file.seekg(0, std::ios::beg);
		size_t a_bytes = X() * Z() * sizeof(T);					//calculate the number of bytes in a single plane of this file
		size_t b_bytes = C->X() * C->Z() * sizeof(T);			//calculate the number of bytes in a single plane of the appending file
		T* a = (T*)malloc(a_bytes);								//allocate space for a plane of the current file
		T* b = (T*)malloc(b_bytes);								//allocate space for a plane of the appended file
		if (PROGRESS) progress = 0;
		for (size_t y = 0; y < Y(); y++) {
			read_plane_xz(a, y);								//read a plane from the current file
			out.write((char*)a, a_bytes);								//write the plane to disk
			C->read_plane_xz(b, y);								//read a plane from the appending file
			out.write((char*)b, b_bytes);
			if (PROGRESS) progress = (double)(y + 1) / (double)(Y()) * 100;
		}
		
		out.close();
	}

	/// Convolve the given band range with a kernel specified by a vector of coefficients.

	/// @param outfile is an already open stream to the output file
	/// @param C is an array of coefficients
	/// @param start is the band to start processing (the first coefficient starts here)
	/// @param nbands is the number of bands to process
	/// @param center is the index for the center coefficient for the kernel (used to set the wavelengths in the output file)

	void convolve(std::string outfile, std::vector<double> C, size_t start, size_t end, unsigned char* mask = NULL, bool PROGRESS = false){
		std::ofstream out(outfile.c_str(), std::ios::binary);	//open the output file for writing
		size_t B = end - start + 1;
		size_t Xb = X() * sizeof(T);								//calculate the size of a band (frame) in bytes
		size_t XBb = Xb * Z();

		size_t N = C.size();										//get the number of bands that the kernel spans
		std::deque<T*> frame(N, NULL);								//create an array to store pointers to each frame
		for(size_t f = 0; f < N; f++)
			frame[f] = (T*)malloc(Xb);								//allocate space for the frame

		T* outline = (T*)malloc(Xb);								//allocate space for the output band

		//Implementation: In order to minimize reads from secondary storage, each band is only loaded once into the 'frame' deque.
		//					When a new band is loaded, the last band is popped, a new frame is copied to the pointer, and it is
		//					re-inserted into the deque.
		for(size_t y = 0; y < Y(); y++){
			file.seekg(y * XBb + start * Xb, std::ios::beg);		//move to the beginning of the 'start' band
			for(size_t f = 0; f < N; f++)								//for each frame
				file.read((char*)frame[f], Xb);							//load the frame
			for(size_t b = 0; b < B; b++){								//for each band
				memset(outline, 0, Xb);									//set the output band to zero (0)
				for(size_t c = 0; c < N; c++){							//for each frame (corresponding to each coefficient)
					for(size_t x = 0; x < X(); x++){					//for each pixel
						if(mask == NULL || mask[y * X() + x]){
							outline[x] += (T)(C[c] * frame[c][x]);		//calculate the contribution of the current frame (scaled by the corresponding coefficient)
						}
					}
				}
				out.write((char*)outline, Xb);							//output the band
				file.read((char*)frame[0], Xb);							//read the next band
				frame.push_back(frame.front());							//put the first element in the back
				frame.pop_front();										//pop the first element
			}
			if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100;
		}
	}

	/// Approximate the spectral derivative of the image
	void deriv(std::string outfile, size_t d, size_t order, const std::vector<double> w = std::vector<double>(), unsigned char* mask = NULL, bool PROGRESS = false){
		std::ofstream out(outfile.c_str(), std::ios::binary);		//open the output file for writing

		size_t Xb = X() * sizeof(T);								//calculate the size of a line (frame) in bytes
		size_t XBb = Xb * Z();
		size_t B = Z();

		file.seekg(0, std::ios::beg);								//seek to the beginning of the file

		size_t N = order + d;										//approximating a derivative requires order + d samples
		std::deque<T*> frame(N, NULL);								//create an array to store pointers to each frame
		for(size_t f = 0; f < N; f++)								//for each frame
			frame[f] = (T*)malloc(Xb);								//allocate space for the frame

		T* outline = (T*)malloc(Xb);								//allocate space for the output band

		//Implementation: In order to minimize reads from secondary storage, each band is only loaded once into the 'frame' deque.
		//					When a new band is loaded, the last band is popped, a new frame is copied to the pointer, and it is
		//					re-inserted into the deque.
		size_t mid = (size_t)(N / 2);								//calculate the mid point of the kernel
		size_t iw;													//index to the first wavelength used to evaluate the derivative at this band

		for(size_t y = 0; y < Y(); y++){
			//file.seekg(y * XBb, std::ios::beg);							//seek to the beginning of the current Y plane
			for(size_t f = 0; f < N; f++)								//for each frame
				file.read((char*)frame[f], Xb);							//load a line

			for(size_t b = 0; b < B; b++){								//for each band
				if(b < mid)												//calculate the first wavelength used to evaluate the derivative at this band
					iw = 0;
				else if(b > B - (N - mid + 1))
					iw = B - N;
				else{
					iw = b - mid;
					file.read((char*)frame[0], Xb);											//read the next band
					frame.push_back(frame.front());											//put the first element in the back
					frame.pop_front();														//pop the first element
				}
				std::vector<double> w_pts(w.begin() + iw, w.begin() + iw + N);			//get the wavelengths corresponding to each sample
				std::vector<double> C = diff_coefficients(w[b], w_pts, d);				//get the optimal sample weights

				memset(outline, 0, Xb);													//set the output band to zero (0)
				for(size_t c = 0; c < N; c++){											//for each frame (corresponding to each coefficient)
					for(size_t x = 0; x < X(); x++){									//for each pixel
						if(mask == NULL || mask[y * X() + x]){
							outline[x] += (T)(C[c] * frame[c][x]);			//calculate the contribution of the current frame (scaled by the corresponding coefficient)
						}
					}
				}
				out.write((char*)outline, Xb);											//output the band

			}
			if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100;
		}

	}

	bool multiply(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){
		unsigned long long B = Z();									//calculate the number of bands
		unsigned long long ZX = Z() * X();
		unsigned long long XY = X() * Y();							//calculate the number of pixels in a band
		unsigned long long S = XY * sizeof(T);						//calculate the number of bytes in a band
		unsigned long long L = ZX * sizeof(T);

		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file

		T * c;														//pointer to the current ZX slice
		c = (T*)malloc( L );										//allocate space for the slice

		for(unsigned long long j = 0; j < Y(); j++){				//for each line
			read_plane_xz(c, j);										//load the line into memory
			for(unsigned long long i = 0; i < B; i++){				//for each band
				for(unsigned long long m = 0; m < X(); m++){		//for each sample
					if( mask == NULL && mask[m + j * X()] )			//if the pixel is masked
						c[m + i * X()] *= (T)v;
				}
			}
			target.write(reinterpret_cast<const char*>(c), L);		//write normalized data into destination

			if(PROGRESS) progress = (double)(j+1) / Y() * 100;		//update the progress
		}

		free(c);													//free the slice memory
		target.close();												//close the output file
		return true;
	}

	bool add(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){
		unsigned long long B = Z();									//calculate the number of bands
		unsigned long long ZX = Z() * X();
		unsigned long long XY = X() * Y();							//calculate the number of pixels in a band
		unsigned long long S = XY * sizeof(T);						//calculate the number of bytes in a band
		unsigned long long L = ZX * sizeof(T);

		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file

		T * c;														//pointer to the current ZX slice
		c = (T*)malloc( L );										//allocate space for the slice

		for(unsigned long long j = 0; j < Y(); j++){				//for each line
			read_plane_xz(c, j);										//load the line into memory
			for(unsigned long long i = 0; i < B; i++){				//for each band
				for(unsigned long long m = 0; m < X(); m++){		//for each sample
					if( mask == NULL && mask[m + j * X()] )			//if the pixel is masked
						c[m + i * X()] += (T)v;
				}
			}
			target.write(reinterpret_cast<const char*>(c), L);		//write normalized data into destination

			if(PROGRESS) progress = (double)(j+1) / Y() * 100;		//update the progress
		}

		free(c);													//free the slice memory
		target.close();												//close the output file
		return true;
	}

	/// Close the file.
	bool close(){
		file.close();
		return true;
	}

	};
}

#endif