Blame view

stim/envi/bip.h 80.8 KB
e8eb202f   David Mayerich   added a new ENVI ...
1
2
3
4
  #ifndef STIM_BIP_H
  #define STIM_BIP_H
  
  #include "../envi/envi_header.h"
6aa04ba2   David Mayerich   interleave types
5
  #include "../envi/bil.h"
9d3ba0b1   David Mayerich   added stim::hsi a...
6
  #include "../envi/hsi.h"
88c3e636   heziqi   Ziqi added big.h
7
  #include <cstring>
9d34229e   David Mayerich   restructured the ...
8
  #include <complex>
88c3e636   heziqi   Ziqi added big.h
9
10
  #include <utility>
  
c6f526c2   David Mayerich   added cuBLAS supp...
11
  //CUDA
2781a48c   David Mayerich   fixed cuBLAS func...
12
13
14
15
  //#ifdef CUDA_FOUND
  #include <stim/cuda/cudatools/error.h>
  #include <cuda_runtime.h>
  #include "cublas_v2.h"
9d34229e   David Mayerich   restructured the ...
16
  #include "cufft.h"
2781a48c   David Mayerich   fixed cuBLAS func...
17
  //#endif
c6f526c2   David Mayerich   added cuBLAS supp...
18
  
8a86bd56   David Mayerich   changed rts names...
19
  namespace stim{
88c3e636   heziqi   Ziqi added big.h
20
  
a23c4132   David Mayerich   Doxygen comments ...
21
22
23
24
25
26
27
28
  /**
  	The BIP class represents a 3-dimensional binary file stored using band interleaved by pixel (BIP) image encoding. The binary file is stored
  	such that Z-X "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 Z, X, 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.
  */
88c3e636   heziqi   Ziqi added big.h
29
30
  template <typename T>
  
9d3ba0b1   David Mayerich   added stim::hsi a...
31
  class bip: public hsi<T> {
88c3e636   heziqi   Ziqi added big.h
32
33
  
  protected:
1a224b6a   David Mayerich   fixed linux compa...
34
35
  
  
9d3ba0b1   David Mayerich   added stim::hsi a...
36
37
  	//std::vector<double> w; //band wavelength
  	unsigned long long offset;		//header offset
88c3e636   heziqi   Ziqi added big.h
38
  
9d3ba0b1   David Mayerich   added stim::hsi a...
39
40
  	using hsi<T>::w;				//use the wavelength array in stim::hsi
  	using hsi<T>::nnz;
63fc1df8   David Mayerich   added an inverse ...
41
  	using binary<T>::progress;
1a224b6a   David Mayerich   fixed linux compa...
42
43
44
45
  	using hsi<T>::X;
  	using hsi<T>::Y;
  	using hsi<T>::Z;
  
88c3e636   heziqi   Ziqi added big.h
46
47
48
49
  public:
  
  	using binary<T>::open;
  	using binary<T>::file;
6708cc25   heziqi   Ziqi added envi c...
50
  	using binary<T>::R;
63fc1df8   David Mayerich   added an inverse ...
51
  	using binary<T>::read_line_0;
88c3e636   heziqi   Ziqi added big.h
52
  
9d3ba0b1   David Mayerich   added stim::hsi a...
53
54
  	bip(){ hsi<T>::init_bip(); }
  
a23c4132   David Mayerich   Doxygen comments ...
55
56
57
58
59
60
61
62
  	/// 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
1a224b6a   David Mayerich   fixed linux compa...
63
64
65
66
67
  	bool open(std::string filename,
  			  unsigned long long X,
  			  unsigned long long Y,
  			  unsigned long long B,
  			  unsigned long long header_offset,
bb457896   David Mayerich   improved error me...
68
69
  			  std::vector<double> wavelengths,
  			  stim::iotype io = stim::io_in){
88c3e636   heziqi   Ziqi added big.h
70
  
6708cc25   heziqi   Ziqi added envi c...
71
72
73
74
  		//copy the wavelengths to the BSQ file structure
  		w = wavelengths;
  		//copy the offset to the structure
  		offset = header_offset;
88c3e636   heziqi   Ziqi added big.h
75
  
bb457896   David Mayerich   improved error me...
76
  		return open(filename, vec<unsigned long long>(B, X, Y), header_offset, io);
1a224b6a   David Mayerich   fixed linux compa...
77
  
88c3e636   heziqi   Ziqi added big.h
78
79
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
80
81
82
83
  	/// 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.
9d3ba0b1   David Mayerich   added stim::hsi a...
84
  	bool band_index( T * p, unsigned long long page, bool PROGRESS = false){
63fc1df8   David Mayerich   added an inverse ...
85
  		return binary<T>::read_plane_0(p, page, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
86
87
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
88
89
90
91
  	/// 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.
63fc1df8   David Mayerich   added an inverse ...
92
  	bool band( T * p, double wavelength, bool PROGRESS = false){
88c3e636   heziqi   Ziqi added big.h
93
  
ee4dea28   David Mayerich   fixed errors in c...
94
95
  		//if there are no wavelengths in the BSQ file
  		if(w.size() == 0)
9d3ba0b1   David Mayerich   added stim::hsi a...
96
  			return band_index(p, (unsigned long long)wavelength, PROGRESS);
ee4dea28   David Mayerich   fixed errors in c...
97
  
9d3ba0b1   David Mayerich   added stim::hsi a...
98
  		unsigned long long XY = X() * Y();	//calculate the number of pixels in a band
88c3e636   heziqi   Ziqi added big.h
99
100
  
  		unsigned page=0;                      //bands around the wavelength
517876d6   heziqi   metrics finished ...
101
  
88c3e636   heziqi   Ziqi added big.h
102
103
104
105
  
  		//get the bands numbers around the wavelength
  
  		//if wavelength is smaller than the first one in header file
6708cc25   heziqi   Ziqi added envi c...
106
  		if ( w[page] > wavelength ){
63fc1df8   David Mayerich   added an inverse ...
107
  			band_index(p, page, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
108
109
110
  			return true;
  		}
  
6708cc25   heziqi   Ziqi added envi c...
111
  		while( w[page] < wavelength )
88c3e636   heziqi   Ziqi added big.h
112
113
114
  		{
  			page++;
  			//if wavelength is larger than the last wavelength in header file
724ec347   David Mayerich   simplified the EN...
115
  			if (page == Z()) {
63fc1df8   David Mayerich   added an inverse ...
116
  				band_index(p, Z()-1, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
117
118
119
  				return true;
  			}
  		}
6708cc25   heziqi   Ziqi added envi c...
120
  		if ( wavelength < w[page] )
88c3e636   heziqi   Ziqi added big.h
121
  		{
517876d6   heziqi   metrics finished ...
122
123
  			T * p1;
  			T * p2;
88c3e636   heziqi   Ziqi added big.h
124
125
126
  			p1=(T*)malloc( XY * sizeof(T));                     //memory allocation
  			p2=(T*)malloc( XY * sizeof(T));
  			band_index(p1, page - 1);
63fc1df8   David Mayerich   added an inverse ...
127
  			band_index(p2, page, PROGRESS);
9d3ba0b1   David Mayerich   added stim::hsi a...
128
  			for(unsigned long long i=0; i < XY; i++){
6708cc25   heziqi   Ziqi added envi c...
129
  				double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
9d3ba0b1   David Mayerich   added stim::hsi a...
130
  				p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]);
88c3e636   heziqi   Ziqi added big.h
131
  			}
517876d6   heziqi   metrics finished ...
132
133
  			free(p1);
  			free(p2);
88c3e636   heziqi   Ziqi added big.h
134
135
136
  		}
  		else                           //if the wavelength is equal to a wavelength in header file
  		{
63fc1df8   David Mayerich   added an inverse ...
137
  			band_index(p, page, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
138
  		}
88c3e636   heziqi   Ziqi added big.h
139
140
  		return true;
  	}
724ec347   David Mayerich   simplified the EN...
141
142
143
144
145
146
  
  	/// Retrieve a single spectrum (Z-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.
9d3ba0b1   David Mayerich   added stim::hsi a...
147
  	bool spectrum(T * p, unsigned long long x, unsigned long long y, bool PROGRESS = false){
63fc1df8   David Mayerich   added an inverse ...
148
  		return read_line_0(p, x, y, PROGRESS);				//read a line in the binary YZ plane (dimension order for BIP is ZXY)
724ec347   David Mayerich   simplified the EN...
149
  	}
735a2a24   David Mayerich   started testing o...
150
151
152
153
154
  	bool spectrum(T* p, size_t n, bool PROGRESS = false){
  		size_t y = n / X();
  		size_t x = n - y * X();
  		return read_line_0(p, x, y, PROGRESS);				//read a line in the binary YZ plane (dimension order for BIP is ZXY)
  	}
724ec347   David Mayerich   simplified the EN...
155
156
157
158
159
160
161
  
  	/// Retrieves a band of x values from a given xz plane.
  
  	/// @param p is a pointer to pre-allocated memory at least X * sizeof(T) in size
  	/// @param c is a pointer to an existing XZ plane (size X*Z*sizeof(T))
  	/// @param wavelength is the wavelength of X values to retrieve
  	bool read_x_from_xz(T* p, T* c, double wavelength)
11f177d5   heziqi   Ziqi added functi...
162
  	{
9d3ba0b1   David Mayerich   added stim::hsi a...
163
  		unsigned long long B = Z();
11f177d5   heziqi   Ziqi added functi...
164
  
9d3ba0b1   David Mayerich   added stim::hsi a...
165
  		unsigned long long page=0;                      //samples around the wavelength
517876d6   heziqi   metrics finished ...
166
  
11f177d5   heziqi   Ziqi added functi...
167
168
169
170
  
  		//get the bands numbers around the wavelength
  
  		//if wavelength is smaller than the first one in header file
6708cc25   heziqi   Ziqi added envi c...
171
  		if ( w[page] > wavelength ){
9d3ba0b1   David Mayerich   added stim::hsi a...
172
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
173
174
  			{
  				p[j] = c[j * B];
1a224b6a   David Mayerich   fixed linux compa...
175
  			}
11f177d5   heziqi   Ziqi added functi...
176
177
178
  			return true;
  		}
  
6708cc25   heziqi   Ziqi added envi c...
179
  		while( w[page] < wavelength )
11f177d5   heziqi   Ziqi added functi...
180
181
182
183
  		{
  			page++;
  			//if wavelength is larger than the last wavelength in header file
  			if (page == B) {
9d3ba0b1   David Mayerich   added stim::hsi a...
184
  				for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
185
186
187
188
189
190
  				{
  					p[j] = c[(j + 1) * B - 1];
  				}
  				return true;
  			}
  		}
6708cc25   heziqi   Ziqi added envi c...
191
  		if ( wavelength < w[page] )
11f177d5   heziqi   Ziqi added functi...
192
  		{
517876d6   heziqi   metrics finished ...
193
194
  			T * p1;
  			T * p2;
724ec347   David Mayerich   simplified the EN...
195
196
  			p1=(T*)malloc( X() * sizeof(T));                     //memory allocation
  			p2=(T*)malloc( X() * sizeof(T));
11f177d5   heziqi   Ziqi added functi...
197
  			//band_index(p1, page - 1);
9d3ba0b1   David Mayerich   added stim::hsi a...
198
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
199
200
201
202
  			{
  				p1[j] = c[j * B + page - 1];
  			}
  			//band_index(p2, page );
9d3ba0b1   David Mayerich   added stim::hsi a...
203
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
204
205
206
  			{
  				p2[j] = c[j * B + page];
  			}
1a224b6a   David Mayerich   fixed linux compa...
207
  
9d3ba0b1   David Mayerich   added stim::hsi a...
208
  			for(unsigned long long i=0; i < X(); i++){
6708cc25   heziqi   Ziqi added envi c...
209
  				double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
11f177d5   heziqi   Ziqi added functi...
210
211
  				p[i] = (p2[i] - p1[i]) * r + p1[i];
  			}
517876d6   heziqi   metrics finished ...
212
213
  			free(p1);
  			free(p2);
11f177d5   heziqi   Ziqi added functi...
214
215
216
217
  		}
  		else                           //if the wavelength is equal to a wavelength in header file
  		{
  			//band_index(p, page);
9d3ba0b1   David Mayerich   added stim::hsi a...
218
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
219
220
221
222
223
  			{
  				p[j] = c[j * B + page];
  			}
  		}
  
1a224b6a   David Mayerich   fixed linux compa...
224
  		return true;
11f177d5   heziqi   Ziqi added functi...
225
  	}
bfe9c6db   heziqi   added feature_mat...
226
  
c6f526c2   David Mayerich   added cuBLAS supp...
227
  	/// Retrieve a single pixel and store it in a pre-allocated double array.
9d3ba0b1   David Mayerich   added stim::hsi a...
228
229
  	bool pixeld(double* p, unsigned long long n){
  		unsigned long long bandnum = X() * Y();		//calculate numbers in one band
c6f526c2   David Mayerich   added cuBLAS supp...
230
231
232
233
  		if ( n >= bandnum){							//make sure the pixel number is right
  			std::cout<<"ERROR: sample or line out of range"<<std::endl;
  			return false;
  		}
9d3ba0b1   David Mayerich   added stim::hsi a...
234
  		unsigned long long B = Z();
c6f526c2   David Mayerich   added cuBLAS supp...
235
236
237
238
239
  
  		T* temp = (T*) malloc(B * sizeof(T));						//allocate space for the raw pixel data
  		file.seekg(n * B * sizeof(T), std::ios::beg);				//point to the certain pixel
  		file.read((char *)temp, sizeof(T) * B);					//read the spectrum from disk to the temp pointer
  
9d3ba0b1   David Mayerich   added stim::hsi a...
240
  		for(unsigned long long i = 0; i < B; i++)						//for each element of the spectrum
c6f526c2   David Mayerich   added cuBLAS supp...
241
  			p[i] = (double) temp[i];							//cast each element to a double value
193bb525   David Mayerich   fixed memory leak...
242
  		free(temp);												//free temporary memory
9d3ba0b1   David Mayerich   added stim::hsi a...
243
  		return true;
c6f526c2   David Mayerich   added cuBLAS supp...
244
245
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
246
247
248
249
  	/// 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.
9d3ba0b1   David Mayerich   added stim::hsi a...
250
  	bool pixel(T * p, unsigned long long n){
bfe9c6db   heziqi   added feature_mat...
251
  
9d3ba0b1   David Mayerich   added stim::hsi a...
252
  		unsigned long long N = X() * Y();					//calculate numbers in one band
63fc1df8   David Mayerich   added an inverse ...
253
  		if ( n >= N){							//make sure the pixel number is right
bfe9c6db   heziqi   added feature_mat...
254
255
256
257
  			std::cout<<"ERROR: sample or line out of range"<<std::endl;
  			return false;
  		}
  
724ec347   David Mayerich   simplified the EN...
258
259
260
  		file.seekg(n * Z() * sizeof(T), std::ios::beg);           //point to the certain pixel
  		file.read((char *)p, sizeof(T) * Z());
  		return true;
bfe9c6db   heziqi   added feature_mat...
261
  	}
1a224b6a   David Mayerich   fixed linux compa...
262
  
11f177d5   heziqi   Ziqi added functi...
263
  	//given a Y ,return a ZX slice
2781a48c   David Mayerich   fixed cuBLAS func...
264
  	bool read_plane_y(T * p, size_t y){
724ec347   David Mayerich   simplified the EN...
265
  		return binary<T>::read_plane_2(p, y);
11f177d5   heziqi   Ziqi added functi...
266
  	}
1a224b6a   David Mayerich   fixed linux compa...
267
  
a23c4132   David Mayerich   Doxygen comments ...
268
269
270
  	/// 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.
1a224b6a   David Mayerich   fixed linux compa...
271
  	/// @param wls is the list of baseline points based on band labels.
798cea97   David Mayerich   added progress ba...
272
  	bool baseline(std::string outname, std::vector<double> base_pts, unsigned char* mask = NULL, bool PROGRESS = false){
1a224b6a   David Mayerich   fixed linux compa...
273
  
88c3e636   heziqi   Ziqi added big.h
274
  		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file
1a224b6a   David Mayerich   fixed linux compa...
275
  
63fc1df8   David Mayerich   added an inverse ...
276
277
278
279
280
281
282
283
  		unsigned long long N = X() * Y();						//calculate the total number of pixels to be processed
  		unsigned long long B = Z();								//get the number of bands
  		T* s = (T*)malloc(sizeof(T) * B);						//allocate memory to store a pixel
  		T* sbc = (T*)malloc(sizeof(T) * B);						//allocate memory to store the baseline corrected spectrum
  
  		std::vector<T> base_vals;								//allocate space for the values at each baseline point
  		double aw, bw;											//surrounding baseline point wavelengths
  		T av, bv;												//surrounding baseline point values
63fc1df8   David Mayerich   added an inverse ...
284
285
  		unsigned long long ai, bi;								//surrounding baseline point band indices
  		for(unsigned long long n = 0; n < N; n++){				//for each pixel in the image
798cea97   David Mayerich   added progress ba...
286
287
288
289
  			if(mask != NULL && !mask[n]){						//if the pixel isn't masked
  				memset(sbc, 0, sizeof(T) * B);					//set the baseline corrected spectrum to zero
  			}
  			else{
1a224b6a   David Mayerich   fixed linux compa...
290
  
798cea97   David Mayerich   added progress ba...
291
  				pixel(s, n);										//retrieve the spectrum s
1a224b6a   David Mayerich   fixed linux compa...
292
  				base_vals = hsi<T>::interp_spectrum(s, base_pts);			//get the values at each baseline point
798cea97   David Mayerich   added progress ba...
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
  
  				ai = bi = 0;
  				aw = w[0];											//initialize the current baseline points (assume the spectrum starts at 0)
  				av = s[0];
  				bw = base_pts[0];
  				for(unsigned long long b = 0; b < B; b++){			//for each band in the spectrum
  					while(bi < base_pts.size() && base_pts[bi] < w[b])	//while the current wavelength is higher than the second baseline point
  						bi++;											//move to the next baseline point
  					if(bi < base_pts.size()){
  						bw = base_pts[bi];							//set the wavelength for the upper bound baseline point
  						bv = base_vals[bi];							//set the value for the upper bound baseline point
  					}
  					if(bi == base_pts.size()){						//if we have passed the last baseline point
  						bw = w[B-1];								//set the outer bound to the last spectral band
  						bv = s[B-1];
  					}
  					if(bi != 0){
  						ai = bi - 1;								//set the lower bound baseline point index
  						aw = base_pts[ai];							//set the wavelength for the lower bound baseline point
  						av = base_vals[ai];							//set the value for the lower bound baseline point
  					}
1a224b6a   David Mayerich   fixed linux compa...
314
  					sbc[b] = s[b] - hsi<T>::lerp(w[b], av, aw, bv, bw);		//perform the baseline correction and save the new value
63fc1df8   David Mayerich   added an inverse ...
315
  				}
63fc1df8   David Mayerich   added an inverse ...
316
  			}
1a224b6a   David Mayerich   fixed linux compa...
317
  
63fc1df8   David Mayerich   added an inverse ...
318
  			if(PROGRESS) progress = (double)(n+1) / N * 100;	//set the current progress
6708cc25   heziqi   Ziqi added envi c...
319
  
63fc1df8   David Mayerich   added an inverse ...
320
321
  			target.write((char*)sbc, sizeof(T) * B);	//write the corrected data into destination
  		}														//end for each pixel
1a224b6a   David Mayerich   fixed linux compa...
322
  
63fc1df8   David Mayerich   added an inverse ...
323
324
  		free(s);
  		free(sbc);
88c3e636   heziqi   Ziqi added big.h
325
  		target.close();
1a224b6a   David Mayerich   fixed linux compa...
326
327
328
  
  		return true;
  
88c3e636   heziqi   Ziqi added big.h
329
  	}
1a224b6a   David Mayerich   fixed linux compa...
330
  
a23c4132   David Mayerich   Doxygen comments ...
331
332
333
334
335
  	/// 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.
c30dd3e6   David Mayerich   added vector norm...
336
  	bool ratio(std::string outname, double w, unsigned char* mask = NULL, bool PROGRESS = false)
88c3e636   heziqi   Ziqi added big.h
337
  	{
88c3e636   heziqi   Ziqi added big.h
338
  		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file
1a224b6a   David Mayerich   fixed linux compa...
339
340
  		std::string headername = outname + ".hdr";              //the header file name
  
63fc1df8   David Mayerich   added an inverse ...
341
342
343
344
345
  		unsigned long long N = X() * Y();						//calculate the total number of pixels to be processed
  		unsigned long long B = Z();								//get the number of bands
  		T* s = (T*)malloc(sizeof(T) * B);						//allocate memory to store a pixel
  		T nv;													//stores the value of the normalized band
  		for(unsigned long long n = 0; n < N; n++){				//for each pixel in the image
798cea97   David Mayerich   added progress ba...
346
  			if(mask != NULL && !mask[n])						//if the normalization band is below threshold
63fc1df8   David Mayerich   added an inverse ...
347
348
  				memset(s, 0, sizeof(T) * B);					//set the output to zero
  			else{
9d3ba0b1   David Mayerich   added stim::hsi a...
349
  				pixel(s, n);										//retrieve the spectrum s
1a224b6a   David Mayerich   fixed linux compa...
350
351
  				nv = hsi<T>::interp_spectrum(s, w);							//find the value of the normalization band
  
9d3ba0b1   David Mayerich   added stim::hsi a...
352
  				for(unsigned long long b = 0; b < B; b++)			//for each band in the spectrum
63fc1df8   David Mayerich   added an inverse ...
353
  					s[b] /= nv;										//divide by the normalization value
88c3e636   heziqi   Ziqi added big.h
354
  			}
1a224b6a   David Mayerich   fixed linux compa...
355
  
63fc1df8   David Mayerich   added an inverse ...
356
  			if(PROGRESS) progress = (double)(n+1) / N * 100;	//set the current progress
6708cc25   heziqi   Ziqi added envi c...
357
  
63fc1df8   David Mayerich   added an inverse ...
358
359
  			target.write((char*)s, sizeof(T) * B);	//write the corrected data into destination
  		}														//end for each pixel
1a224b6a   David Mayerich   fixed linux compa...
360
  
63fc1df8   David Mayerich   added an inverse ...
361
  		free(s);
88c3e636   heziqi   Ziqi added big.h
362
  		target.close();
88c3e636   heziqi   Ziqi added big.h
363
364
  		return true;
  	}
c30dd3e6   David Mayerich   added vector norm...
365
  
904614c3   David Mayerich   added normalizati...
366
367
368
369
370
371
372
373
374
375
376
377
378
379
  	void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){
  
  		std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file
  		file.seekg(0, std::ios::beg);								//move the pointer to the current file to the beginning
  
  		size_t B = Z();												//number of spectral components
  		size_t XY = X() * Y();										//calculate the number of pixels
  		size_t Bb = B * sizeof(T);									//number of bytes in a spectrum
  
  		T* spec = (T*) malloc(Bb);									//allocate space for the spectrum
  		T len;
  		for(size_t xy = 0; xy < XY; xy++){							//for each pixel
  			memset(spec, 0, Bb);									//set the spectrum to zero
  			if(mask == NULL || mask[xy]){							//if the pixel is masked
334547d8   sam   changes for mnf: ...
380
  				len = 0;											//initialize the
904614c3   David Mayerich   added normalizati...
381
382
383
384
385
386
387
388
389
390
391
  				file.read((char*)spec, Bb);							//read a spectrum
  				for(size_t b = 0; b < B; b++)						//for each band
  					len += spec[b]*spec[b];							//add the square of the spectral band
  				len = sqrt(len);									//calculate the square of the sum of squared components
  				for(size_t b = 0; b < B; b++)						//for each band
  					spec[b] /= len;									//divide by the length
  			}
  			else
  				file.seekg(Bb, std::ios::cur);						//otherwise skip a spectrum
  			target.write((char*)spec, Bb);							//output the normalized spectrum
  			if(PROGRESS) progress = (double)(xy + 1) / (double)XY * 100;		//update the progress
334547d8   sam   changes for mnf: ...
392
  		}
904614c3   David Mayerich   added normalizati...
393
  	}
1a224b6a   David Mayerich   fixed linux compa...
394
  
b31c02eb   David Mayerich   added a --select ...
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
  	/// This function loads a specified set of bands and saves them into a new output file
  	bool select(std::string outfile, std::vector<double> bandlist, unsigned char* mask = NULL, bool PROGRESS = false) {
  		std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file
  		if (!target) {
  			std::cout << "ERROR opening output file: " << outfile << std::endl;
  			return false;
  		}
  		file.seekg(0, std::ios::beg);								//move the pointer to the current file to the beginning
  
  		size_t B = Z();												//number of spectral components
  		size_t XY = X() * Y();										//calculate the number of pixels
  		size_t Bout = bandlist.size();
  		size_t in_bytes = B * sizeof(T);							//number of bytes in a spectrum
  		size_t out_bytes = Bout * sizeof(T);						//number of bytes in an output spectrum
  
  		T* in = (T*)malloc(in_bytes);								//allocate space for the input spectrum
  		T* out = (T*)malloc(out_bytes);								//allocate space for the output spectrum
  
  		double wc;													//register to store the desired wavelength
  		double w0, w1;												//registers to store the wavelengths surrounding the given band
  		size_t b0, b1;												//indices of the bands surrounding the specified wavelength
  		for (size_t xy = 0; xy < XY; xy++) {						//for each pixel
  			//memset(out, 0, out_bytes);								//set the spectrum to zero
  			if (mask == NULL || mask[xy]) {							//if the pixel is masked
  				file.read((char*)in, in_bytes);						//read an input spectrum
  				for (size_t b = 0; b < Bout; b++) {					//for each band
  					wc = bandlist[b];									//set the desired wavelength
  					band_bounds(wc, b0, b1);								//get the surrounding bands
  					w0 = w[b0];											//get the wavelength for the lower band
  					w1 = w[b1];											//get the wavelength for the higher band
  					out[b] = lerp(wc, in[b0], w0, in[b1], w1);			//interpolate the spectral values to get the desired output value
  				}
  			}
  			else
  				file.seekg(Bout, std::ios::cur);						//otherwise skip a spectrum
  			target.write((char*)out, out_bytes);							//output the normalized spectrum
  			if (PROGRESS) progress = (double)(xy + 1) / (double)XY * 100;		//update the progress
  		}
  
  		free(in);
  		free(out);
  		return true;
  	}
  
f6169dea   heziqi   Ziqi completed bi...
439
  
a23c4132   David Mayerich   Doxygen comments ...
440
441
442
  	/// Convert the current BIP file to a BIL file with the specified file name.
  
  	/// @param outname is the name of the output BIL file to be saved to disk.
63fc1df8   David Mayerich   added an inverse ...
443
  	bool bil(std::string outname, bool PROGRESS = false)
11f177d5   heziqi   Ziqi added functi...
444
  	{
9d3ba0b1   David Mayerich   added stim::hsi a...
445
  		unsigned long long S = X() * Z() * sizeof(T);		//calculate the number of bytes in a ZX slice
1a224b6a   David Mayerich   fixed linux compa...
446
  
11f177d5   heziqi   Ziqi added functi...
447
  		std::ofstream target(outname.c_str(), std::ios::binary);
9d3ba0b1   David Mayerich   added stim::hsi a...
448
  		//std::string headername = outname + ".hdr";
1a224b6a   David Mayerich   fixed linux compa...
449
  
f6169dea   heziqi   Ziqi completed bi...
450
  		T * p;			//pointer to the current ZX slice for bip file
11f177d5   heziqi   Ziqi added functi...
451
  		p = (T*)malloc(S);
f6169dea   heziqi   Ziqi completed bi...
452
453
  		T * q;			//pointer to the current XZ slice for bil file
  		q = (T*)malloc(S);
1a224b6a   David Mayerich   fixed linux compa...
454
  
9d3ba0b1   David Mayerich   added stim::hsi a...
455
  		for ( unsigned long long i = 0; i < Y(); i++)
1a224b6a   David Mayerich   fixed linux compa...
456
  		{
724ec347   David Mayerich   simplified the EN...
457
  			read_plane_y(p, i);
9d3ba0b1   David Mayerich   added stim::hsi a...
458
  			for ( unsigned long long k = 0; k < Z(); k++)
f6169dea   heziqi   Ziqi completed bi...
459
  			{
9d3ba0b1   David Mayerich   added stim::hsi a...
460
461
  				unsigned long long ks = k * X();
  				for ( unsigned long long j = 0; j < X(); j++)
724ec347   David Mayerich   simplified the EN...
462
  					q[ks + j] = p[k + j * Z()];
cbce396e   David Mayerich   added support for...
463
  
63fc1df8   David Mayerich   added an inverse ...
464
  				if(PROGRESS) progress = (double)(i * Z() + k+1) / (Y() * Z()) * 100;
f6169dea   heziqi   Ziqi completed bi...
465
  			}
1a224b6a   David Mayerich   fixed linux compa...
466
  			target.write(reinterpret_cast<const char*>(q), S);   //write a band data into target file
11f177d5   heziqi   Ziqi added functi...
467
  		}
f6169dea   heziqi   Ziqi completed bi...
468
  
11f177d5   heziqi   Ziqi added functi...
469
  		free(p);
f6169dea   heziqi   Ziqi completed bi...
470
  		free(q);
11f177d5   heziqi   Ziqi added functi...
471
472
473
  		target.close();
  		return true;
  	}
88c3e636   heziqi   Ziqi added big.h
474
  
a23c4132   David Mayerich   Doxygen comments ...
475
476
477
478
479
480
481
482
  	/// 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.
20c212c0   heziqi   Ziqi added functi...
483
484
  	bool baseline_band(double lb, double rb, T* lp, T* rp, double wavelength, T* result){
  
9d3ba0b1   David Mayerich   added stim::hsi a...
485
  		unsigned long long XY = X() * Y();
20c212c0   heziqi   Ziqi added functi...
486
487
488
489
  		band(result, wavelength);		//get band
  
  		//perform the baseline correction
  		double r = (double) (wavelength - lb) / (double) (rb - lb);
9d3ba0b1   David Mayerich   added stim::hsi a...
490
  		for(unsigned long long i=0; i < XY; i++){
20c212c0   heziqi   Ziqi added functi...
491
492
493
494
  			result[i] =(T) (result[i] - (rp[i] - lp[i]) * r - lp[i] );
  		}
  		return true;
  	}
1a224b6a   David Mayerich   fixed linux compa...
495
  
a23c4132   David Mayerich   Doxygen comments ...
496
497
498
499
500
501
  	/// 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.
70407ea9   heziqi   Ziqi completed he...
502
  	bool height(double lb, double rb, double bandwavelength, T* result){
20c212c0   heziqi   Ziqi added functi...
503
504
505
  
  		T* lp;
  		T* rp;
9d3ba0b1   David Mayerich   added stim::hsi a...
506
507
  		unsigned long long XY = X() * Y();
  		unsigned long long S = XY * sizeof(T);
70407ea9   heziqi   Ziqi completed he...
508
509
  		lp = (T*) malloc(S);			//memory allocation
  		rp = (T*) malloc(S);
20c212c0   heziqi   Ziqi added functi...
510
  
70407ea9   heziqi   Ziqi completed he...
511
  		band(lp, lb);
1a224b6a   David Mayerich   fixed linux compa...
512
  		band(rp, rb);
20c212c0   heziqi   Ziqi added functi...
513
514
515
  
  		baseline_band(lb, rb, lp, rp, bandwavelength, result);
  
517876d6   heziqi   metrics finished ...
516
517
  		free(lp);
  		free(rp);
20c212c0   heziqi   Ziqi added functi...
518
519
520
  		return true;
  	}
  
c359422d   heziqi   Ziqi added functi...
521
  
a23c4132   David Mayerich   Doxygen comments ...
522
523
524
525
526
527
528
  	/// 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
70407ea9   heziqi   Ziqi completed he...
529
  	bool area(double lb, double rb, double lab, double rab, T* result){
20c212c0   heziqi   Ziqi added functi...
530
531
532
533
534
  
  		T* lp;	//left band pointer
  		T* rp;	//right band pointer
  		T* cur;		//current band 1
  		T* cur2;	//current band 2
20c212c0   heziqi   Ziqi added functi...
535
  
9d3ba0b1   David Mayerich   added stim::hsi a...
536
537
  		unsigned long long XY = X() * Y();
  		unsigned long long S = XY * sizeof(T);
20c212c0   heziqi   Ziqi added functi...
538
539
540
541
542
  
  		lp = (T*) malloc(S);			//memory allocation
  		rp = (T*) malloc(S);
  		cur = (T*) malloc(S);
  		cur2 = (T*) malloc(S);
20c212c0   heziqi   Ziqi added functi...
543
544
545
546
  
  		memset(result, (char)0, S);
  
  		//find the wavelenght position in the whole band
9d3ba0b1   David Mayerich   added stim::hsi a...
547
548
549
  		unsigned long long n = w.size();
  		unsigned long long ai = 0;		//left bound position
  		unsigned long long bi = n - 1;		//right bound position
20c212c0   heziqi   Ziqi added functi...
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
  
  
  
  		//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
70407ea9   heziqi   Ziqi completed he...
565
  		while (lab >= w[ai]){
20c212c0   heziqi   Ziqi added functi...
566
567
  			ai++;
  		}
70407ea9   heziqi   Ziqi completed he...
568
  		while (rab <= w[bi]){
20c212c0   heziqi   Ziqi added functi...
569
570
571
572
573
574
575
  			bi--;
  		}
  
  		band(lp, lb);
  		band(rp, rb);
  
  		//calculate the beginning and the ending part
70407ea9   heziqi   Ziqi completed he...
576
577
  		baseline_band(lb, rb, lp, rp, rab, cur2);		//ending part
  		baseline_band(lb, rb, lp, rp, w[bi], cur);
9d3ba0b1   David Mayerich   added stim::hsi a...
578
579
  		for(unsigned long long j = 0; j < XY; j++){
  			result[j] += (T)((rab - w[bi]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
20c212c0   heziqi   Ziqi added functi...
580
  		}
70407ea9   heziqi   Ziqi completed he...
581
582
  		baseline_band(lb, rb, lp, rp, lab, cur2);		//beginnning part
  		baseline_band(lb, rb, lp, rp, w[ai], cur);
1a224b6a   David Mayerich   fixed linux compa...
583
  		for(unsigned long long j = 0; j < XY; j++){
9d3ba0b1   David Mayerich   added stim::hsi a...
584
  			result[j] += (T)((w[ai] - lab) * ((double)cur[j] + (double)cur2[j]) / 2.0);
20c212c0   heziqi   Ziqi added functi...
585
  		}
517876d6   heziqi   metrics finished ...
586
  
20c212c0   heziqi   Ziqi added functi...
587
588
  		//calculate the area
  		ai++;
9d3ba0b1   David Mayerich   added stim::hsi a...
589
  		for(unsigned long long i = ai; i <= bi ;i++)
20c212c0   heziqi   Ziqi added functi...
590
591
  		{
  			baseline_band(lb, rb, lp, rp, w[ai], cur2);
9d3ba0b1   David Mayerich   added stim::hsi a...
592
  			for(unsigned long long j = 0; j < XY; j++)
20c212c0   heziqi   Ziqi added functi...
593
  			{
1a224b6a   David Mayerich   fixed linux compa...
594
  				result[j] += (T)((w[ai] - w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
20c212c0   heziqi   Ziqi added functi...
595
596
597
  			}
  			std::swap(cur,cur2);		//swap the band pointers
  		}
70407ea9   heziqi   Ziqi completed he...
598
  
517876d6   heziqi   metrics finished ...
599
600
601
602
  		free(lp);
  		free(rp);
  		free(cur);
  		free(cur2);
20c212c0   heziqi   Ziqi added functi...
603
604
605
  		return true;
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
606
607
608
609
610
611
612
613
614
  	/// 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
ba51ae6a   David Mayerich   fixed metric calc...
615
  	bool ph_to_ph(T* result, double lb1, double rb1, double pos1, double lb2, double rb2, double pos2, unsigned char* mask = NULL){
70407ea9   heziqi   Ziqi completed he...
616
  
1a224b6a   David Mayerich   fixed linux compa...
617
618
619
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
70407ea9   heziqi   Ziqi completed he...
620
621
622
623
  		//get the two peak band
  		height(lb1, rb1, pos1, p1);
  		height(lb2, rb2, pos2, p2);
  		//calculate the ratio in result
9d3ba0b1   David Mayerich   added stim::hsi a...
624
  		for(unsigned long long i = 0; i < X() * Y(); i++){
70407ea9   heziqi   Ziqi completed he...
625
626
627
628
629
  			if(p1[i] == 0 && p2[i] ==0)
  				result[i] = 1;
  			else
  				result[i] = p1[i] / p2[i];
  		}
517876d6   heziqi   metrics finished ...
630
631
632
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
633
  		return true;
70407ea9   heziqi   Ziqi completed he...
634
  	}
1a224b6a   David Mayerich   fixed linux compa...
635
  
a23c4132   David Mayerich   Doxygen comments ...
636
637
638
639
640
641
642
643
644
  	/// 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
ba51ae6a   David Mayerich   fixed metric calc...
645
646
  	bool pa_to_ph(T* result, double lb1, double rb1, double lab1, double rab1,
  					double lb2, double rb2, double pos, unsigned char* mask = NULL){
1a224b6a   David Mayerich   fixed linux compa...
647
648
649
650
  
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
70407ea9   heziqi   Ziqi completed he...
651
652
653
654
  		//get the area and the peak band
  		area(lb1, rb1, lab1, rab1, p1);
  		height(lb2, rb2, pos, p2);
  		//calculate the ratio in result
9d3ba0b1   David Mayerich   added stim::hsi a...
655
  		for(unsigned long long i = 0; i < X() * Y(); i++){
70407ea9   heziqi   Ziqi completed he...
656
657
658
659
660
  			if(p1[i] == 0 && p2[i] ==0)
  				result[i] = 1;
  			else
  				result[i] = p1[i] / p2[i];
  		}
517876d6   heziqi   metrics finished ...
661
662
663
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
664
665
666
  		return true;
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
667
668
669
670
671
672
673
674
675
  	/// 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)
1a224b6a   David Mayerich   fixed linux compa...
676
  	/// @param rab2 is the label value for the right bound (end of the integration) of the second peak (denominator)
a23c4132   David Mayerich   Doxygen comments ...
677
  	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
ba51ae6a   David Mayerich   fixed metric calc...
678
679
  	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){
1a224b6a   David Mayerich   fixed linux compa...
680
681
682
683
  
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
70407ea9   heziqi   Ziqi completed he...
684
685
686
687
  		//get the area and the peak band
  		area(lb1, rb1, lab1, rab1, p1);
  		area(lb2, rb2, lab2, rab2, p2);
  		//calculate the ratio in result
9d3ba0b1   David Mayerich   added stim::hsi a...
688
  		for(unsigned long long i = 0; i < X() * Y(); i++){
70407ea9   heziqi   Ziqi completed he...
689
690
691
692
693
  			if(p1[i] == 0 && p2[i] ==0)
  				result[i] = 1;
  			else
  				result[i] = p1[i] / p2[i];
  		}
517876d6   heziqi   metrics finished ...
694
695
696
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
697
698
  		return true;
  	}
517876d6   heziqi   metrics finished ...
699
  
a23c4132   David Mayerich   Doxygen comments ...
700
701
702
703
704
705
706
  	/// 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
517876d6   heziqi   metrics finished ...
707
708
709
710
711
712
  	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
  
9d3ba0b1   David Mayerich   added stim::hsi a...
713
714
  		unsigned long long XY = X() * Y();
  		unsigned long long S = XY * sizeof(T);
517876d6   heziqi   metrics finished ...
715
716
717
718
719
720
721
722
723
  
  		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
9d3ba0b1   David Mayerich   added stim::hsi a...
724
725
726
  		unsigned long long n = w.size();
  		unsigned long long ai = 0;		//left bound position
  		unsigned long long bi = n - 1;		//right bound position
517876d6   heziqi   metrics finished ...
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
  
  		//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);
9d3ba0b1   David Mayerich   added stim::hsi a...
753
754
  		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);
517876d6   heziqi   metrics finished ...
755
756
757
  		}
  		baseline_band(lb, rb, lp, rp, lab, cur2);		//beginnning part
  		baseline_band(lb, rb, lp, rp, w[ai], cur);
1a224b6a   David Mayerich   fixed linux compa...
758
  		for(unsigned long long j = 0; j < XY; j++){
9d3ba0b1   David Mayerich   added stim::hsi a...
759
  			result[j] += (T)((w[ai] - lab) * (w[ai] + lab) * ((double)cur[j] + (double)cur2[j]) / 4.0);
517876d6   heziqi   metrics finished ...
760
761
762
763
  		}
  
  		//calculate f(x) times x
  		ai++;
9d3ba0b1   David Mayerich   added stim::hsi a...
764
  		for(unsigned long long i = ai; i <= bi ;i++)
517876d6   heziqi   metrics finished ...
765
766
  		{
  			baseline_band(lb, rb, lp, rp, w[ai], cur2);
9d3ba0b1   David Mayerich   added stim::hsi a...
767
  			for(unsigned long long j = 0; j < XY; j++)
517876d6   heziqi   metrics finished ...
768
  			{
1a224b6a   David Mayerich   fixed linux compa...
769
  				result[j] += (T)((w[ai] - w[ai-1]) * (w[ai] + w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 4.0);
517876d6   heziqi   metrics finished ...
770
771
772
773
774
775
776
777
778
779
780
  			}
  			std::swap(cur,cur2);		//swap the band pointers
  		}
  
  		free(lp);
  		free(rp);
  		free(cur);
  		free(cur2);
  		return true;
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
781
782
783
784
785
786
787
  	/// 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
ba51ae6a   David Mayerich   fixed metric calc...
788
  	bool centroid(T* result, double lb, double rb, double lab, double rab, unsigned char* mask = NULL){
1a224b6a   David Mayerich   fixed linux compa...
789
790
791
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
517876d6   heziqi   metrics finished ...
792
793
794
795
  		//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
9d3ba0b1   David Mayerich   added stim::hsi a...
796
  		for(unsigned long long i = 0; i < X() * Y(); i++){
ba51ae6a   David Mayerich   fixed metric calc...
797
  			if(mask == NULL || mask[i])
517876d6   heziqi   metrics finished ...
798
799
800
801
802
  				result[i] = p1[i] / p2[i];
  		}
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
803
  		return true;
517876d6   heziqi   metrics finished ...
804
805
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
806
807
808
809
810
811
812
  	/// 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
e37ce7bf   David Mayerich   added the ability...
813
  	bool build_mask(unsigned char* out_mask, double mask_band, double lower, double upper, unsigned char* mask = NULL, bool PROGRESS = false){
4a6f666c   heziqi   Added mask method
814
  
724ec347   David Mayerich   simplified the EN...
815
  		T* temp = (T*)malloc(X() * Y() * sizeof(T));		//allocate memory for the certain band
798cea97   David Mayerich   added progress ba...
816
  		band(temp, mask_band, PROGRESS);
4a6f666c   heziqi   Added mask method
817
  
9d3ba0b1   David Mayerich   added stim::hsi a...
818
  		for (unsigned long long i = 0; i < X() * Y();i++) {
e37ce7bf   David Mayerich   added the ability...
819
820
821
822
823
824
825
  			if(mask == NULL || mask[i] != 0){
  				if(temp[i] > lower && temp[i] < upper){
  					out_mask[i] = 255;
  				}
  				else
  					out_mask[i] = 0;
  			}
4a6f666c   heziqi   Added mask method
826
827
  		}
  
517876d6   heziqi   metrics finished ...
828
  		free(temp);
4a6f666c   heziqi   Added mask method
829
830
831
832
  		return true;
  
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
833
834
835
836
  	/// 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.
6a46c8ff   David Mayerich   fixed bugs in app...
837
  	bool apply_mask(std::string outfile, unsigned char* p, bool PROGRESS = false){
740f8cd2   heziqi   added apply_mask
838
839
840
  
  		std::ofstream target(outfile.c_str(), std::ios::binary);
  
9d3ba0b1   David Mayerich   added stim::hsi a...
841
842
  		unsigned long long ZX = Z() * X();		//calculate the number of values in a page (XZ in BIP)
  		unsigned long long L = ZX * sizeof(T);	//calculate the number of bytes in a page
740f8cd2   heziqi   added apply_mask
843
  
1a55a328   David Mayerich   Added comments fo...
844
  		T * temp = (T*)malloc(L);		//allocate space for that page
740f8cd2   heziqi   added apply_mask
845
  
9d3ba0b1   David Mayerich   added stim::hsi a...
846
  		for (unsigned long long i = 0; i < Y(); i++)			//for each page (Y in BIP)
740f8cd2   heziqi   added apply_mask
847
  		{
724ec347   David Mayerich   simplified the EN...
848
  			read_plane_y(temp, i);							//load that page (it's pointed to by temp)
9d3ba0b1   David Mayerich   added stim::hsi a...
849
  			for ( unsigned long long j = 0; j < X(); j++)	//for each X value
740f8cd2   heziqi   added apply_mask
850
  			{
9d3ba0b1   David Mayerich   added stim::hsi a...
851
  				for (unsigned long long k = 0; k < Z(); k++)	//for each B value (band)
740f8cd2   heziqi   added apply_mask
852
  				{
724ec347   David Mayerich   simplified the EN...
853
854
  					if (p[i * X() + j] == 0)	//if the mask value is zero
  					temp[j * Z() + k] = 0;			//set the pixel value to zero
1a55a328   David Mayerich   Added comments fo...
855
  				else								//otherwise just continue
740f8cd2   heziqi   added apply_mask
856
857
858
  					continue;
  				}
  			}
1a55a328   David Mayerich   Added comments fo...
859
  			target.write(reinterpret_cast<const char*>(temp), L);   //write the edited band data into target file
6a46c8ff   David Mayerich   fixed bugs in app...
860
  			if(PROGRESS) progress = (double)(i+1) / (double)Y() * 100;
740f8cd2   heziqi   added apply_mask
861
  		}
1a55a328   David Mayerich   Added comments fo...
862
863
864
  		target.close();						//close the target file
  		free(temp);							//free allocated memory
  		return true;						//return true
740f8cd2   heziqi   added apply_mask
865
866
  	}
  
ba51ae6a   David Mayerich   fixed metric calc...
867
868
869
870
871
  	/// 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
814eb271   David Mayerich   fixed problems wi...
872
  	bool sift(T* matrix, unsigned char* mask = NULL, bool PROGRESS = false){
ba51ae6a   David Mayerich   fixed metric calc...
873
874
875
876
877
  		size_t Bbytes = sizeof(T) * Z();
  		size_t XY = X() * Y();
  		T* band = (T*) malloc( Bbytes );					//allocate space for a line
  
  		file.seekg(0, std::ios::beg);	//seek to the beginning of the file
1a224b6a   David Mayerich   fixed linux compa...
878
  
ba51ae6a   David Mayerich   fixed metric calc...
879
880
  		size_t p = 0;										//create counter variables
  		for(size_t xy = 0; xy < XY; xy++){					//for each pixel
d4aa4cf7   David Mayerich   fixed the in-memo...
881
  			if(mask == NULL || mask[xy]){									//if the current pixel is masked
ba51ae6a   David Mayerich   fixed metric calc...
882
883
884
885
886
887
888
889
  				file.read( (char*)band, Bbytes );			//read the current line
  				for(size_t b = 0; b < Z(); b++){			//copy each band value to the sifted matrix
  					size_t i = p * Z() + b;					//calculate the index in the sifted matrix
  					matrix[i] = band[b];					//store the current value in the line at the correct matrix location
  				}
  				p++;									//increment the pixel pointer
  			}
  			else
1a224b6a   David Mayerich   fixed linux compa...
890
  				file.seekg(Bbytes, std::ios::cur);			//otherwise skip this band
814eb271   David Mayerich   fixed problems wi...
891
  			if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
ba51ae6a   David Mayerich   fixed metric calc...
892
893
894
  		}
  		return true;
  	}
e843658b   Brad Deutsch   Previous push did...
895
  
2fcab4f0   David Mayerich   added unsifting f...
896
  	/// Saves to disk only those spectra corresponding to mask values != 0
c1078e19   David Mayerich   HSIproc bug fixes...
897
  	bool sift(std::string outfile, unsigned char* mask, bool PROGRESS = false){
8550e94f   David Mayerich   added sifting and...
898
899
  
  		//reset the file pointer to the beginning of the file
1a224b6a   David Mayerich   fixed linux compa...
900
  		file.seekg(0, std::ios::beg);
8550e94f   David Mayerich   added sifting and...
901
902
  
  		// open an output stream
e843658b   Brad Deutsch   Previous push did...
903
904
  		std::ofstream target(outfile.c_str(), std::ios::binary);
  
8550e94f   David Mayerich   added sifting and...
905
  		//allocate space for a single spectrum
9d3ba0b1   David Mayerich   added stim::hsi a...
906
  		unsigned long long B = Z();
8550e94f   David Mayerich   added sifting and...
907
  		T* spectrum = (T*) malloc(B * sizeof(T));
e843658b   Brad Deutsch   Previous push did...
908
  
8550e94f   David Mayerich   added sifting and...
909
  		//calculate the number of pixels in a band
9d3ba0b1   David Mayerich   added stim::hsi a...
910
  		unsigned long long XY = X() * Y();
e843658b   Brad Deutsch   Previous push did...
911
  
8550e94f   David Mayerich   added sifting and...
912
  		//for each pixel
9d3ba0b1   David Mayerich   added stim::hsi a...
913
914
  		unsigned long long skip = 0;					//number of spectra to skip
  		for(unsigned long long x = 0; x < XY; x++){
8550e94f   David Mayerich   added sifting and...
915
916
917
918
919
  
  			//if the current pixel isn't masked
  			if( mask[x] == 0){
  				//increment the number of skipped pixels
  				skip++;
e843658b   Brad Deutsch   Previous push did...
920
  			}
8550e94f   David Mayerich   added sifting and...
921
922
923
924
925
926
927
928
929
930
931
932
933
934
  			//if the current pixel is masked
  			else{
  
  				//skip the intermediate pixels
  				file.seekg(skip * B * sizeof(T), std::ios::cur);
  
  				//set the skip value to zero
  				skip = 0;
  
  				//read this pixel into memory
  				file.read((char *)spectrum, B * sizeof(T));
  
  				//write this pixel out
  				target.write((char *)spectrum, B * sizeof(T));
8550e94f   David Mayerich   added sifting and...
935
  			}
9d3ba0b1   David Mayerich   added stim::hsi a...
936
  			if(PROGRESS) progress = (double) (x+1) / XY * 100;
8550e94f   David Mayerich   added sifting and...
937
  
e843658b   Brad Deutsch   Previous push did...
938
  		}
8550e94f   David Mayerich   added sifting and...
939
940
  
  		//close the output file
e843658b   Brad Deutsch   Previous push did...
941
  		target.close();
8550e94f   David Mayerich   added sifting and...
942
  		free(spectrum);
cbce396e   David Mayerich   added support for...
943
  
cbce396e   David Mayerich   added support for...
944
  		return true;
8550e94f   David Mayerich   added sifting and...
945
946
  	}
  
a2bf1d08   David Mayerich   general bug fixes...
947
  	bool unsift(std::string outfile, unsigned char* mask, unsigned long long samples, unsigned long long lines, bool PROGRESS = false){
8550e94f   David Mayerich   added sifting and...
948
949
950
951
952
  
  		// open an output stream
  		std::ofstream target(outfile.c_str(), std::ios::binary);
  
  		//reset the file pointer to the beginning of the file
1a224b6a   David Mayerich   fixed linux compa...
953
  		file.seekg(0, std::ios::beg);
8550e94f   David Mayerich   added sifting and...
954
955
  
  		//allocate space for a single spectrum
9d3ba0b1   David Mayerich   added stim::hsi a...
956
  		unsigned long long B = Z();
8550e94f   David Mayerich   added sifting and...
957
958
959
960
961
962
963
  		T* spectrum = (T*) malloc(B * sizeof(T));
  
  		//allocate space for a spectrum of zeros
  		T* zeros = (T*) malloc(B * sizeof(T));
  		memset(zeros, 0, B * sizeof(T));
  
  		//calculate the number of pixels in a band
9d3ba0b1   David Mayerich   added stim::hsi a...
964
  		unsigned long long XY = samples * lines;
8550e94f   David Mayerich   added sifting and...
965
966
  
  		//for each pixel
9d3ba0b1   David Mayerich   added stim::hsi a...
967
968
  		unsigned long long skip = 0;					//number of spectra to skip
  		for(unsigned long long x = 0; x < XY; x++){
8550e94f   David Mayerich   added sifting and...
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
  
  			//if the current pixel isn't masked
  			if( mask[x] == 0){
  
  				//write a bunch of zeros
  				target.write((char *)zeros, B * sizeof(T));
  			}
  			//if the current pixel is masked
  			else{
  
  				//read a pixel into memory
  				file.read((char *)spectrum, B * sizeof(T));
  
  				//write this pixel out
  				target.write((char *)spectrum, B * sizeof(T));
  			}
  
a2bf1d08   David Mayerich   general bug fixes...
986
  			if(PROGRESS) progress = (double)(x + 1) / XY * 100;
3150c537   David Mayerich   added further pro...
987
  
8550e94f   David Mayerich   added sifting and...
988
989
990
991
992
993
  		}
  
  		//close the output file
  		target.close();
  		free(spectrum);
  
a2bf1d08   David Mayerich   general bug fixes...
994
  		//progress = 100;
3150c537   David Mayerich   added further pro...
995
  
e843658b   Brad Deutsch   Previous push did...
996
  		return true;
8550e94f   David Mayerich   added sifting and...
997
998
  
  
e843658b   Brad Deutsch   Previous push did...
999
1000
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
1001
1002
1003
1004
  	/// 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
dd5aab2f   David Mayerich   fixed errors in s...
1005
  	bool mean_spectrum(double* m, double* std = NULL, unsigned char* mask = NULL, bool PROGRESS = false){
9d3ba0b1   David Mayerich   added stim::hsi a...
1006
1007
  		unsigned long long XY = X() * Y();							//calculate the total number of pixels in the HSI
  		T* temp = (T*)malloc(sizeof(T) * Z());						//allocate space for the current spectrum to be read
35a7195f   David Mayerich   added median spec...
1008
1009
1010
  		memset(m, 0, Z() * sizeof(double));							//set the mean spectrum to zero
  		double* e_x2 = (double*)malloc(Z() * sizeof(double));		//allocate space for E[x^2]
  		memset(e_x2, 0, Z() * sizeof(double));						//set all values for E[x^2] to zero
9d3ba0b1   David Mayerich   added stim::hsi a...
1011
1012
  
  		unsigned long long count = nnz(mask);									//calculate the number of masked pixels
35a7195f   David Mayerich   added median spec...
1013
  		double x;
9d3ba0b1   David Mayerich   added stim::hsi a...
1014
1015
1016
1017
  		for (unsigned long long i = 0; i < XY; i++){							//for each pixel in the HSI
  			if (mask == NULL || mask[i] != 0){						//if the pixel is masked
  				pixel(temp, i);										//get the spectrum
  				for (unsigned long long j = 0; j < Z(); j++){					//for each spectral component
35a7195f   David Mayerich   added median spec...
1018
1019
1020
  					x = temp[j];
  					m[j] += x / (double)count;		//add the weighted value to the average
  					e_x2[j] += x*x / (double)count;
a43c4fe1   heziqi   Added crop in env...
1021
1022
  				}
  			}
9d3ba0b1   David Mayerich   added stim::hsi a...
1023
  			if(PROGRESS) progress = (double)(i+1) / XY * 100;		//increment the progress
a43c4fe1   heziqi   Added crop in env...
1024
  		}
2aa68315   David Mayerich   major bug fixes f...
1025
  
35a7195f   David Mayerich   added median spec...
1026
  		//calculate the standard deviation
dd5aab2f   David Mayerich   fixed errors in s...
1027
1028
1029
1030
  		if (std != NULL) {
  			for (size_t i = 0; i < Z(); i++)
  				std[i] = sqrt(e_x2[i] - m[i] * m[i]);
  		}
35a7195f   David Mayerich   added median spec...
1031
  
a43c4fe1   heziqi   Added crop in env...
1032
  		free(temp);
a43c4fe1   heziqi   Added crop in env...
1033
1034
  		return true;
  	}
2781a48c   David Mayerich   fixed cuBLAS func...
1035
  //#ifdef CUDA_FOUND
c6f526c2   David Mayerich   added cuBLAS supp...
1036
  	/// Calculate the covariance matrix for masked pixels using cuBLAS
9d3ba0b1   David Mayerich   added stim::hsi a...
1037
  	/// Note that cuBLAS only supports integer-sized arrays, so there may be issues with large spectra
2781a48c   David Mayerich   fixed cuBLAS func...
1038
  	int co_matrix_cublas(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
c6f526c2   David Mayerich   added cuBLAS supp...
1039
  
1a224b6a   David Mayerich   fixed linux compa...
1040
  		cudaError_t cudaStat;
c6f526c2   David Mayerich   added cuBLAS supp...
1041
1042
  		cublasStatus_t stat;
  		cublasHandle_t handle;
1a224b6a   David Mayerich   fixed linux compa...
1043
  
c6f526c2   David Mayerich   added cuBLAS supp...
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
  		unsigned long long XY = X() * Y();									//calculate the number of elements in a band image
  		unsigned long long B = Z();											//calculate the number of spectral elements
  
  		double* s = (double*)malloc(sizeof(double) * B);					//allocate space for the spectrum that will be pulled from the file
  		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
  		cudaStat = cudaMalloc(&s_dev, B * sizeof(double));					//allocate space on the CUDA device for the spectrum
  		cudaStat = cudaMalloc(&A_dev, B * B * sizeof(double));				//allocate space on the CUDA device for the covariance matrix
  		cudaStat = cudaMemset(A_dev, 0, B * B * sizeof(double));			//initialize the covariance matrix to zero (0)
  		cudaStat = cudaMalloc(&avg_dev, B * sizeof(double));				//allocate space on the CUDA device for the average spectrum
9d3ba0b1   David Mayerich   added stim::hsi a...
1055
  		stat = cublasSetVector((int)B, sizeof(double), avg, 1, avg_dev, 1);		//copy the average spectrum to the CUDA device
1a224b6a   David Mayerich   fixed linux compa...
1056
  
c6f526c2   David Mayerich   added cuBLAS supp...
1057
1058
1059
1060
  		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)
  
  		stat = cublasCreate(&handle);										//create a cuBLAS instance
2781a48c   David Mayerich   fixed cuBLAS func...
1061
1062
  		if (stat != CUBLAS_STATUS_SUCCESS) return 1;						//test the cuBLAS instance to make sure it is valid
  
1c92dffe   David Mayerich   simplified CUDA i...
1063
  		//else std::cout<<"Using cuBLAS to calculate the mean covariance matrix..."<<std::endl;
c6f526c2   David Mayerich   added cuBLAS supp...
1064
  		for (unsigned long long xy = 0; xy < XY; xy++){										//for each pixel
5a5359da   David Mayerich   fixed - allows th...
1065
1066
  			if (mask == NULL || mask[xy] != 0){
  				pixeld(s, xy);																	//retreive the spectrum at the current xy pixel location
9d3ba0b1   David Mayerich   added stim::hsi a...
1067
1068
1069
  				stat = cublasSetVector((int)B, sizeof(double), s, 1, s_dev, 1);						//copy the spectrum from the host to the device
  				stat = cublasDaxpy(handle, (int)B, &axpy_alpha, avg_dev, 1, s_dev, 1);				//subtract the average spectrum
  				stat = cublasDsyr(handle, CUBLAS_FILL_MODE_UPPER, (int)B, &ger_alpha, s_dev, 1, A_dev, (int)B);	//calculate the covariance matrix (symmetric outer product)
5a5359da   David Mayerich   fixed - allows th...
1070
  			}
63fc1df8   David Mayerich   added an inverse ...
1071
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;													//record the current progress
5a5359da   David Mayerich   fixed - allows th...
1072
  
c6f526c2   David Mayerich   added cuBLAS supp...
1073
1074
  		}
  
9d3ba0b1   David Mayerich   added stim::hsi a...
1075
  		cublasGetMatrix((int)B, (int)B, sizeof(double), A_dev, (int)B, co, (int)B);					//copy the result from the GPU to the CPU
c6f526c2   David Mayerich   added cuBLAS supp...
1076
1077
1078
1079
1080
  
  		cudaFree(A_dev);														//clean up allocated device memory
  		cudaFree(s_dev);
  		cudaFree(avg_dev);
  
9d3ba0b1   David Mayerich   added stim::hsi a...
1081
1082
  		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++){
c6f526c2   David Mayerich   added cuBLAS supp...
1083
1084
1085
1086
  				co[B * i + j] = co[B * j + i];
  			}
  		}
  
2781a48c   David Mayerich   fixed cuBLAS func...
1087
  		return 0;
c6f526c2   David Mayerich   added cuBLAS supp...
1088
  	}
2781a48c   David Mayerich   fixed cuBLAS func...
1089
  //#endif
1a224b6a   David Mayerich   fixed linux compa...
1090
  
2aa68315   David Mayerich   major bug fixes f...
1091
  	/// Calculate the covariance matrix for all masked pixels in the image with 64-bit floating point precision.
a23c4132   David Mayerich   Doxygen comments ...
1092
1093
1094
1095
  
  	/// @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
dd5aab2f   David Mayerich   fixed errors in s...
1096
  	bool co_matrix(double* co, double* avg, unsigned char *mask, int cuda_device = 0, bool PROGRESS = false){
2781a48c   David Mayerich   fixed cuBLAS func...
1097
  		progress = 0;
c6f526c2   David Mayerich   added cuBLAS supp...
1098
  
dd5aab2f   David Mayerich   fixed errors in s...
1099
  		if(cuda_device >= 0){													//if a CUDA device is specified
dc99699e   David Mayerich   implemented a no-...
1100
1101
  			int dev_count;
  			HANDLE_ERROR(cudaGetDeviceCount(&dev_count));						//get the number of CUDA devices
dd5aab2f   David Mayerich   fixed errors in s...
1102
  			if(dev_count > 0 && dev_count > cuda_device){							//if the first device is not an emulator
1c92dffe   David Mayerich   simplified CUDA i...
1103
  				cudaDeviceProp prop;
dd5aab2f   David Mayerich   fixed errors in s...
1104
1105
  				cudaGetDeviceProperties(&prop, cuda_device);									//get the property of the requested CUDA device
  				if (prop.major != 9999) {
1c92dffe   David Mayerich   simplified CUDA i...
1106
  					std::cout << "Using CUDA device [" << cuda_device << "] to calculate the mean covariance matrix..."<<std::endl;
dd5aab2f   David Mayerich   fixed errors in s...
1107
1108
1109
1110
  					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
  				}
1c92dffe   David Mayerich   simplified CUDA i...
1111
1112
  			}																	//otherwise continue using the CPU		
  			std::cout<<"WARNING: cuBLAS failed, using CPU"<<std::endl;
dc99699e   David Mayerich   implemented a no-...
1113
  		}
a43c4fe1   heziqi   Added crop in env...
1114
  		//memory allocation
2aa68315   David Mayerich   major bug fixes f...
1115
1116
  		unsigned long long XY = X() * Y();
  		unsigned long long B = Z();
a43c4fe1   heziqi   Added crop in env...
1117
  		T* temp = (T*)malloc(sizeof(T) * B);
9d3ba0b1   David Mayerich   added stim::hsi a...
1118
1119
  
  		unsigned long long count = nnz(mask);								//count the number of masked pixels
63fc1df8   David Mayerich   added an inverse ...
1120
  
2aa68315   David Mayerich   major bug fixes f...
1121
1122
1123
1124
  		//initialize covariance matrix
  		memset(co, 0, B * B * sizeof(double));
  
  		//calculate covariance matrix
2aa68315   David Mayerich   major bug fixes f...
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
  		double* co_half = (double*) malloc(B * B * sizeof(double));			//allocate space for a higher-precision intermediate matrix
  		double* temp_precise = (double*) malloc(B * sizeof(double));
  		memset(co_half, 0, B * B * sizeof(double));							//initialize the high-precision matrix with zeros
  		unsigned long long idx;													//stores i*B to speed indexing
  		for (unsigned long long xy = 0; xy < XY; xy++){
  			if (mask == NULL || mask[xy] != 0){
  				pixel(temp, xy);												//retreive the spectrum at the current xy pixel location
  				for(unsigned long long b = 0; b < B; b++)									//subtract the mean from this spectrum and increase the precision
  					temp_precise[b] = (double)temp[b] - (double)avg[b];
  				idx = 0;
  				for (unsigned long long b0 = 0; b0 < B; b0++){								//for each band
  					for (unsigned long long b1 = b0; b1 < B; b1++)
  						co_half[idx++] += temp_precise[b0] * temp_precise[b1];
a43c4fe1   heziqi   Added crop in env...
1138
1139
  				}
  			}
63fc1df8   David Mayerich   added an inverse ...
1140
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;
a43c4fe1   heziqi   Added crop in env...
1141
  		}
2aa68315   David Mayerich   major bug fixes f...
1142
1143
1144
1145
  		idx = 0;
  		for (unsigned long long i = 0; i < B; i++){										//copy the precision matrix to both halves of the output matrix
  			for (unsigned long long j = i; j < B; j++){
  				co[j * B + i] = co[i * B + j] = co_half[idx++] / (double) count;
a43c4fe1   heziqi   Added crop in env...
1146
1147
1148
  			}
  		}
  
a43c4fe1   heziqi   Added crop in env...
1149
  		free(temp);
2aa68315   David Mayerich   major bug fixes f...
1150
  		free(temp_precise);
a43c4fe1   heziqi   Added crop in env...
1151
1152
1153
  		return true;
  	}
  
334547d8   sam   changes for mnf: ...
1154
  
1ee79b84   David Mayerich   more accurate noi...
1155
  	int coNoise_matrix_cublas(double* coN, double* avg, unsigned char *mask, bool PROGRESS = false) {
334547d8   sam   changes for mnf: ...
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
  
  		cudaError_t cudaStat;
  		cublasStatus_t stat;
  		cublasHandle_t handle;
  
  		progress = 0;													    //initialize the progress to zero (0)
  		unsigned long long XY = X() * Y();									//calculate the number of elements in a band image
  		unsigned long long B = Z();											//calculate the number of spectral elements
  
  		double* s = (double*)malloc(sizeof(double) * B);					//allocate space for the spectrum that will be pulled from the file
  		double* s_dev;														//declare a device pointer that will store the spectrum on the GPU
  
1ee79b84   David Mayerich   more accurate noi...
1168
1169
1170
1171
  		double* s2 = (double*)malloc(sizeof(double) * B);					//allocate space for the spectrum of second pixel that will be pulled from the file
  		double* s2_dev;														//  device pointer on the GPU
  		cudaStat = cudaMalloc(&s2_dev, B * sizeof(double));					//  allocate space on the CUDA device
  		cudaStat = cudaMemset(s2_dev, 0, B * sizeof(double));               //  initialize s2_dev to zero (0)
334547d8   sam   changes for mnf: ...
1172
1173
1174
1175
1176
1177
1178
1179
1180
  
  		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
  		cudaStat = cudaMalloc(&s_dev, B * sizeof(double));					//allocate space on the CUDA device for the spectrum
  		cudaStat = cudaMalloc(&A_dev, B * B * sizeof(double));				//allocate space on the CUDA device for the covariance matrix
  		cudaStat = cudaMemset(A_dev, 0, B * B * sizeof(double));			//initialize the covariance matrix to zero (0)
  		cudaStat = cudaMalloc(&avg_dev, B * sizeof(double));				//allocate space on the CUDA device for the average spectrum
  		stat = cublasSetVector((int)B, sizeof(double), avg, 1, avg_dev, 1);		//copy the average spectrum to the CUDA device
  
1ee79b84   David Mayerich   more accurate noi...
1181
  		double ger_alpha = 1.0 / (double)XY;									//scale the outer product by the inverse of the number of samples (mean outer product)
334547d8   sam   changes for mnf: ...
1182
1183
  		double axpy_alpha = -1;												//multiplication factor for the average spectrum (in order to perform a subtraction)
  
2781a48c   David Mayerich   fixed cuBLAS func...
1184
1185
1186
  		CUBLAS_HANDLE_ERROR(cublasCreate(&handle));							//create a cuBLAS instance
  		if (stat != CUBLAS_STATUS_SUCCESS) return 1;						//test the cuBLAS instance to make sure it is valid
  
1ee79b84   David Mayerich   more accurate noi...
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
  		for (unsigned long long xy = 0; xy < XY; xy++) {										//for each pixel
  			if (mask == NULL || mask[xy] != 0) {
  				pixeld(s, xy);                                        //retreive the spectrum at the current xy pixel location
  				if (xy < XY - X()) {
  					pixeld(s2, xy + X());                              //retreive the spectrum at the current xy+X pixel location, which is adjacent (bellow) to the pixel at xy location (in y direction)
  				}
  				else {
  					pixeld(s2, xy - X());                               //for the last row we consider the the adjacent pixel which is located above pixel xy
  				}
  				stat = cublasSetVector((int)B, sizeof(double), s, 1, s_dev, 1);						//copy the spectrum of first pixel from the host to the device
334547d8   sam   changes for mnf: ...
1197
1198
  				stat = cublasDaxpy(handle, (int)B, &axpy_alpha, avg_dev, 1, s_dev, 1);				//subtract the average spectrum
  
1ee79b84   David Mayerich   more accurate noi...
1199
1200
  				stat = cublasSetVector((int)B, sizeof(double), s2, 1, s2_dev, 1);					//copy the spectrum of second pixel from the host to the device
  				stat = cublasDaxpy(handle, (int)B, &axpy_alpha, avg_dev, 1, s2_dev, 1);				//subtract the average spectrum
334547d8   sam   changes for mnf: ...
1201
  
1ee79b84   David Mayerich   more accurate noi...
1202
  				stat = cublasDaxpy(handle, (int)B, &axpy_alpha, s2_dev, 1, s_dev, 1);	   //Minimum/Maximum Autocorrelation Factors (MAF) method : subtranct each pixel from adjacent pixel (in y direction)
334547d8   sam   changes for mnf: ...
1203
1204
1205
  
  				stat = cublasDsyr(handle, CUBLAS_FILL_MODE_UPPER, (int)B, &ger_alpha, s_dev, 1, A_dev, (int)B);	//calculate the covariance matrix (symmetric outer product)
  			}
1ee79b84   David Mayerich   more accurate noi...
1206
  			if (PROGRESS) progress = (double)(xy + 1) / XY * 100;													//record the current progress
334547d8   sam   changes for mnf: ...
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
  
  		}
  
  		cublasGetMatrix((int)B, (int)B, sizeof(double), A_dev, (int)B, coN, (int)B);					//copy the result from the GPU to the CPU
  
  		cudaFree(A_dev);														//clean up allocated device memory
  		cudaFree(s_dev);
  		cudaFree(s2_dev);
  		cudaFree(avg_dev);
  
1ee79b84   David Mayerich   more accurate noi...
1217
1218
  		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++) {
334547d8   sam   changes for mnf: ...
1219
1220
1221
1222
  				coN[B * i + j] = coN[B * j + i];
  			}
  		}
  
dc99699e   David Mayerich   implemented a no-...
1223
  		return 0;
334547d8   sam   changes for mnf: ...
1224
  	}
1ee79b84   David Mayerich   more accurate noi...
1225
  	//#endif
334547d8   sam   changes for mnf: ...
1226
1227
1228
1229
1230
1231
  
  	/// Calculate the covariance of noise matrix for all masked pixels in the image with 64-bit floating point precision.
  
  	/// @param coN 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
1ee79b84   David Mayerich   more accurate noi...
1232
  	bool coNoise_matrix(double* coN, double* avg, unsigned char *mask, int cuda_device = 0, bool PROGRESS = false) {
dc99699e   David Mayerich   implemented a no-...
1233
  
1c92dffe   David Mayerich   simplified CUDA i...
1234
  		if (cuda_device >= 0) {													//if a CUDA device is specified
dc99699e   David Mayerich   implemented a no-...
1235
1236
  			int dev_count;
  			HANDLE_ERROR(cudaGetDeviceCount(&dev_count));						//get the number of CUDA devices
1c92dffe   David Mayerich   simplified CUDA i...
1237
1238
1239
1240
1241
1242
1243
1244
1245
  			if (dev_count > 0 && dev_count > cuda_device) {							//if the first device is not an emulator
  				cudaDeviceProp prop;
  				cudaGetDeviceProperties(&prop, cuda_device);									//get the property of the requested CUDA device
  				if (prop.major != 9999) {
  					std::cout << "Using CUDA device [" << cuda_device << "] to calculate the noise covariance matrix..." << std::endl;
  					HANDLE_ERROR(cudaSetDevice(cuda_device));
  					int status = coNoise_matrix_cublas(coN, avg, mask, PROGRESS);			//use cuBLAS to calculate the covariance matrix
  					if (status == 0) return true;									//if the cuBLAS function returned correctly, we're done
  				}
1ee79b84   David Mayerich   more accurate noi...
1246
  			}																	//otherwise continue using the CPU
1c92dffe   David Mayerich   simplified CUDA i...
1247
  			std::cout << "WARNING: cuBLAS failed, using CPU" << std::endl;
dc99699e   David Mayerich   implemented a no-...
1248
  		}
334547d8   sam   changes for mnf: ...
1249
1250
1251
1252
1253
1254
  
  		progress = 0;
  		//memory allocation
  		unsigned long long XY = X() * Y();
  		unsigned long long B = Z();
  		T* temp = (T*)malloc(sizeof(T) * B);
1ee79b84   David Mayerich   more accurate noi...
1255
  		T* temp2 = (T*)malloc(sizeof(T) * B);
334547d8   sam   changes for mnf: ...
1256
1257
1258
  
  		unsigned long long count = nnz(mask);								//count the number of masked pixels
  
1ee79b84   David Mayerich   more accurate noi...
1259
  																			//initialize covariance matrix of noise
334547d8   sam   changes for mnf: ...
1260
1261
1262
  		memset(coN, 0, B * B * sizeof(double));
  
  		//calculate covariance matrix
1ee79b84   David Mayerich   more accurate noi...
1263
1264
1265
  		double* coN_half = (double*)malloc(B * B * sizeof(double));			//allocate space for a higher-precision intermediate matrix
  		double* temp_precise = (double*)malloc(B * sizeof(double));
  		double* temp_precise2 = (double*)malloc(B * sizeof(double));
334547d8   sam   changes for mnf: ...
1266
1267
  		memset(coN_half, 0, B * B * sizeof(double));							//initialize the high-precision matrix with zeros
  		unsigned long long idx;													//stores i*B to speed indexing
1ee79b84   David Mayerich   more accurate noi...
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
  		for (unsigned long long xy = 0; xy < XY; xy++) {
  			if (mask == NULL || mask[xy] != 0) {
  				pixel(temp, xy);                                        //retreive the spectrum at the current xy pixel location
  				if (xy < XY - X()) {
  					pixel(temp2, xy + X());                              //retreive the spectrum at the current xy+X pixel location, which is adjacent (bellow) to the pixel at xy location (in y direction)
  				}
  				else {
  					pixel(temp2, xy - X());                               //for the last row we consider the the adjacent pixel which is located above pixel xy
  				}
  				for (unsigned long long b = 0; b < B; b++) {									//subtract the mean from this spectrum and increase the precision
334547d8   sam   changes for mnf: ...
1278
  					temp_precise[b] = (double)temp[b] - (double)avg[b];
1ee79b84   David Mayerich   more accurate noi...
1279
1280
  					temp_precise2[b] = (double)temp2[b] - (double)avg[b];
  				}
334547d8   sam   changes for mnf: ...
1281
  
1ee79b84   David Mayerich   more accurate noi...
1282
1283
  				for (unsigned long long b2 = 0; b2 < B; b2++)	    //Minimum/Maximum Autocorrelation Factors (MAF) method : subtranct each pixel from adjacent pixel (in y direction)
  					temp_precise[b2] -= temp_precise2[b2];
334547d8   sam   changes for mnf: ...
1284
1285
  
  				idx = 0;
1ee79b84   David Mayerich   more accurate noi...
1286
  				for (unsigned long long b0 = 0; b0 < B; b0++) {								//for each band
334547d8   sam   changes for mnf: ...
1287
1288
1289
1290
  					for (unsigned long long b1 = b0; b1 < B; b1++)
  						coN_half[idx++] += temp_precise[b0] * temp_precise[b1];
  				}
  			}
1ee79b84   David Mayerich   more accurate noi...
1291
  			if (PROGRESS) progress = (double)(xy + 1) / XY * 100;
334547d8   sam   changes for mnf: ...
1292
1293
  		}
  		idx = 0;
1ee79b84   David Mayerich   more accurate noi...
1294
1295
1296
  		for (unsigned long long i = 0; i < B; i++) {										//copy the precision matrix to both halves of the output matrix
  			for (unsigned long long j = i; j < B; j++) {
  				coN[j * B + i] = coN[i * B + j] = coN_half[idx++] / (double)count;
334547d8   sam   changes for mnf: ...
1297
1298
1299
1300
  			}
  		}
  
  		free(temp);
1ee79b84   David Mayerich   more accurate noi...
1301
  		free(temp2);
334547d8   sam   changes for mnf: ...
1302
  		free(temp_precise);
1ee79b84   David Mayerich   more accurate noi...
1303
  		free(temp_precise2);
334547d8   sam   changes for mnf: ...
1304
1305
1306
  		return true;
  	}
  
334547d8   sam   changes for mnf: ...
1307
1308
1309
1310
1311
1312
1313
1314
      /// Project the spectra onto a set of basis functions
  	/// @param outfile is the name of the new binary output file that will be created
  	/// @param center is a spectrum about which the data set will be rotated (ex. when performing mean centering)
  	/// @param basis a set of basis vectors that the data set will be projected onto (after centering)
  	/// @param M is the number of basis vectors
  	/// @param mask is a character mask used to limit processing to valid pixels
  	bool project_cublas(std::string outfile, double* center, double* basis, unsigned long long M, unsigned char* mask = NULL, bool PROGRESS = false){
  
1c92dffe   David Mayerich   simplified CUDA i...
1315
1316
  		//cudaError_t cudaStat;
  		//cublasStatus_t stat;
334547d8   sam   changes for mnf: ...
1317
1318
1319
  		cublasHandle_t handle;
  
  		std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file
334547d8   sam   changes for mnf: ...
1320
1321
1322
1323
1324
1325
1326
  
  		progress = 0;													    //initialize the progress to zero (0)
  		unsigned long long XY = X() * Y();									//calculate the number of elements in a band image
  		unsigned long long B = Z();											//calculate the number of spectral elements
  
  		double* s = (double*)malloc(sizeof(double) * B);					//allocate space for the spectrum that will be pulled from the file
  		double* s_dev;														//declare a device pointer that will store the spectrum on the GPU
1c92dffe   David Mayerich   simplified CUDA i...
1327
  		HANDLE_ERROR(cudaMalloc(&s_dev, B * sizeof(double)));				//allocate space on the CUDA device for the spectrum
334547d8   sam   changes for mnf: ...
1328
1329
1330
  
  
          double* basis_dev;														//  device pointer on the GPU
1c92dffe   David Mayerich   simplified CUDA i...
1331
1332
          HANDLE_ERROR(cudaMalloc(&basis_dev, M * B * sizeof(double)));					//  allocate space on the CUDA device
          HANDLE_ERROR(cudaMemset(basis_dev, 0, M * B * sizeof(double)));               //  initialize basis_dev to zero (0)
334547d8   sam   changes for mnf: ...
1333
1334
1335
1336
1337
1338
1339
1340
  
  
          /// transposing basis matrix (because cuBLAS is column-major)
          double *basis_Transposed = (double*)malloc(M * B * sizeof(double));
          memset(basis_Transposed, 0, M * B * sizeof(double));
          for (int i = 0; i<M; i++)
              for (int j = 0; j<B; j++)
              basis_Transposed[i+j*M] = basis[i*B+j];
1c92dffe   David Mayerich   simplified CUDA i...
1341
1342
  		//copy the basis_Transposed matrix to the CUDA device (both matrices are stored in column-major format)
  		CUBLAS_HANDLE_ERROR(cublasSetMatrix((int)M, (int)B, sizeof(double),basis_Transposed, (int)M, basis_dev, (int)M));
334547d8   sam   changes for mnf: ...
1343
  
1c92dffe   David Mayerich   simplified CUDA i...
1344
1345
1346
  		double* center_dev;																						//declare a device pointer that will store the center (average)
  		HANDLE_ERROR(cudaMalloc(&center_dev, B * sizeof(double)));									//allocate space on the CUDA device for the center (average)
  		CUBLAS_HANDLE_ERROR(cublasSetVector((int)B, sizeof(double), center, 1, center_dev, 1));			//copy the center vector (average) to the CUDA device (from host to device)
334547d8   sam   changes for mnf: ...
1347
1348
  
  
1c92dffe   David Mayerich   simplified CUDA i...
1349
1350
1351
1352
          double* A = (double*)malloc(sizeof(double) * M);								//allocate space for the projected pixel on the host
          double* A_dev;																	//declare a device pointer that will store the projected pixel on the GPU
  		HANDLE_ERROR(cudaMalloc(&A_dev,M * sizeof(double)));				    //allocate space on the CUDA device for the projected pixel
  		HANDLE_ERROR(cudaMemset(A_dev, 0,M * sizeof(double)));		        //initialize the projected pixel to zero (0)
334547d8   sam   changes for mnf: ...
1353
1354
1355
1356
1357
  
  		double axpy_alpha = -1;												//multiplication factor for the center (in order to perform a subtraction)
  		double axpy_alpha2 = 1;												//multiplication factor for the matrix-vector multiplication
          double axpy_beta = 0;												//multiplication factor for the matrix-vector multiplication (there is no second scalor)
  
1c92dffe   David Mayerich   simplified CUDA i...
1358
  		CUBLAS_HANDLE_ERROR(cublasCreate(&handle));					//create a cuBLAS instance
334547d8   sam   changes for mnf: ...
1359
  
08d1c3b9   David Mayerich   fixed casting war...
1360
1361
1362
          T* temp = (T*)malloc(sizeof(T) * M);													//allocate space for the projected pixel to be written on the disc
  		size_t i;
  		for (unsigned long long xy = 0; xy < XY; xy++){											//for each pixel
334547d8   sam   changes for mnf: ...
1363
  			if (mask == NULL || mask[xy] != 0){
08d1c3b9   David Mayerich   fixed casting war...
1364
  				pixeld(s, xy);																	//retreive the spectrum at the current xy pixel location
334547d8   sam   changes for mnf: ...
1365
  
1c92dffe   David Mayerich   simplified CUDA i...
1366
1367
1368
1369
1370
  				CUBLAS_HANDLE_ERROR(cublasSetVector((int)B, sizeof(double), s, 1, s_dev, 1));						    //copy the spectrum from the host to the device
                  CUBLAS_HANDLE_ERROR(cublasDaxpy(handle, (int)B, &axpy_alpha, center_dev, 1, s_dev, 1));					//subtract the center (average)
                  CUBLAS_HANDLE_ERROR(cublasDgemv(handle,CUBLAS_OP_N,(int)M,(int)B,&axpy_alpha2,basis_dev,(int)M,s_dev,1,&axpy_beta,A_dev,1));         //performs the matrix-vector multiplication
                  CUBLAS_HANDLE_ERROR(cublasGetVector((int)M, sizeof(double), A_dev, 1, A, 1));							//copy the projected pixel to the host (from GPU to CPU)
  						
08d1c3b9   David Mayerich   fixed casting war...
1371
  				for(i = 0; i < M; i++)	temp[i] = (T)A[i];										//casting projected pixel from double to whatever T is
334547d8   sam   changes for mnf: ...
1372
  			}
3a880d1c   David Mayerich   fixed PCA and BIP...
1373
1374
  			else
  				memset(temp, 0, sizeof(T) * M);
334547d8   sam   changes for mnf: ...
1375
1376
1377
1378
1379
1380
1381
  
  			target.write(reinterpret_cast<const char*>(temp), sizeof(T) * M);					  //write the projected vector
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;									    //record the current progress
  
  		}
  
          //clean up allocated device memory
1c92dffe   David Mayerich   simplified CUDA i...
1382
1383
1384
1385
1386
  		HANDLE_ERROR(cudaFree(A_dev));
  		HANDLE_ERROR(cudaFree(s_dev));
  		HANDLE_ERROR(cudaFree(basis_dev));
  		HANDLE_ERROR(cudaFree(center_dev));
  		CUBLAS_HANDLE_ERROR(cublasDestroy(handle));
334547d8   sam   changes for mnf: ...
1387
1388
1389
1390
  		free(A);
  		free(s);
  		free(temp);
  		target.close();												//close the output file
08d1c3b9   David Mayerich   fixed casting war...
1391
  		
334547d8   sam   changes for mnf: ...
1392
1393
  		return true;
  	}
334547d8   sam   changes for mnf: ...
1394
  
812a2741   David Mayerich   added documentati...
1395
1396
1397
1398
1399
1400
  	/// Project the spectra onto a set of basis functions
  	/// @param outfile is the name of the new binary output file that will be created
  	/// @param center is a spectrum about which the data set will be rotated (ex. when performing mean centering)
  	/// @param basis a set of basis vectors that the data set will be projected onto (after centering)
  	/// @param M is the number of basis vectors
  	/// @param mask is a character mask used to limit processing to valid pixels
dd5aab2f   David Mayerich   fixed errors in s...
1401
  	bool project(std::string outfile, double* center, double* basis, unsigned long long M, unsigned char* mask = NULL, int cuda_device = 0, bool PROGRESS = false){
dd5aab2f   David Mayerich   fixed errors in s...
1402
1403
1404
  		if (cuda_device >= 0) {													//if a CUDA device is specified
  			int dev_count;
  			HANDLE_ERROR(cudaGetDeviceCount(&dev_count));						//get the number of CUDA devices
dd5aab2f   David Mayerich   fixed errors in s...
1405
  			if (dev_count > 0 && dev_count > cuda_device) {							//if the first device is not an emulator
1c92dffe   David Mayerich   simplified CUDA i...
1406
  				cudaDeviceProp prop;
dd5aab2f   David Mayerich   fixed errors in s...
1407
1408
  				cudaGetDeviceProperties(&prop, cuda_device);									//get the property of the requested CUDA device
  				if (prop.major != 9999) {
1c92dffe   David Mayerich   simplified CUDA i...
1409
  					std::cout << "Using CUDA device [" << cuda_device << "] to perform a basis projection..." << std::endl;
dd5aab2f   David Mayerich   fixed errors in s...
1410
  					HANDLE_ERROR(cudaSetDevice(cuda_device));
1c92dffe   David Mayerich   simplified CUDA i...
1411
  					return project_cublas(outfile, center, basis, M, mask, PROGRESS);
dd5aab2f   David Mayerich   fixed errors in s...
1412
  				}
1c92dffe   David Mayerich   simplified CUDA i...
1413
1414
  			}																	//otherwise continue using the CPU		
  			std::cout << "WARNING: cuBLAS failed, using CPU" << std::endl;
dd5aab2f   David Mayerich   fixed errors in s...
1415
1416
  		}
  		
63fc1df8   David Mayerich   added an inverse ...
1417
  		std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file
9d3ba0b1   David Mayerich   added stim::hsi a...
1418
  		//std::string headername = outfile + ".hdr";					//the header file name
63fc1df8   David Mayerich   added an inverse ...
1419
1420
1421
1422
1423
1424
1425
1426
1427
  
  		//memory allocation
  		unsigned long long XY = X() * Y();
  		unsigned long long B = Z();
  
  		T* s = (T*)malloc(sizeof(T) * B);							//allocate space for the spectrum
  		T* rs = (T*)malloc(sizeof(T) * M);							//allocate space for the projected spectrum
  		double* bv;													//pointer to the current projection vector
  		for(unsigned long long xy = 0; xy < XY; xy++){				//for each spectrum in the image
9d3ba0b1   David Mayerich   added stim::hsi a...
1428
1429
1430
1431
1432
1433
1434
1435
  			memset(rs, 0, sizeof(T) * M);
  			if(mask == NULL || mask[xy]){
  				pixel(s, xy);											//load the spectrum
  				for(unsigned long long m = 0; m < M; m++){				//for each basis vector
  					bv = &basis[m * B];									//assign 'bv' to the beginning of the basis vector
  					for(unsigned long long b = 0; b < B; b++){			//for each band
  						rs[m] += (T)(((double)s[b] - center[b]) * bv[b]);		//center the spectrum and perform the projection
  					}
63fc1df8   David Mayerich   added an inverse ...
1436
1437
  				}
  			}
1a224b6a   David Mayerich   fixed linux compa...
1438
  
63fc1df8   David Mayerich   added an inverse ...
1439
1440
1441
1442
1443
1444
1445
  			target.write(reinterpret_cast<const char*>(rs), sizeof(T) * M);					//write the projected vector
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;
  		}
  
  		free(s);													//free temporary storage arrays
  		free(rs);
  		target.close();												//close the output file
08d1c3b9   David Mayerich   fixed casting war...
1446
  		
63fc1df8   David Mayerich   added an inverse ...
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
  		return true;
  	}
  
  	bool inverse(std::string outfile, double* center, double* basis, unsigned long long B, unsigned long long C = 0, bool PROGRESS = false){
  
  		std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file
  		std::string headername = outfile + ".hdr";					//the header file name
  
  		//memory allocation
  		unsigned long long XY = X() * Y();
  		if(C == 0) C = Z();											//if no coefficient number is given, assume all are used
  		C = std::min<unsigned long long>(C, Z());					//set the number of coefficients (the user can specify fewer)
  
  		T* coeff = (T*)malloc(sizeof(T) * Z());						//allocate space for the coefficients
  		T* s = (T*)malloc(sizeof(T) * B);							//allocate space for the spectrum
  		double* bv;													//pointer to the current projection vector
  		for(unsigned long long xy = 0; xy < XY; xy++){				//for each pixel in the image (expressed as a set of coefficients)
  			pixel(coeff, xy);										//load the coefficients
  			memset(s, 0, sizeof(T) * B);							//initialize the spectrum to zero (0)
  			for(unsigned long long c = 0; c < C; c++){				//for each basis vector coefficient
  				bv = &basis[c * B];									//assign 'bv' to the beginning of the basis vector
1a224b6a   David Mayerich   fixed linux compa...
1468
  				for(unsigned long long b = 0; b < B; b++){			//for each component of the basis vector
9d3ba0b1   David Mayerich   added stim::hsi a...
1469
  					s[b] += (T)((double)coeff[c] * bv[b] + center[b]);			//calculate the contribution of each element of the basis vector in the final spectrum
63fc1df8   David Mayerich   added an inverse ...
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
  				}
  			}
  
  			target.write(reinterpret_cast<const char*>(s), sizeof(T) * B);					//write the projected vector
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;
  		}
  
  		free(coeff);												//free temporary storage arrays
  		free(s);
  		target.close();												//close the output file
  
  		return true;
  	}
  
0df38ff3   heziqi   fixed head detached
1484
  
a23c4132   David Mayerich   Doxygen comments ...
1485
1486
1487
1488
1489
1490
1491
  	/// 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
e9bddc57   David Mayerich   added band croppi...
1492
1493
1494
1495
1496
  	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,
63fc1df8   David Mayerich   added an inverse ...
1497
1498
  								   unsigned long long b1,
  								   bool PROGRESS = false){
a43c4fe1   heziqi   Added crop in env...
1499
  
3033b886   David Mayerich   implemented spect...
1500
  		//calculate the new number of samples, lines, and bands
f78ceb64   David Mayerich   fixed BIP croppin...
1501
1502
1503
  		unsigned long long samples = x1 - x0 + 1;
  		unsigned long long lines = y1 - y0 + 1;
  		unsigned long long bands = b1 - b0 + 1;
3033b886   David Mayerich   implemented spect...
1504
1505
1506
1507
1508
1509
1510
  
  		//calculate the length of one cropped spectrum
  		unsigned long long L = bands * sizeof(T);
  
  		//unsigned long long L = Z() * sizeof(T);
  
  		//allocate space for the spectrum
f78ceb64   David Mayerich   fixed BIP croppin...
1511
  		char* temp = (char*)malloc(L);
3033b886   David Mayerich   implemented spect...
1512
1513
  
  		//open an output file for binary writing
a43c4fe1   heziqi   Added crop in env...
1514
  		std::ofstream out(outfile.c_str(), std::ios::binary);
3033b886   David Mayerich   implemented spect...
1515
1516
  
  		//seek to the first pixel in the cropped image
f78ceb64   David Mayerich   fixed BIP croppin...
1517
1518
1519
1520
  		size_t startx = x0 * Z();
  		size_t starty = y0 * X() * Z();
  		size_t startb = b0;
  		file.seekg( (starty + startx + startb) * sizeof(T), std::ios::beg);
3033b886   David Mayerich   implemented spect...
1521
1522
  
  		//distance between sample spectra in the same line
f78ceb64   David Mayerich   fixed BIP croppin...
1523
1524
  		size_t dist_between_samples = Z() - bands;
  		size_t jump_sample = dist_between_samples * sizeof(T);
3033b886   David Mayerich   implemented spect...
1525
1526
  
  		//distance between sample spectra in adjacent lines
f78ceb64   David Mayerich   fixed BIP croppin...
1527
1528
1529
  		//unsigned long long jump_line = ( X() - x1 + x0 ) * Z() * sizeof(T);
  		size_t dist_between_lines = X() - samples;
  		size_t jump_line = dist_between_lines * Z() * sizeof(T);
3033b886   David Mayerich   implemented spect...
1530
1531
1532
1533
1534
  
  
  		//unsigned long long sp = y0 * X() + x0;		//start pixel
  
  		//for each pixel in the image
f78ceb64   David Mayerich   fixed BIP croppin...
1535
1536
  		for (unsigned y = 0; y < lines; y++) {
  			for (unsigned x = 0; x < samples; x++) {
3033b886   David Mayerich   implemented spect...
1537
  				//read the cropped spectral region
f78ceb64   David Mayerich   fixed BIP croppin...
1538
  				file.read(temp, L );
3033b886   David Mayerich   implemented spect...
1539
  				//pixel(temp, sp + x + y * X());
f78ceb64   David Mayerich   fixed BIP croppin...
1540
  				out.write(temp, L);   //write slice data into target file
3150c537   David Mayerich   added further pro...
1541
  
3033b886   David Mayerich   implemented spect...
1542
1543
  				file.seekg(jump_sample, std::ios::cur);
  
dd5aab2f   David Mayerich   fixed errors in s...
1544
  				if(PROGRESS) progress = (double)(y * samples + x + 1) / (lines * samples) * 100;
a43c4fe1   heziqi   Added crop in env...
1545
  			}
3033b886   David Mayerich   implemented spect...
1546
1547
  
  			file.seekg(jump_line, std::ios::cur);
a43c4fe1   heziqi   Added crop in env...
1548
1549
  		}
  		free(temp);
f78ceb64   David Mayerich   fixed BIP croppin...
1550
  		out.close();
3150c537   David Mayerich   added further pro...
1551
  
a43c4fe1   heziqi   Added crop in env...
1552
1553
1554
  		return true;
  	}
  
dc8eb8aa   David Mayerich   added band trimmi...
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
  	/// 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 B = Z();								//calculate the number of elements in a spectrum
  		size_t Bdst = Z() - band_array.size();		//calculate the number of elements in an output spectrum
  		size_t Bb = B * sizeof(T);					//calculate the number of bytes in a spectrum
  		size_t XY = X() * Y();						//calculate the number of pixels in the image
  		T* src = (T*)malloc(Bb);					//allocate space to store an input spectrum
  		T* dst = (T*)malloc(Bdst * sizeof(T));		//allocate space to store an output spectrum
  
  		size_t i;									//index into the band array
  		size_t bdst;								//index into the output array
  		for(size_t xy = 0; xy < XY; xy++){			//for each pixel
  			i = 0;
  			bdst = 0;
  			file.read((char*)src, Bb);				//read a spectrum
  			for(size_t b = 0; b < B; b++){			//for each band
  				if(b != band_array[i]){				//if the band isn't trimmed
  					dst[bdst] = src[b];				//copy the band value to the output spectrum
  					bdst++;
  				}
  				else i++;							//otherwise increment i
  			}
  			out.write((char*)dst, Bdst * sizeof(T));	//write the output spectrum
  			if(PROGRESS) progress = (double)(xy + 1) / (double) XY * 100;
  		}
  		free(src);
  		free(dst);
  	}
  
2ce6954b   David Mayerich   added the ability...
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
  	/// Combine two BIP 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, bip<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 spec_bytes = Z() * sizeof(T);						//calculate the number of bytes in a spectrum
  		T* spec = (T*)malloc(spec_bytes);							//allocate space for a spectrum
  
  		for(size_t y = 0; y < S[1]; y++){							//for each pixel in the destination image
  			for(size_t x = 0; x < S[0]; x++){
  				if(x >= p0[0] && x < p0[0] + X() && y >= p0[1] && y < p0[1] + Y())	//if this pixel is in the current image
  					file.read((char*)spec, spec_bytes);
  				else if(x >= p1[0] && x < p1[0] + C->X() && y >= p1[1] && y < p1[1] + C->Y())	//if this pixel is in the source image
  					C->file.read((char*)spec, spec_bytes);
  				else
  					memset(spec, 0, spec_bytes);
  				out.write((char*)spec, spec_bytes);					//write to the output file
  			}
1a224b6a   David Mayerich   fixed linux compa...
1622
  			if(PROGRESS) progress = (double)( (y+1) * S[0] + 1) / (double) (S[0] * S[1]) * 100;
2ce6954b   David Mayerich   added the ability...
1623
1624
1625
  		}
  	}
  
82a7c08c   David Mayerich   added an append c...
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
  	///Append two files together along the band dimension
  	void append(std::string outfile, bip<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 = Z() * sizeof(T);					//calculate the number of bytes in a single plane of this file
  		size_t b_bytes = 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 xy = 0; xy < X()*Y(); xy++) {
  			spectrum(a, xy);								//read a plane from the current file
  			out.write((char*)a, a_bytes);								//write the plane to disk
  			C->spectrum(b, xy);								//read a plane from the appending file
  			out.write((char*)b, b_bytes);
  			if (PROGRESS) progress = (double)(xy + 1) / (double)(X() * Y()) * 100;
  		}
  
  		out.close();
  	}
9b2a5d71   David Mayerich   implemented finit...
1646
1647
1648
1649
1650
1651
1652
1653
  	/// 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)
  
9d77bbd9   David Mayerich   updates for HSIvi...
1654
  	void convolve(std::string outfile, std::vector<double> C, size_t start, size_t end, unsigned char* mask = NULL, bool PROGRESS = false){
9b2a5d71   David Mayerich   implemented finit...
1655
  		std::ofstream out(outfile.c_str(), std::ios::binary);		//open the output file for writing
1a224b6a   David Mayerich   fixed linux compa...
1656
  
9b2a5d71   David Mayerich   implemented finit...
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
  		size_t N = end - start + 1;									//number of bands in the output spectrum
  		size_t Nb = N * sizeof(T);									//size of the output spectrum in bytes
  		size_t B = Z();												//calculate the number of values in a spectrum
  		size_t Bb = B * sizeof(T);									//calculate the size of a spectrum in bytes
  
  		file.seekg(0, std::ios::beg);								//move to the beginning of the input file
  
  		size_t nC = C.size();										//get the number of bands that the kernel spans
  		T* inspec = (T*)malloc(Bb);									//allocate space for the input spectrum
  		T* outspec = (T*)malloc(Nb);								//allocate space for the output spectrum
  
  		size_t XY = X() * Y();										//number of pixels in the image
  		for(size_t xy = 0; xy < XY; xy++){							//for each pixel
  			file.read((char*)inspec, Bb);							//read an input spectrum
  			memset(outspec, 0, Nb);									//set the output spectrum to zero (0)
9d77bbd9   David Mayerich   updates for HSIvi...
1672
1673
1674
1675
1676
  			if(mask == NULL || mask[xy]){
  				for(size_t b = 0; b < N; b++){							//for each component of the spectrum
  					for(size_t c = 0; c < nC; c++){						//for each coefficient in the kernel
  						outspec[b] += (T)(inspec[b + start + c] * C[c]);		//perform the sum/multiply part of the convolution
  					}
9b2a5d71   David Mayerich   implemented finit...
1677
1678
1679
1680
1681
1682
1683
  				}
  			}
  			out.write((char*)outspec, Nb);							//output the band
  			if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
  		}
  	}
  
9d77bbd9   David Mayerich   updates for HSIvi...
1684
  	void deriv(std::string outfile, size_t d, size_t order, const std::vector<double> w, unsigned char* mask = NULL, bool PROGRESS = false){
9b2a5d71   David Mayerich   implemented finit...
1685
  		std::ofstream out(outfile.c_str(), std::ios::binary);		//open the output file for writing
1a224b6a   David Mayerich   fixed linux compa...
1686
1687
  
  
9b2a5d71   David Mayerich   implemented finit...
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
  		size_t B = Z();												//calculate the number of values in a spectrum
  		size_t Bb = B * sizeof(T);									//calculate the size of a spectrum in bytes
  
  		bool UNIFORM = true;
  		double ds = w[1] - w[0];									//initialize ds
  		for(size_t b = 1; b < B; b++)								//test to see if the spectral spacing is uniform (if it is, convolution is much faster)
  			if(w[b] - w[b-1] != ds) UNIFORM = false;
  
  		size_t nC = order + d;									//approximating a derivative requires order + d samples
  
  		file.seekg(0, std::ios::beg);								//move to the beginning of the input file
  
  		T* inspec = (T*)malloc(Bb);									//allocate space for the input spectrum
  		T* outspec = (T*)malloc(Bb);								//allocate space for the output spectrum
  
  		size_t XY = X() * Y();										//number of pixels in the image
  		size_t mid = (size_t)(nC / 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 xy = 0; xy < XY; xy++){							//for each pixel
  			file.read((char*)inspec, Bb);							//read an input spectrum
  			memset(outspec, 0, Bb);									//set the output spectrum to zero (0)
9d77bbd9   David Mayerich   updates for HSIvi...
1709
1710
1711
  			if(mask == NULL || mask[xy]){
  				iw = 0;
  				for(size_t b = 0; b < mid; b++){							//for each component of the spectrum
9b2a5d71   David Mayerich   implemented finit...
1712
1713
  					std::vector<double> w_pts(w.begin() + iw, w.begin() + iw + nC);			//get the wavelengths corresponding to each sample
  					std::vector<double> C = diff_coefficients(w[b], w_pts, d);					//get the optimal sample weights
9d77bbd9   David Mayerich   updates for HSIvi...
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
  					for(size_t c = 0; c < nC; c++)						//for each coefficient in the kernel
  						outspec[b] += (T)(inspec[iw + c] * C[c]);		//perform the sum/multiply part of the convolution
  				}
  				std::vector<double> w_pts(w.begin(), w.begin() + nC);			//get the wavelengths corresponding to each sample
  				std::vector<double> C = diff_coefficients(w[0], w_pts, d);					//get the optimal sample weights
  				for(size_t b = mid; b <= B - (nC - mid); b++){
  					iw = b - mid;
  					if(!UNIFORM){																//if the spacing is non-uniform, we have to re-calculate these points every iteration
  						std::vector<double> w_pts(w.begin() + iw, w.begin() + iw + nC);			//get the wavelengths corresponding to each sample
  						std::vector<double> C = diff_coefficients(w[b], w_pts, d);					//get the optimal sample weights
  					}
  					for(size_t c = 0; c < nC; c++)						//for each coefficient in the kernel
  						outspec[b] += (T)(inspec[iw + c] * C[c]);		//perform the sum/multiply part of the convolution
  				}
  				iw = B - nC;
  				for(size_t b = B - (nC - mid) + 1; b < B; b++){
  					std::vector<double> w_pts(w.begin() + iw, w.begin() + iw + nC);			//get the wavelengths corresponding to each sample
  					std::vector<double> C = diff_coefficients(w[b], w_pts, d);					//get the optimal sample weights
  					for(size_t c = 0; c < nC; c++)						//for each coefficient in the kernel
  						outspec[b] += (T)(inspec[iw + c] * C[c]);		//perform the sum/multiply part of the convolution
9b2a5d71   David Mayerich   implemented finit...
1734
  				}
9b2a5d71   David Mayerich   implemented finit...
1735
  			}
9b2a5d71   David Mayerich   implemented finit...
1736
1737
1738
1739
1740
  			out.write((char*)outspec, Bb);							//output the band
  			if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
  		}
  	}
  
278c622b   David Mayerich   finished Windows ...
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
  	bool multiply(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){
  		std::ofstream target(outname.c_str(), std::ios::binary);		//open the target binary file
  		std::string headername = outname + ".hdr";						//the header file name
  
  		unsigned long long N = X() * Y();								//calculate the total number of pixels to be processed
  		unsigned long long B = Z();										//get the number of bands
  		T* s = (T*)malloc(sizeof(T) * B);								//allocate memory to store a pixel
  		for(unsigned long long n = 0; n < N; n++){						//for each pixel in the image
  			if(mask == NULL || mask[n]){								//if the pixel is masked
  				for(size_t b = 0; b < B; b++)							//for each band in the spectrum
  					s[b] *= (T)v;											//multiply
  			}
  
  			if(PROGRESS) progress = (double)(n+1) / N * 100;			//set the current progress
  
  			target.write((char*)s, sizeof(T) * B);						//write the corrected data into destination
  		}																//end for each pixel
  
  		free(s);														//free the spectrum
  		target.close();													//close the output file
  		return true;
  	}
  
  	bool add(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){
  		std::ofstream target(outname.c_str(), std::ios::binary);		//open the target binary file
  		std::string headername = outname + ".hdr";						//the header file name
  
  		unsigned long long N = X() * Y();								//calculate the total number of pixels to be processed
  		unsigned long long B = Z();										//get the number of bands
  		T* s = (T*)malloc(sizeof(T) * B);								//allocate memory to store a pixel
  		for(unsigned long long n = 0; n < N; n++){						//for each pixel in the image
  			if(mask == NULL || mask[n]){								//if the pixel is masked
  				for(size_t b = 0; b < B; b++)							//for each band in the spectrum
  					s[b] += (T)v;											//multiply
  			}
  
  			if(PROGRESS) progress = (double)(n+1) / N * 100;			//set the current progress
  
  			target.write((char*)s, sizeof(T) * B);						//write the corrected data into destination
  		}																//end for each pixel
  
  		free(s);														//free the spectrum
  		target.close();													//close the output file
  		return true;
  	}
  
9d34229e   David Mayerich   restructured the ...
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
  	int fft(std::string outname, size_t bandmin, size_t bandmax, size_t samples = 0, T* ratio = NULL, size_t rx = 0, size_t ry = 0, bool PROGRESS = false, int device = 0){
  		if(device == -1){
  			std::cout<<"ERROR: GPU required for FFT (uses cuFFT)."<<std::endl;
  			exit(1);
  		}
  		if(samples == 0) samples = Z();								//if samples are specified, use all of them
  		if(samples > Z()){
  			std::cout<<"ERROR: stim::envi doesn't support FFT padding just yet."<<std::endl;
  			exit(1);
  		}
  		int nd;															//stores the number of CUDA devices
  		HANDLE_ERROR(cudaGetDeviceCount(&nd));							//get the number of CUDA devices
  		if(device >= nd){												//test for the existence of the requested device
  			std::cout<<"ERROR: requested CUDA device for stim::envi::fft() doesn't exist"<<std::endl;
  			exit(1);
  		}
  		HANDLE_ERROR(cudaSetDevice(device));							//set the CUDA device
  		cudaDeviceProp prop;
  		HANDLE_ERROR(cudaGetDeviceProperties(&prop, device));			//get the CUDA device properties
  
  		size_t B = Z();
  		size_t S = samples;
  		size_t fft_size = S * sizeof(T);								//number of bytes for each FFT
  		size_t cuda_bytes = prop.totalGlobalMem;						//get the number of bytes of global memory available
  		size_t cuda_use = (size_t)floor(cuda_bytes * 0.2);								//only use 80%
  		size_t nS = cuda_use / fft_size;								//calculate the number of spectra that can be loaded onto the GPU as a single batch
  		size_t batch_bytes = nS * fft_size;								//calculate the size of a batch (in bytes)
  		size_t fft_bytes = nS * (S/2 + 1) * sizeof(cufftComplex);
  		T* batch = (T*) malloc(batch_bytes);							//allocate space in host memory to store a batch
  		memset(batch, 0, batch_bytes);
  		std::complex<T>* batch_fft = (std::complex<T>*) malloc(fft_bytes);
  		T* gpu_batch;													//device pointer to the batch
  		HANDLE_ERROR(cudaMalloc(&gpu_batch, batch_bytes));				//allocate space on the device for the FFT batch
  		cufftComplex* gpu_batch_fft;												//allocate space for the FFT result
  		HANDLE_ERROR(cudaMalloc(&gpu_batch_fft, fft_bytes));
  		int N[1];														//create an array with the interferogram size (required for cuFFT input)
  		N[0] = (int)S;													//set the only array value to the interferogram size
  
  		//if a background is provided for a ratio
  		std::complex<T>* ratio_fft = NULL;											//create a pointer for the FFT of the ratio image (if it exists)
  		if(ratio){
  			size_t bkg_bytes = rx * ry * S * sizeof(T);								//calculate the total number of bytes in the background image
  			T* bkg_copy = (T*) malloc(bkg_bytes);									//allocate space to copy the background
  			if(S == Z()) memcpy(bkg_copy, ratio, bkg_bytes);						//if the number of samples used in processing equals the number of available samples
  			else{
  				for(size_t xyi = 0; xyi < rx*ry; xyi++)
  					memcpy(&bkg_copy[xyi * S], &ratio[xyi * B], S * sizeof(T));
  			}
  			T* gpu_ratio;
  			HANDLE_ERROR(cudaMalloc(&gpu_ratio, bkg_bytes));
  			HANDLE_ERROR(cudaMemcpy(gpu_ratio, bkg_copy, bkg_bytes, cudaMemcpyHostToDevice));
  			cufftHandle bkg_plan;
  			CUFFT_HANDLE_ERROR(cufftPlanMany(&bkg_plan, 1, N, NULL, 1, N[0], NULL, 1, N[0], CUFFT_R2C, (int)(rx * ry)));
  			size_t bkg_fft_bytes = rx * ry * (S / 2 + 1) * sizeof(cufftComplex);
  			T* gpu_ratio_fft;
  			HANDLE_ERROR(cudaMalloc(&gpu_ratio_fft, bkg_fft_bytes));
  			CUFFT_HANDLE_ERROR(cufftExecR2C(bkg_plan, (cufftReal*)gpu_ratio, (cufftComplex*)gpu_ratio_fft));
  			ratio_fft = (std::complex<T>*) malloc(bkg_fft_bytes);
  			HANDLE_ERROR(cudaMemcpy(ratio_fft, gpu_ratio_fft, bkg_fft_bytes, cudaMemcpyDeviceToHost));
  			HANDLE_ERROR(cudaFree(gpu_ratio));
  			HANDLE_ERROR(cudaFree(gpu_ratio_fft));
  			CUFFT_HANDLE_ERROR(cufftDestroy(bkg_plan));
  		}
dc8eb8aa   David Mayerich   added band trimmi...
1850
  
9d34229e   David Mayerich   restructured the ...
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
  		cufftHandle plan;												//create a CUFFT plan
  		CUFFT_HANDLE_ERROR(cufftPlanMany(&plan, 1, N, NULL, 1, N[0], NULL, 1, N[0], CUFFT_R2C, (int)nS));
  
  		std::ofstream outfile(outname, std::ios::binary);				//open a file for writing
  
  		size_t XY = X() * Y();											//calculate the number of spectra
  		size_t xy = 0;
  		size_t bs;														//stores the number of spectra in the current batch
  		size_t s, b;
  		size_t S_fft = S/2 + 1;
  		size_t bandkeep = bandmax - bandmin + 1;
  		size_t x, y;
  		size_t ratio_i;
  		T* temp_spec = (T*) malloc(Z() * sizeof(T));					//allocate space to hold a single pixel
  		while(xy < XY){													//while there are unprocessed spectra
  			bs = min(XY - xy, nS);										//calculate the number of spectra to include in the batch
  			for(s = 0; s < bs; s++){									//for each spectrum in the batch
  				pixel(temp_spec, xy + s);						//read a pixel from disk
  				memcpy(&batch[s * S], temp_spec, S * sizeof(T));
  				//pixel(&batch[s * S], xy + s);							//read the next spectrum
  			}
  			HANDLE_ERROR(cudaMemcpy(gpu_batch, batch, batch_bytes, cudaMemcpyHostToDevice));
  			CUFFT_HANDLE_ERROR(cufftExecR2C(plan, (cufftReal*)gpu_batch, gpu_batch_fft));			//execute the (implicitly forward) transform
  			HANDLE_ERROR(cudaMemcpy(batch_fft, gpu_batch_fft, fft_bytes, cudaMemcpyDeviceToHost));	//copy the data back to the GPU
  			for(s = 0; s < bs; s++){															//for each spectrum in the batch
  				y = (xy + s)/X();
  				x = xy + s - y * X();
  				if(ratio_fft)	ratio_i = (y % ry) * rx + (x % rx);								//if a background is used, calculate the coordinates into it
  				for(b = 0; b < S/2 + 1; b++){														//for each sample
  					if(ratio_fft)						
  						batch[s * S + b] = -log(abs(batch_fft[s * S_fft + b]) / abs(ratio_fft[ratio_i * S_fft + b]));
  					else
  						batch[s * S + b] = abs(batch_fft[s * S_fft + b]);		//calculate the magnitude of the spectrum					
  				}
  				outfile.write((char*)&batch[s * S + bandmin], bandkeep * sizeof(T));							//save the resulting spectrum
  			}
  			xy += bs;													//increment xy by the number of spectra processed
  			if(PROGRESS) progress = (double)xy / (double)XY * 100;
  		}
  		outfile.close();
  		free(ratio_fft);
  		free(batch_fft);
  		free(batch);
  		HANDLE_ERROR(cudaFree(gpu_batch));
  		HANDLE_ERROR(cudaFree(gpu_batch_fft));
  		return 0;
  	}
0df38ff3   heziqi   fixed head detached
1898
  
a23c4132   David Mayerich   Doxygen comments ...
1899
  	/// Close the file.
f6169dea   heziqi   Ziqi completed bi...
1900
1901
1902
1903
1904
  	bool close(){
  		file.close();
  		return true;
  	}
  
88c3e636   heziqi   Ziqi added big.h
1905
  	};
6aa04ba2   David Mayerich   interleave types
1906
  }
e8eb202f   David Mayerich   added a new ENVI ...
1907
1908
  
  #endif