Blame view

stim/envi/bip.h 52.2 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
8
9
  #include <cstring>
  #include <utility>
  
c6f526c2   David Mayerich   added cuBLAS supp...
10
  //CUDA
a2bf1d08   David Mayerich   general bug fixes...
11
12
13
14
  #ifdef CUDA_FOUND
  	#include <cuda_runtime.h>
  	#include "cublas_v2.h"
  #endif
c6f526c2   David Mayerich   added cuBLAS supp...
15
  
8a86bd56   David Mayerich   changed rts names...
16
  namespace stim{
88c3e636   heziqi   Ziqi added big.h
17
  
a23c4132   David Mayerich   Doxygen comments ...
18
19
20
21
22
23
24
25
  /**
  	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
26
27
  template <typename T>
  
9d3ba0b1   David Mayerich   added stim::hsi a...
28
  class bip: public hsi<T> {
88c3e636   heziqi   Ziqi added big.h
29
30
  
  protected:
1a224b6a   David Mayerich   fixed linux compa...
31
32
  
  
9d3ba0b1   David Mayerich   added stim::hsi a...
33
34
  	//std::vector<double> w; //band wavelength
  	unsigned long long offset;		//header offset
88c3e636   heziqi   Ziqi added big.h
35
  
9d3ba0b1   David Mayerich   added stim::hsi a...
36
37
  	using hsi<T>::w;				//use the wavelength array in stim::hsi
  	using hsi<T>::nnz;
63fc1df8   David Mayerich   added an inverse ...
38
  	using binary<T>::progress;
1a224b6a   David Mayerich   fixed linux compa...
39
40
41
42
  	using hsi<T>::X;
  	using hsi<T>::Y;
  	using hsi<T>::Z;
  
88c3e636   heziqi   Ziqi added big.h
43
44
45
46
  public:
  
  	using binary<T>::open;
  	using binary<T>::file;
6708cc25   heziqi   Ziqi added envi c...
47
  	using binary<T>::R;
63fc1df8   David Mayerich   added an inverse ...
48
  	using binary<T>::read_line_0;
88c3e636   heziqi   Ziqi added big.h
49
  
9d3ba0b1   David Mayerich   added stim::hsi a...
50
51
  	bip(){ hsi<T>::init_bip(); }
  
a23c4132   David Mayerich   Doxygen comments ...
52
53
54
55
56
57
58
59
  	/// 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...
60
61
62
63
64
  	bool open(std::string filename,
  			  unsigned long long X,
  			  unsigned long long Y,
  			  unsigned long long B,
  			  unsigned long long header_offset,
9d3ba0b1   David Mayerich   added stim::hsi a...
65
  			  std::vector<double> wavelengths){
88c3e636   heziqi   Ziqi added big.h
66
  
6708cc25   heziqi   Ziqi added envi c...
67
68
69
70
  		//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
71
  
9d3ba0b1   David Mayerich   added stim::hsi a...
72
  		return open(filename, vec<unsigned long long>(B, X, Y), header_offset);
1a224b6a   David Mayerich   fixed linux compa...
73
  
88c3e636   heziqi   Ziqi added big.h
74
75
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
76
77
78
79
  	/// 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...
80
  	bool band_index( T * p, unsigned long long page, bool PROGRESS = false){
63fc1df8   David Mayerich   added an inverse ...
81
  		return binary<T>::read_plane_0(p, page, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
82
83
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
84
85
86
87
  	/// 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 ...
88
  	bool band( T * p, double wavelength, bool PROGRESS = false){
88c3e636   heziqi   Ziqi added big.h
89
  
ee4dea28   David Mayerich   fixed errors in c...
90
91
  		//if there are no wavelengths in the BSQ file
  		if(w.size() == 0)
9d3ba0b1   David Mayerich   added stim::hsi a...
92
  			return band_index(p, (unsigned long long)wavelength, PROGRESS);
ee4dea28   David Mayerich   fixed errors in c...
93
  
9d3ba0b1   David Mayerich   added stim::hsi a...
94
  		unsigned long long XY = X() * Y();	//calculate the number of pixels in a band
88c3e636   heziqi   Ziqi added big.h
95
96
  
  		unsigned page=0;                      //bands around the wavelength
517876d6   heziqi   metrics finished ...
97
  
88c3e636   heziqi   Ziqi added big.h
98
99
100
101
  
  		//get the bands numbers around the wavelength
  
  		//if wavelength is smaller than the first one in header file
6708cc25   heziqi   Ziqi added envi c...
102
  		if ( w[page] > wavelength ){
63fc1df8   David Mayerich   added an inverse ...
103
  			band_index(p, page, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
104
105
106
  			return true;
  		}
  
6708cc25   heziqi   Ziqi added envi c...
107
  		while( w[page] < wavelength )
88c3e636   heziqi   Ziqi added big.h
108
109
110
  		{
  			page++;
  			//if wavelength is larger than the last wavelength in header file
724ec347   David Mayerich   simplified the EN...
111
  			if (page == Z()) {
63fc1df8   David Mayerich   added an inverse ...
112
  				band_index(p, Z()-1, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
113
114
115
  				return true;
  			}
  		}
6708cc25   heziqi   Ziqi added envi c...
116
  		if ( wavelength < w[page] )
88c3e636   heziqi   Ziqi added big.h
117
  		{
517876d6   heziqi   metrics finished ...
118
119
  			T * p1;
  			T * p2;
88c3e636   heziqi   Ziqi added big.h
120
121
122
  			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 ...
123
  			band_index(p2, page, PROGRESS);
9d3ba0b1   David Mayerich   added stim::hsi a...
124
  			for(unsigned long long i=0; i < XY; i++){
6708cc25   heziqi   Ziqi added envi c...
125
  				double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
9d3ba0b1   David Mayerich   added stim::hsi a...
126
  				p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]);
88c3e636   heziqi   Ziqi added big.h
127
  			}
517876d6   heziqi   metrics finished ...
128
129
  			free(p1);
  			free(p2);
88c3e636   heziqi   Ziqi added big.h
130
131
132
  		}
  		else                           //if the wavelength is equal to a wavelength in header file
  		{
63fc1df8   David Mayerich   added an inverse ...
133
  			band_index(p, page, PROGRESS);
88c3e636   heziqi   Ziqi added big.h
134
  		}
88c3e636   heziqi   Ziqi added big.h
135
136
  		return true;
  	}
724ec347   David Mayerich   simplified the EN...
137
138
139
140
141
142
  
  	/// 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...
143
  	bool spectrum(T * p, unsigned long long x, unsigned long long y, bool PROGRESS = false){
63fc1df8   David Mayerich   added an inverse ...
144
  		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...
145
146
147
148
149
150
151
152
  	}
  
  	/// 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...
153
  	{
9d3ba0b1   David Mayerich   added stim::hsi a...
154
  		unsigned long long B = Z();
11f177d5   heziqi   Ziqi added functi...
155
  
9d3ba0b1   David Mayerich   added stim::hsi a...
156
  		unsigned long long page=0;                      //samples around the wavelength
517876d6   heziqi   metrics finished ...
157
  
11f177d5   heziqi   Ziqi added functi...
158
159
160
161
  
  		//get the bands numbers around the wavelength
  
  		//if wavelength is smaller than the first one in header file
6708cc25   heziqi   Ziqi added envi c...
162
  		if ( w[page] > wavelength ){
9d3ba0b1   David Mayerich   added stim::hsi a...
163
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
164
165
  			{
  				p[j] = c[j * B];
1a224b6a   David Mayerich   fixed linux compa...
166
  			}
11f177d5   heziqi   Ziqi added functi...
167
168
169
  			return true;
  		}
  
6708cc25   heziqi   Ziqi added envi c...
170
  		while( w[page] < wavelength )
11f177d5   heziqi   Ziqi added functi...
171
172
173
174
  		{
  			page++;
  			//if wavelength is larger than the last wavelength in header file
  			if (page == B) {
9d3ba0b1   David Mayerich   added stim::hsi a...
175
  				for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
176
177
178
179
180
181
  				{
  					p[j] = c[(j + 1) * B - 1];
  				}
  				return true;
  			}
  		}
6708cc25   heziqi   Ziqi added envi c...
182
  		if ( wavelength < w[page] )
11f177d5   heziqi   Ziqi added functi...
183
  		{
517876d6   heziqi   metrics finished ...
184
185
  			T * p1;
  			T * p2;
724ec347   David Mayerich   simplified the EN...
186
187
  			p1=(T*)malloc( X() * sizeof(T));                     //memory allocation
  			p2=(T*)malloc( X() * sizeof(T));
11f177d5   heziqi   Ziqi added functi...
188
  			//band_index(p1, page - 1);
9d3ba0b1   David Mayerich   added stim::hsi a...
189
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
190
191
192
193
  			{
  				p1[j] = c[j * B + page - 1];
  			}
  			//band_index(p2, page );
9d3ba0b1   David Mayerich   added stim::hsi a...
194
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
195
196
197
  			{
  				p2[j] = c[j * B + page];
  			}
1a224b6a   David Mayerich   fixed linux compa...
198
  
9d3ba0b1   David Mayerich   added stim::hsi a...
199
  			for(unsigned long long i=0; i < X(); i++){
6708cc25   heziqi   Ziqi added envi c...
200
  				double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
11f177d5   heziqi   Ziqi added functi...
201
202
  				p[i] = (p2[i] - p1[i]) * r + p1[i];
  			}
517876d6   heziqi   metrics finished ...
203
204
  			free(p1);
  			free(p2);
11f177d5   heziqi   Ziqi added functi...
205
206
207
208
  		}
  		else                           //if the wavelength is equal to a wavelength in header file
  		{
  			//band_index(p, page);
9d3ba0b1   David Mayerich   added stim::hsi a...
209
  			for(unsigned long long j = 0; j < X(); j++)
11f177d5   heziqi   Ziqi added functi...
210
211
212
213
214
  			{
  				p[j] = c[j * B + page];
  			}
  		}
  
1a224b6a   David Mayerich   fixed linux compa...
215
  		return true;
11f177d5   heziqi   Ziqi added functi...
216
  	}
bfe9c6db   heziqi   added feature_mat...
217
  
c6f526c2   David Mayerich   added cuBLAS supp...
218
  	/// Retrieve a single pixel and store it in a pre-allocated double array.
9d3ba0b1   David Mayerich   added stim::hsi a...
219
220
  	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...
221
222
223
224
  		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...
225
  		unsigned long long B = Z();
c6f526c2   David Mayerich   added cuBLAS supp...
226
227
228
229
230
  
  		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...
231
  		for(unsigned long long i = 0; i < B; i++)						//for each element of the spectrum
c6f526c2   David Mayerich   added cuBLAS supp...
232
  			p[i] = (double) temp[i];							//cast each element to a double value
9d3ba0b1   David Mayerich   added stim::hsi a...
233
  		return true;
c6f526c2   David Mayerich   added cuBLAS supp...
234
235
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
236
237
238
239
  	/// 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...
240
  	bool pixel(T * p, unsigned long long n){
bfe9c6db   heziqi   added feature_mat...
241
  
9d3ba0b1   David Mayerich   added stim::hsi a...
242
  		unsigned long long N = X() * Y();					//calculate numbers in one band
63fc1df8   David Mayerich   added an inverse ...
243
  		if ( n >= N){							//make sure the pixel number is right
bfe9c6db   heziqi   added feature_mat...
244
245
246
247
  			std::cout<<"ERROR: sample or line out of range"<<std::endl;
  			return false;
  		}
  
724ec347   David Mayerich   simplified the EN...
248
249
250
  		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...
251
  	}
1a224b6a   David Mayerich   fixed linux compa...
252
  
11f177d5   heziqi   Ziqi added functi...
253
  	//given a Y ,return a ZX slice
9d3ba0b1   David Mayerich   added stim::hsi a...
254
  	bool read_plane_y(T * p, unsigned long long y){
724ec347   David Mayerich   simplified the EN...
255
  		return binary<T>::read_plane_2(p, y);
11f177d5   heziqi   Ziqi added functi...
256
  	}
1a224b6a   David Mayerich   fixed linux compa...
257
  
a23c4132   David Mayerich   Doxygen comments ...
258
259
260
  	/// 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...
261
  	/// @param wls is the list of baseline points based on band labels.
798cea97   David Mayerich   added progress ba...
262
  	bool baseline(std::string outname, std::vector<double> base_pts, unsigned char* mask = NULL, bool PROGRESS = false){
1a224b6a   David Mayerich   fixed linux compa...
263
  
88c3e636   heziqi   Ziqi added big.h
264
  		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file
1a224b6a   David Mayerich   fixed linux compa...
265
  
63fc1df8   David Mayerich   added an inverse ...
266
267
268
269
270
271
272
273
  		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 ...
274
275
  		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...
276
277
278
279
  			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...
280
  
798cea97   David Mayerich   added progress ba...
281
  				pixel(s, n);										//retrieve the spectrum s
1a224b6a   David Mayerich   fixed linux compa...
282
  				base_vals = hsi<T>::interp_spectrum(s, base_pts);			//get the values at each baseline point
798cea97   David Mayerich   added progress ba...
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
  
  				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...
304
  					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 ...
305
  				}
63fc1df8   David Mayerich   added an inverse ...
306
  			}
1a224b6a   David Mayerich   fixed linux compa...
307
  
63fc1df8   David Mayerich   added an inverse ...
308
  			if(PROGRESS) progress = (double)(n+1) / N * 100;	//set the current progress
6708cc25   heziqi   Ziqi added envi c...
309
  
63fc1df8   David Mayerich   added an inverse ...
310
311
  			target.write((char*)sbc, sizeof(T) * B);	//write the corrected data into destination
  		}														//end for each pixel
1a224b6a   David Mayerich   fixed linux compa...
312
  
63fc1df8   David Mayerich   added an inverse ...
313
314
  		free(s);
  		free(sbc);
88c3e636   heziqi   Ziqi added big.h
315
  		target.close();
1a224b6a   David Mayerich   fixed linux compa...
316
317
318
  
  		return true;
  
88c3e636   heziqi   Ziqi added big.h
319
  	}
1a224b6a   David Mayerich   fixed linux compa...
320
  
a23c4132   David Mayerich   Doxygen comments ...
321
322
323
324
325
  	/// 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...
326
  	bool ratio(std::string outname, double w, unsigned char* mask = NULL, bool PROGRESS = false)
88c3e636   heziqi   Ziqi added big.h
327
  	{
88c3e636   heziqi   Ziqi added big.h
328
  		std::ofstream target(outname.c_str(), std::ios::binary);	//open the target binary file
1a224b6a   David Mayerich   fixed linux compa...
329
330
  		std::string headername = outname + ".hdr";              //the header file name
  
63fc1df8   David Mayerich   added an inverse ...
331
332
333
334
335
  		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...
336
  			if(mask != NULL && !mask[n])						//if the normalization band is below threshold
63fc1df8   David Mayerich   added an inverse ...
337
338
  				memset(s, 0, sizeof(T) * B);					//set the output to zero
  			else{
9d3ba0b1   David Mayerich   added stim::hsi a...
339
  				pixel(s, n);										//retrieve the spectrum s
1a224b6a   David Mayerich   fixed linux compa...
340
341
  				nv = hsi<T>::interp_spectrum(s, w);							//find the value of the normalization band
  
9d3ba0b1   David Mayerich   added stim::hsi a...
342
  				for(unsigned long long b = 0; b < B; b++)			//for each band in the spectrum
63fc1df8   David Mayerich   added an inverse ...
343
  					s[b] /= nv;										//divide by the normalization value
88c3e636   heziqi   Ziqi added big.h
344
  			}
1a224b6a   David Mayerich   fixed linux compa...
345
  
63fc1df8   David Mayerich   added an inverse ...
346
  			if(PROGRESS) progress = (double)(n+1) / N * 100;	//set the current progress
6708cc25   heziqi   Ziqi added envi c...
347
  
63fc1df8   David Mayerich   added an inverse ...
348
349
  			target.write((char*)s, sizeof(T) * B);	//write the corrected data into destination
  		}														//end for each pixel
1a224b6a   David Mayerich   fixed linux compa...
350
  
63fc1df8   David Mayerich   added an inverse ...
351
  		free(s);
88c3e636   heziqi   Ziqi added big.h
352
  		target.close();
88c3e636   heziqi   Ziqi added big.h
353
354
  		return true;
  	}
c30dd3e6   David Mayerich   added vector norm...
355
356
  
  	void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){}
1a224b6a   David Mayerich   fixed linux compa...
357
  
f6169dea   heziqi   Ziqi completed bi...
358
  
a23c4132   David Mayerich   Doxygen comments ...
359
360
361
  	/// 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 ...
362
  	bool bil(std::string outname, bool PROGRESS = false)
11f177d5   heziqi   Ziqi added functi...
363
  	{
9d3ba0b1   David Mayerich   added stim::hsi a...
364
  		unsigned long long S = X() * Z() * sizeof(T);		//calculate the number of bytes in a ZX slice
1a224b6a   David Mayerich   fixed linux compa...
365
  
11f177d5   heziqi   Ziqi added functi...
366
  		std::ofstream target(outname.c_str(), std::ios::binary);
9d3ba0b1   David Mayerich   added stim::hsi a...
367
  		//std::string headername = outname + ".hdr";
1a224b6a   David Mayerich   fixed linux compa...
368
  
f6169dea   heziqi   Ziqi completed bi...
369
  		T * p;			//pointer to the current ZX slice for bip file
11f177d5   heziqi   Ziqi added functi...
370
  		p = (T*)malloc(S);
f6169dea   heziqi   Ziqi completed bi...
371
372
  		T * q;			//pointer to the current XZ slice for bil file
  		q = (T*)malloc(S);
1a224b6a   David Mayerich   fixed linux compa...
373
  
9d3ba0b1   David Mayerich   added stim::hsi a...
374
  		for ( unsigned long long i = 0; i < Y(); i++)
1a224b6a   David Mayerich   fixed linux compa...
375
  		{
724ec347   David Mayerich   simplified the EN...
376
  			read_plane_y(p, i);
9d3ba0b1   David Mayerich   added stim::hsi a...
377
  			for ( unsigned long long k = 0; k < Z(); k++)
f6169dea   heziqi   Ziqi completed bi...
378
  			{
9d3ba0b1   David Mayerich   added stim::hsi a...
379
380
  				unsigned long long ks = k * X();
  				for ( unsigned long long j = 0; j < X(); j++)
724ec347   David Mayerich   simplified the EN...
381
  					q[ks + j] = p[k + j * Z()];
cbce396e   David Mayerich   added support for...
382
  
63fc1df8   David Mayerich   added an inverse ...
383
  				if(PROGRESS) progress = (double)(i * Z() + k+1) / (Y() * Z()) * 100;
f6169dea   heziqi   Ziqi completed bi...
384
  			}
1a224b6a   David Mayerich   fixed linux compa...
385
  			target.write(reinterpret_cast<const char*>(q), S);   //write a band data into target file
11f177d5   heziqi   Ziqi added functi...
386
  		}
f6169dea   heziqi   Ziqi completed bi...
387
  
11f177d5   heziqi   Ziqi added functi...
388
  		free(p);
f6169dea   heziqi   Ziqi completed bi...
389
  		free(q);
11f177d5   heziqi   Ziqi added functi...
390
391
392
  		target.close();
  		return true;
  	}
88c3e636   heziqi   Ziqi added big.h
393
  
a23c4132   David Mayerich   Doxygen comments ...
394
395
396
397
398
399
400
401
  	/// 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...
402
403
  	bool baseline_band(double lb, double rb, T* lp, T* rp, double wavelength, T* result){
  
9d3ba0b1   David Mayerich   added stim::hsi a...
404
  		unsigned long long XY = X() * Y();
20c212c0   heziqi   Ziqi added functi...
405
406
407
408
  		band(result, wavelength);		//get band
  
  		//perform the baseline correction
  		double r = (double) (wavelength - lb) / (double) (rb - lb);
9d3ba0b1   David Mayerich   added stim::hsi a...
409
  		for(unsigned long long i=0; i < XY; i++){
20c212c0   heziqi   Ziqi added functi...
410
411
412
413
  			result[i] =(T) (result[i] - (rp[i] - lp[i]) * r - lp[i] );
  		}
  		return true;
  	}
1a224b6a   David Mayerich   fixed linux compa...
414
  
a23c4132   David Mayerich   Doxygen comments ...
415
416
417
418
419
420
  	/// 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...
421
  	bool height(double lb, double rb, double bandwavelength, T* result){
20c212c0   heziqi   Ziqi added functi...
422
423
424
  
  		T* lp;
  		T* rp;
9d3ba0b1   David Mayerich   added stim::hsi a...
425
426
  		unsigned long long XY = X() * Y();
  		unsigned long long S = XY * sizeof(T);
70407ea9   heziqi   Ziqi completed he...
427
428
  		lp = (T*) malloc(S);			//memory allocation
  		rp = (T*) malloc(S);
20c212c0   heziqi   Ziqi added functi...
429
  
70407ea9   heziqi   Ziqi completed he...
430
  		band(lp, lb);
1a224b6a   David Mayerich   fixed linux compa...
431
  		band(rp, rb);
20c212c0   heziqi   Ziqi added functi...
432
433
434
  
  		baseline_band(lb, rb, lp, rp, bandwavelength, result);
  
517876d6   heziqi   metrics finished ...
435
436
  		free(lp);
  		free(rp);
20c212c0   heziqi   Ziqi added functi...
437
438
439
  		return true;
  	}
  
c359422d   heziqi   Ziqi added functi...
440
  
a23c4132   David Mayerich   Doxygen comments ...
441
442
443
444
445
446
447
  	/// 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...
448
  	bool area(double lb, double rb, double lab, double rab, T* result){
20c212c0   heziqi   Ziqi added functi...
449
450
451
452
453
  
  		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...
454
  
9d3ba0b1   David Mayerich   added stim::hsi a...
455
456
  		unsigned long long XY = X() * Y();
  		unsigned long long S = XY * sizeof(T);
20c212c0   heziqi   Ziqi added functi...
457
458
459
460
461
  
  		lp = (T*) malloc(S);			//memory allocation
  		rp = (T*) malloc(S);
  		cur = (T*) malloc(S);
  		cur2 = (T*) malloc(S);
20c212c0   heziqi   Ziqi added functi...
462
463
464
465
  
  		memset(result, (char)0, S);
  
  		//find the wavelenght position in the whole band
9d3ba0b1   David Mayerich   added stim::hsi a...
466
467
468
  		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...
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
  
  
  
  		//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...
484
  		while (lab >= w[ai]){
20c212c0   heziqi   Ziqi added functi...
485
486
  			ai++;
  		}
70407ea9   heziqi   Ziqi completed he...
487
  		while (rab <= w[bi]){
20c212c0   heziqi   Ziqi added functi...
488
489
490
491
492
493
494
  			bi--;
  		}
  
  		band(lp, lb);
  		band(rp, rb);
  
  		//calculate the beginning and the ending part
70407ea9   heziqi   Ziqi completed he...
495
496
  		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...
497
498
  		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...
499
  		}
70407ea9   heziqi   Ziqi completed he...
500
501
  		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...
502
  		for(unsigned long long j = 0; j < XY; j++){
9d3ba0b1   David Mayerich   added stim::hsi a...
503
  			result[j] += (T)((w[ai] - lab) * ((double)cur[j] + (double)cur2[j]) / 2.0);
20c212c0   heziqi   Ziqi added functi...
504
  		}
517876d6   heziqi   metrics finished ...
505
  
20c212c0   heziqi   Ziqi added functi...
506
507
  		//calculate the area
  		ai++;
9d3ba0b1   David Mayerich   added stim::hsi a...
508
  		for(unsigned long long i = ai; i <= bi ;i++)
20c212c0   heziqi   Ziqi added functi...
509
510
  		{
  			baseline_band(lb, rb, lp, rp, w[ai], cur2);
9d3ba0b1   David Mayerich   added stim::hsi a...
511
  			for(unsigned long long j = 0; j < XY; j++)
20c212c0   heziqi   Ziqi added functi...
512
  			{
1a224b6a   David Mayerich   fixed linux compa...
513
  				result[j] += (T)((w[ai] - w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
20c212c0   heziqi   Ziqi added functi...
514
515
516
  			}
  			std::swap(cur,cur2);		//swap the band pointers
  		}
70407ea9   heziqi   Ziqi completed he...
517
  
517876d6   heziqi   metrics finished ...
518
519
520
521
  		free(lp);
  		free(rp);
  		free(cur);
  		free(cur2);
20c212c0   heziqi   Ziqi added functi...
522
523
524
  		return true;
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
525
526
527
528
529
530
531
532
533
  	/// 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...
534
  	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...
535
  
1a224b6a   David Mayerich   fixed linux compa...
536
537
538
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
70407ea9   heziqi   Ziqi completed he...
539
540
541
542
  		//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...
543
  		for(unsigned long long i = 0; i < X() * Y(); i++){
70407ea9   heziqi   Ziqi completed he...
544
545
546
547
548
  			if(p1[i] == 0 && p2[i] ==0)
  				result[i] = 1;
  			else
  				result[i] = p1[i] / p2[i];
  		}
517876d6   heziqi   metrics finished ...
549
550
551
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
552
  		return true;
70407ea9   heziqi   Ziqi completed he...
553
  	}
1a224b6a   David Mayerich   fixed linux compa...
554
  
a23c4132   David Mayerich   Doxygen comments ...
555
556
557
558
559
560
561
562
563
  	/// 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...
564
565
  	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...
566
567
568
569
  
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
70407ea9   heziqi   Ziqi completed he...
570
571
572
573
  		//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...
574
  		for(unsigned long long i = 0; i < X() * Y(); i++){
70407ea9   heziqi   Ziqi completed he...
575
576
577
578
579
  			if(p1[i] == 0 && p2[i] ==0)
  				result[i] = 1;
  			else
  				result[i] = p1[i] / p2[i];
  		}
517876d6   heziqi   metrics finished ...
580
581
582
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
583
584
585
  		return true;
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
586
587
588
589
590
591
592
593
594
  	/// 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...
595
  	/// @param rab2 is the label value for the right bound (end of the integration) of the second peak (denominator)
a23c4132   David Mayerich   Doxygen comments ...
596
  	/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
ba51ae6a   David Mayerich   fixed metric calc...
597
598
  	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...
599
600
601
602
  
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
70407ea9   heziqi   Ziqi completed he...
603
604
605
606
  		//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...
607
  		for(unsigned long long i = 0; i < X() * Y(); i++){
70407ea9   heziqi   Ziqi completed he...
608
609
610
611
612
  			if(p1[i] == 0 && p2[i] ==0)
  				result[i] = 1;
  			else
  				result[i] = p1[i] / p2[i];
  		}
517876d6   heziqi   metrics finished ...
613
614
615
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
616
617
  		return true;
  	}
517876d6   heziqi   metrics finished ...
618
  
a23c4132   David Mayerich   Doxygen comments ...
619
620
621
622
623
624
625
  	/// 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 ...
626
627
628
629
630
631
  	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...
632
633
  		unsigned long long XY = X() * Y();
  		unsigned long long S = XY * sizeof(T);
517876d6   heziqi   metrics finished ...
634
635
636
637
638
639
640
641
642
  
  		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...
643
644
645
  		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 ...
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
  
  		//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...
672
673
  		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 ...
674
675
676
  		}
  		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...
677
  		for(unsigned long long j = 0; j < XY; j++){
9d3ba0b1   David Mayerich   added stim::hsi a...
678
  			result[j] += (T)((w[ai] - lab) * (w[ai] + lab) * ((double)cur[j] + (double)cur2[j]) / 4.0);
517876d6   heziqi   metrics finished ...
679
680
681
682
  		}
  
  		//calculate f(x) times x
  		ai++;
9d3ba0b1   David Mayerich   added stim::hsi a...
683
  		for(unsigned long long i = ai; i <= bi ;i++)
517876d6   heziqi   metrics finished ...
684
685
  		{
  			baseline_band(lb, rb, lp, rp, w[ai], cur2);
9d3ba0b1   David Mayerich   added stim::hsi a...
686
  			for(unsigned long long j = 0; j < XY; j++)
517876d6   heziqi   metrics finished ...
687
  			{
1a224b6a   David Mayerich   fixed linux compa...
688
  				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 ...
689
690
691
692
693
694
695
696
697
698
699
  			}
  			std::swap(cur,cur2);		//swap the band pointers
  		}
  
  		free(lp);
  		free(rp);
  		free(cur);
  		free(cur2);
  		return true;
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
700
701
702
703
704
705
706
  	/// 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...
707
  	bool centroid(T* result, double lb, double rb, double lab, double rab, unsigned char* mask = NULL){
1a224b6a   David Mayerich   fixed linux compa...
708
709
710
  		T* p1 = (T*)malloc(X() * Y() * sizeof(T));
  		T* p2 = (T*)malloc(X() * Y() * sizeof(T));
  
517876d6   heziqi   metrics finished ...
711
712
713
714
  		//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...
715
  		for(unsigned long long i = 0; i < X() * Y(); i++){
ba51ae6a   David Mayerich   fixed metric calc...
716
  			if(mask == NULL || mask[i])
517876d6   heziqi   metrics finished ...
717
718
719
720
721
  				result[i] = p1[i] / p2[i];
  		}
  
  		free(p1);
  		free(p2);
1a224b6a   David Mayerich   fixed linux compa...
722
  		return true;
517876d6   heziqi   metrics finished ...
723
724
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
725
726
727
728
729
730
731
  	/// 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
ba51ae6a   David Mayerich   fixed metric calc...
732
  	bool build_mask(unsigned char* mask, double mask_band, double threshold, bool PROGRESS = false){
4a6f666c   heziqi   Added mask method
733
  
724ec347   David Mayerich   simplified the EN...
734
  		T* temp = (T*)malloc(X() * Y() * sizeof(T));		//allocate memory for the certain band
798cea97   David Mayerich   added progress ba...
735
  		band(temp, mask_band, PROGRESS);
4a6f666c   heziqi   Added mask method
736
  
9d3ba0b1   David Mayerich   added stim::hsi a...
737
  		for (unsigned long long i = 0; i < X() * Y();i++) {
798cea97   David Mayerich   added progress ba...
738
  			if (temp[i] < threshold)
ba51ae6a   David Mayerich   fixed metric calc...
739
  				mask[i] = 0;
798cea97   David Mayerich   added progress ba...
740
  			else
ba51ae6a   David Mayerich   fixed metric calc...
741
  				mask[i] = 255;
4a6f666c   heziqi   Added mask method
742
743
  		}
  
517876d6   heziqi   metrics finished ...
744
  		free(temp);
4a6f666c   heziqi   Added mask method
745
746
747
748
  		return true;
  
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
749
750
751
752
  	/// 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...
753
  	bool apply_mask(std::string outfile, unsigned char* p, bool PROGRESS = false){
740f8cd2   heziqi   added apply_mask
754
755
756
  
  		std::ofstream target(outfile.c_str(), std::ios::binary);
  
9d3ba0b1   David Mayerich   added stim::hsi a...
757
758
  		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
759
  
1a55a328   David Mayerich   Added comments fo...
760
  		T * temp = (T*)malloc(L);		//allocate space for that page
740f8cd2   heziqi   added apply_mask
761
  
9d3ba0b1   David Mayerich   added stim::hsi a...
762
  		for (unsigned long long i = 0; i < Y(); i++)			//for each page (Y in BIP)
740f8cd2   heziqi   added apply_mask
763
  		{
724ec347   David Mayerich   simplified the EN...
764
  			read_plane_y(temp, i);							//load that page (it's pointed to by temp)
9d3ba0b1   David Mayerich   added stim::hsi a...
765
  			for ( unsigned long long j = 0; j < X(); j++)	//for each X value
740f8cd2   heziqi   added apply_mask
766
  			{
9d3ba0b1   David Mayerich   added stim::hsi a...
767
  				for (unsigned long long k = 0; k < Z(); k++)	//for each B value (band)
740f8cd2   heziqi   added apply_mask
768
  				{
724ec347   David Mayerich   simplified the EN...
769
770
  					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...
771
  				else								//otherwise just continue
740f8cd2   heziqi   added apply_mask
772
773
774
  					continue;
  				}
  			}
1a55a328   David Mayerich   Added comments fo...
775
  			target.write(reinterpret_cast<const char*>(temp), L);   //write the edited band data into target file
6a46c8ff   David Mayerich   fixed bugs in app...
776
  			if(PROGRESS) progress = (double)(i+1) / (double)Y() * 100;
740f8cd2   heziqi   added apply_mask
777
  		}
1a55a328   David Mayerich   Added comments fo...
778
779
780
  		target.close();						//close the target file
  		free(temp);							//free allocated memory
  		return true;						//return true
740f8cd2   heziqi   added apply_mask
781
782
  	}
  
ba51ae6a   David Mayerich   fixed metric calc...
783
784
785
786
787
788
789
790
791
792
793
  	/// Copies all spectra corresponding to nonzero values of a mask into a pre-allocated matrix of size (B x P)
  	///		where P is the number of masked pixels and B is the number of bands. The allocated memory can be accessed
  	///		using the following indexing: i = p*B + b
  	/// @param matrix is the destination for the pixel data
  	/// @param mask is the mask
  	bool sift(T* matrix, unsigned char* mask){
  		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...
794
  
ba51ae6a   David Mayerich   fixed metric calc...
795
796
797
798
799
800
801
802
803
804
805
  		size_t p = 0;										//create counter variables
  		for(size_t xy = 0; xy < XY; xy++){					//for each pixel
  			if(mask[xy]){									//if the current pixel is masked
  				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...
806
  				file.seekg(Bbytes, std::ios::cur);			//otherwise skip this band
ba51ae6a   David Mayerich   fixed metric calc...
807
808
809
  		}
  		return true;
  	}
e843658b   Brad Deutsch   Previous push did...
810
  
2fcab4f0   David Mayerich   added unsifting f...
811
  	/// Saves to disk only those spectra corresponding to mask values != 0
c1078e19   David Mayerich   HSIproc bug fixes...
812
  	bool sift(std::string outfile, unsigned char* mask, bool PROGRESS = false){
8550e94f   David Mayerich   added sifting and...
813
814
  
  		//reset the file pointer to the beginning of the file
1a224b6a   David Mayerich   fixed linux compa...
815
  		file.seekg(0, std::ios::beg);
8550e94f   David Mayerich   added sifting and...
816
817
  
  		// open an output stream
e843658b   Brad Deutsch   Previous push did...
818
819
  		std::ofstream target(outfile.c_str(), std::ios::binary);
  
8550e94f   David Mayerich   added sifting and...
820
  		//allocate space for a single spectrum
9d3ba0b1   David Mayerich   added stim::hsi a...
821
  		unsigned long long B = Z();
8550e94f   David Mayerich   added sifting and...
822
  		T* spectrum = (T*) malloc(B * sizeof(T));
e843658b   Brad Deutsch   Previous push did...
823
  
8550e94f   David Mayerich   added sifting and...
824
  		//calculate the number of pixels in a band
9d3ba0b1   David Mayerich   added stim::hsi a...
825
  		unsigned long long XY = X() * Y();
e843658b   Brad Deutsch   Previous push did...
826
  
8550e94f   David Mayerich   added sifting and...
827
  		//for each pixel
9d3ba0b1   David Mayerich   added stim::hsi a...
828
829
  		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...
830
831
832
833
834
  
  			//if the current pixel isn't masked
  			if( mask[x] == 0){
  				//increment the number of skipped pixels
  				skip++;
e843658b   Brad Deutsch   Previous push did...
835
  			}
8550e94f   David Mayerich   added sifting and...
836
837
838
839
840
841
842
843
844
845
846
847
848
849
  			//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...
850
  			}
9d3ba0b1   David Mayerich   added stim::hsi a...
851
  			if(PROGRESS) progress = (double) (x+1) / XY * 100;
8550e94f   David Mayerich   added sifting and...
852
  
e843658b   Brad Deutsch   Previous push did...
853
  		}
8550e94f   David Mayerich   added sifting and...
854
855
  
  		//close the output file
e843658b   Brad Deutsch   Previous push did...
856
  		target.close();
8550e94f   David Mayerich   added sifting and...
857
  		free(spectrum);
cbce396e   David Mayerich   added support for...
858
  
cbce396e   David Mayerich   added support for...
859
  		return true;
8550e94f   David Mayerich   added sifting and...
860
861
  	}
  
a2bf1d08   David Mayerich   general bug fixes...
862
  	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...
863
864
865
866
867
  
  		// 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...
868
  		file.seekg(0, std::ios::beg);
8550e94f   David Mayerich   added sifting and...
869
870
  
  		//allocate space for a single spectrum
9d3ba0b1   David Mayerich   added stim::hsi a...
871
  		unsigned long long B = Z();
8550e94f   David Mayerich   added sifting and...
872
873
874
875
876
877
878
  		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...
879
  		unsigned long long XY = samples * lines;
8550e94f   David Mayerich   added sifting and...
880
881
  
  		//for each pixel
9d3ba0b1   David Mayerich   added stim::hsi a...
882
883
  		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...
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
  
  			//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...
901
  			if(PROGRESS) progress = (double)(x + 1) / XY * 100;
3150c537   David Mayerich   added further pro...
902
  
8550e94f   David Mayerich   added sifting and...
903
904
905
906
907
908
  		}
  
  		//close the output file
  		target.close();
  		free(spectrum);
  
a2bf1d08   David Mayerich   general bug fixes...
909
  		//progress = 100;
3150c537   David Mayerich   added further pro...
910
  
e843658b   Brad Deutsch   Previous push did...
911
  		return true;
8550e94f   David Mayerich   added sifting and...
912
913
  
  
e843658b   Brad Deutsch   Previous push did...
914
915
  	}
  
a23c4132   David Mayerich   Doxygen comments ...
916
917
918
919
  	/// 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
63fc1df8   David Mayerich   added an inverse ...
920
  	bool avg_band(double* p, unsigned char* mask = NULL, bool PROGRESS = false){
9d3ba0b1   David Mayerich   added stim::hsi a...
921
922
923
924
925
926
927
928
929
930
931
932
933
934
  		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
  		memset(p, 0, sizeof(double) * Z());							//initialize the average spectrum to zero (0)
  		//for (unsigned j = 0; j < Z(); j++){
  		//	p[j] = 0;
  		//}
  
  		unsigned long long count = nnz(mask);									//calculate the number of masked pixels
  
  		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
  					p[j] += (double)temp[j] / (double)count;		//add the weighted value to the average
a43c4fe1   heziqi   Added crop in env...
935
936
  				}
  			}
9d3ba0b1   David Mayerich   added stim::hsi a...
937
  			if(PROGRESS) progress = (double)(i+1) / XY * 100;		//increment the progress
a43c4fe1   heziqi   Added crop in env...
938
  		}
2aa68315   David Mayerich   major bug fixes f...
939
  
a43c4fe1   heziqi   Added crop in env...
940
  		free(temp);
a43c4fe1   heziqi   Added crop in env...
941
942
  		return true;
  	}
c6f526c2   David Mayerich   added cuBLAS supp...
943
944
  #ifdef CUDA_FOUND
  	/// Calculate the covariance matrix for masked pixels using cuBLAS
9d3ba0b1   David Mayerich   added stim::hsi a...
945
  	/// Note that cuBLAS only supports integer-sized arrays, so there may be issues with large spectra
63fc1df8   David Mayerich   added an inverse ...
946
  	bool co_matrix_cublas(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
c6f526c2   David Mayerich   added cuBLAS supp...
947
  
1a224b6a   David Mayerich   fixed linux compa...
948
  		cudaError_t cudaStat;
c6f526c2   David Mayerich   added cuBLAS supp...
949
950
  		cublasStatus_t stat;
  		cublasHandle_t handle;
1a224b6a   David Mayerich   fixed linux compa...
951
  
63fc1df8   David Mayerich   added an inverse ...
952
  		progress = 0;													//initialize the progress to zero (0)
c6f526c2   David Mayerich   added cuBLAS supp...
953
954
955
956
957
958
959
960
961
962
963
  		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...
964
  		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...
965
  
c6f526c2   David Mayerich   added cuBLAS supp...
966
967
968
969
970
971
972
973
974
  		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
  		if (stat != CUBLAS_STATUS_SUCCESS) {								//test the cuBLAS instance to make sure it is valid
  			printf ("CUBLAS initialization failed\n");
  			return EXIT_FAILURE;
  		}
  		for (unsigned long long xy = 0; xy < XY; xy++){										//for each pixel
5a5359da   David Mayerich   fixed - allows th...
975
976
  			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...
977
978
979
  				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...
980
  			}
63fc1df8   David Mayerich   added an inverse ...
981
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;													//record the current progress
5a5359da   David Mayerich   fixed - allows th...
982
  
c6f526c2   David Mayerich   added cuBLAS supp...
983
984
  		}
  
9d3ba0b1   David Mayerich   added stim::hsi a...
985
  		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...
986
987
988
989
990
  
  		cudaFree(A_dev);														//clean up allocated device memory
  		cudaFree(s_dev);
  		cudaFree(avg_dev);
  
9d3ba0b1   David Mayerich   added stim::hsi a...
991
992
  		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...
993
994
995
996
  				co[B * i + j] = co[B * j + i];
  			}
  		}
  
c6f526c2   David Mayerich   added cuBLAS supp...
997
998
999
  		return true;
  	}
  #endif
1a224b6a   David Mayerich   fixed linux compa...
1000
  
2aa68315   David Mayerich   major bug fixes f...
1001
  	/// Calculate the covariance matrix for all masked pixels in the image with 64-bit floating point precision.
a23c4132   David Mayerich   Doxygen comments ...
1002
1003
1004
1005
  
  	/// @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
63fc1df8   David Mayerich   added an inverse ...
1006
  	bool co_matrix(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
c6f526c2   David Mayerich   added cuBLAS supp...
1007
1008
  
  #ifdef CUDA_FOUND
b37665a9   David Mayerich   cuBLAS won't run ...
1009
1010
1011
1012
1013
1014
  		int dev_count;
  		cudaGetDeviceCount(&dev_count);									//get the number of CUDA devices
  		cudaDeviceProp prop;
  		cudaGetDeviceProperties(&prop, 0);								//get the property of the first device
  		if(dev_count > 0 && prop.major != 9999)							//if the first device is not an emulator
  			return co_matrix_cublas(co, avg, mask, PROGRESS);			//use cuBLAS to calculate the covariance matrix
c6f526c2   David Mayerich   added cuBLAS supp...
1015
  #endif
63fc1df8   David Mayerich   added an inverse ...
1016
  		progress = 0;
a43c4fe1   heziqi   Added crop in env...
1017
  		//memory allocation
2aa68315   David Mayerich   major bug fixes f...
1018
1019
  		unsigned long long XY = X() * Y();
  		unsigned long long B = Z();
a43c4fe1   heziqi   Added crop in env...
1020
  		T* temp = (T*)malloc(sizeof(T) * B);
9d3ba0b1   David Mayerich   added stim::hsi a...
1021
1022
  
  		unsigned long long count = nnz(mask);								//count the number of masked pixels
63fc1df8   David Mayerich   added an inverse ...
1023
  
2aa68315   David Mayerich   major bug fixes f...
1024
1025
1026
1027
  		//initialize covariance matrix
  		memset(co, 0, B * B * sizeof(double));
  
  		//calculate covariance matrix
2aa68315   David Mayerich   major bug fixes f...
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
  		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...
1041
1042
  				}
  			}
63fc1df8   David Mayerich   added an inverse ...
1043
  			if(PROGRESS) progress = (double)(xy+1) / XY * 100;
a43c4fe1   heziqi   Added crop in env...
1044
  		}
2aa68315   David Mayerich   major bug fixes f...
1045
1046
1047
1048
  		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...
1049
1050
1051
  			}
  		}
  
a43c4fe1   heziqi   Added crop in env...
1052
  		free(temp);
2aa68315   David Mayerich   major bug fixes f...
1053
  		free(temp_precise);
a43c4fe1   heziqi   Added crop in env...
1054
1055
1056
  		return true;
  	}
  
9d3ba0b1   David Mayerich   added stim::hsi a...
1057
  	bool project(std::string outfile, double* center, double* basis, unsigned long long M, unsigned char* mask = NULL, bool PROGRESS = false){
63fc1df8   David Mayerich   added an inverse ...
1058
1059
  
  		std::ofstream target(outfile.c_str(), std::ios::binary);	//open the target binary file
9d3ba0b1   David Mayerich   added stim::hsi a...
1060
  		//std::string headername = outfile + ".hdr";					//the header file name
63fc1df8   David Mayerich   added an inverse ...
1061
1062
1063
1064
1065
1066
1067
1068
1069
  
  		//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...
1070
1071
1072
1073
1074
1075
1076
1077
  			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 ...
1078
1079
  				}
  			}
1a224b6a   David Mayerich   fixed linux compa...
1080
  
63fc1df8   David Mayerich   added an inverse ...
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
  			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
  
  		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...
1110
  				for(unsigned long long b = 0; b < B; b++){			//for each component of the basis vector
9d3ba0b1   David Mayerich   added stim::hsi a...
1111
  					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 ...
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
  				}
  			}
  
  			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
1126
  
a23c4132   David Mayerich   Doxygen comments ...
1127
1128
1129
1130
1131
1132
1133
  	/// 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...
1134
1135
1136
1137
1138
  	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 ...
1139
1140
  								   unsigned long long b1,
  								   bool PROGRESS = false){
a43c4fe1   heziqi   Added crop in env...
1141
  
3033b886   David Mayerich   implemented spect...
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
  		//calculate the new number of samples, lines, and bands
  		unsigned long long samples = x1 - x0;
  		unsigned long long lines = y1 - y0;
  		unsigned long long bands = b1 - b0;
  
  		//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
a43c4fe1   heziqi   Added crop in env...
1153
  		T* temp = (T*)malloc(L);
3033b886   David Mayerich   implemented spect...
1154
1155
  
  		//open an output file for binary writing
a43c4fe1   heziqi   Added crop in env...
1156
  		std::ofstream out(outfile.c_str(), std::ios::binary);
3033b886   David Mayerich   implemented spect...
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
  
  		//seek to the first pixel in the cropped image
  		file.seekg( (y0 * X() * Z() + x0 * Z() + b0) * sizeof(T), std::ios::beg);
  
  		//distance between sample spectra in the same line
  		unsigned long long jump_sample = ( (Z() - b1) + b0 ) * sizeof(T);
  
  		//distance between sample spectra in adjacent lines
  		unsigned long long jump_line = (X() - x1) * Z() * sizeof(T);
  
  
  		//unsigned long long sp = y0 * X() + x0;		//start pixel
  
  		//for each pixel in the image
  		for (unsigned y = 0; y < lines; y++)
a43c4fe1   heziqi   Added crop in env...
1172
  		{
3033b886   David Mayerich   implemented spect...
1173
  			for (unsigned x = 0; x < samples; x++)
a43c4fe1   heziqi   Added crop in env...
1174
  			{
3033b886   David Mayerich   implemented spect...
1175
1176
1177
  				//read the cropped spectral region
  				file.read( (char*) temp, L );
  				//pixel(temp, sp + x + y * X());
1a224b6a   David Mayerich   fixed linux compa...
1178
  				out.write(reinterpret_cast<const char*>(temp), L);   //write slice data into target file
3150c537   David Mayerich   added further pro...
1179
  
3033b886   David Mayerich   implemented spect...
1180
1181
  				file.seekg(jump_sample, std::ios::cur);
  
63fc1df8   David Mayerich   added an inverse ...
1182
  				if(PROGRESS) progress = (double)((y+1) * samples + x + 1) / (lines * samples) * 100;
a43c4fe1   heziqi   Added crop in env...
1183
  			}
3033b886   David Mayerich   implemented spect...
1184
1185
  
  			file.seekg(jump_line, std::ios::cur);
a43c4fe1   heziqi   Added crop in env...
1186
1187
  		}
  		free(temp);
3150c537   David Mayerich   added further pro...
1188
  
a43c4fe1   heziqi   Added crop in env...
1189
1190
1191
  		return true;
  	}
  
dc8eb8aa   David Mayerich   added band trimmi...
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
  	/// 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...
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
  	/// 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...
1259
  			if(PROGRESS) progress = (double)( (y+1) * S[0] + 1) / (double) (S[0] * S[1]) * 100;
2ce6954b   David Mayerich   added the ability...
1260
1261
1262
  		}
  	}
  
9b2a5d71   David Mayerich   implemented finit...
1263
1264
1265
1266
1267
1268
1269
1270
  	/// 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...
1271
  	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...
1272
  		std::ofstream out(outfile.c_str(), std::ios::binary);		//open the output file for writing
1a224b6a   David Mayerich   fixed linux compa...
1273
  
9b2a5d71   David Mayerich   implemented finit...
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
  		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...
1289
1290
1291
1292
1293
  			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...
1294
1295
1296
1297
1298
1299
1300
  				}
  			}
  			out.write((char*)outspec, Nb);							//output the band
  			if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
  		}
  	}
  
9d77bbd9   David Mayerich   updates for HSIvi...
1301
  	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...
1302
  		std::ofstream out(outfile.c_str(), std::ios::binary);		//open the output file for writing
1a224b6a   David Mayerich   fixed linux compa...
1303
1304
  
  
9b2a5d71   David Mayerich   implemented finit...
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
  		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...
1326
1327
1328
  			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...
1329
1330
  					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...
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
  					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...
1351
  				}
9b2a5d71   David Mayerich   implemented finit...
1352
  			}
9b2a5d71   David Mayerich   implemented finit...
1353
1354
1355
1356
1357
  			out.write((char*)outspec, Bb);							//output the band
  			if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
  		}
  	}
  
dc8eb8aa   David Mayerich   added band trimmi...
1358
  
0df38ff3   heziqi   fixed head detached
1359
  
a23c4132   David Mayerich   Doxygen comments ...
1360
  	/// Close the file.
f6169dea   heziqi   Ziqi completed bi...
1361
1362
1363
1364
1365
  	bool close(){
  		file.close();
  		return true;
  	}
  
88c3e636   heziqi   Ziqi added big.h
1366
  	};
6aa04ba2   David Mayerich   interleave types
1367
  }
e8eb202f   David Mayerich   added a new ENVI ...
1368
1369
  
  #endif