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