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
|
}
|
735a2a24
David Mayerich
started testing o...
|
146
147
148
149
150
|
bool spectrum(T* p, size_t n, bool PROGRESS = false){
size_t y = n / X();
size_t x = n - y * X();
return read_line_0(p, x, y, PROGRESS); //read a line in the binary YZ plane (dimension order for BIP is ZXY)
}
|
724ec347
David Mayerich
simplified the EN...
|
151
152
153
154
155
156
157
|
/// 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...
|
158
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
159
|
unsigned long long B = Z();
|
11f177d5
heziqi
Ziqi added functi...
|
160
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
161
|
unsigned long long page=0; //samples around the wavelength
|
517876d6
heziqi
metrics finished ...
|
162
|
|
11f177d5
heziqi
Ziqi added functi...
|
163
164
165
166
|
//get the bands numbers around the wavelength
//if wavelength is smaller than the first one in header file
|
6708cc25
heziqi
Ziqi added envi c...
|
167
|
if ( w[page] > wavelength ){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
168
|
for(unsigned long long j = 0; j < X(); j++)
|
11f177d5
heziqi
Ziqi added functi...
|
169
170
|
{
p[j] = c[j * B];
|
1a224b6a
David Mayerich
fixed linux compa...
|
171
|
}
|
11f177d5
heziqi
Ziqi added functi...
|
172
173
174
|
return true;
}
|
6708cc25
heziqi
Ziqi added envi c...
|
175
|
while( w[page] < wavelength )
|
11f177d5
heziqi
Ziqi added functi...
|
176
177
178
179
|
{
page++;
//if wavelength is larger than the last wavelength in header file
if (page == B) {
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
180
|
for(unsigned long long j = 0; j < X(); j++)
|
11f177d5
heziqi
Ziqi added functi...
|
181
182
183
184
185
186
|
{
p[j] = c[(j + 1) * B - 1];
}
return true;
}
}
|
6708cc25
heziqi
Ziqi added envi c...
|
187
|
if ( wavelength < w[page] )
|
11f177d5
heziqi
Ziqi added functi...
|
188
|
{
|
517876d6
heziqi
metrics finished ...
|
189
190
|
T * p1;
T * p2;
|
724ec347
David Mayerich
simplified the EN...
|
191
192
|
p1=(T*)malloc( X() * sizeof(T)); //memory allocation
p2=(T*)malloc( X() * sizeof(T));
|
11f177d5
heziqi
Ziqi added functi...
|
193
|
//band_index(p1, page - 1);
|
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
198
|
{
p1[j] = c[j * B + page - 1];
}
//band_index(p2, page );
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
199
|
for(unsigned long long j = 0; j < X(); j++)
|
11f177d5
heziqi
Ziqi added functi...
|
200
201
202
|
{
p2[j] = c[j * B + page];
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
203
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
204
|
for(unsigned long long i=0; i < X(); i++){
|
6708cc25
heziqi
Ziqi added envi c...
|
205
|
double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
|
11f177d5
heziqi
Ziqi added functi...
|
206
207
|
p[i] = (p2[i] - p1[i]) * r + p1[i];
}
|
517876d6
heziqi
metrics finished ...
|
208
209
|
free(p1);
free(p2);
|
11f177d5
heziqi
Ziqi added functi...
|
210
211
212
213
|
}
else //if the wavelength is equal to a wavelength in header file
{
//band_index(p, page);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
214
|
for(unsigned long long j = 0; j < X(); j++)
|
11f177d5
heziqi
Ziqi added functi...
|
215
216
217
218
219
|
{
p[j] = c[j * B + page];
}
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
220
|
return true;
|
11f177d5
heziqi
Ziqi added functi...
|
221
|
}
|
bfe9c6db
heziqi
added feature_mat...
|
222
|
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
223
|
/// Retrieve a single pixel and store it in a pre-allocated double array.
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
224
225
|
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...
|
226
227
228
229
|
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...
|
230
|
unsigned long long B = Z();
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
231
232
233
234
235
|
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...
|
236
|
for(unsigned long long i = 0; i < B; i++) //for each element of the spectrum
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
237
|
p[i] = (double) temp[i]; //cast each element to a double value
|
193bb525
David Mayerich
fixed memory leak...
|
238
|
free(temp); //free temporary memory
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
239
|
return true;
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
240
241
|
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
242
243
244
245
|
/// 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...
|
246
|
bool pixel(T * p, unsigned long long n){
|
bfe9c6db
heziqi
added feature_mat...
|
247
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
248
|
unsigned long long N = X() * Y(); //calculate numbers in one band
|
63fc1df8
David Mayerich
added an inverse ...
|
249
|
if ( n >= N){ //make sure the pixel number is right
|
bfe9c6db
heziqi
added feature_mat...
|
250
251
252
253
|
std::cout<<"ERROR: sample or line out of range"<<std::endl;
return false;
}
|
724ec347
David Mayerich
simplified the EN...
|
254
255
256
|
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...
|
257
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
258
|
|
11f177d5
heziqi
Ziqi added functi...
|
259
|
//given a Y ,return a ZX slice
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
260
|
bool read_plane_y(T * p, unsigned long long y){
|
724ec347
David Mayerich
simplified the EN...
|
261
|
return binary<T>::read_plane_2(p, y);
|
11f177d5
heziqi
Ziqi added functi...
|
262
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
263
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
264
265
266
|
/// 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...
|
267
|
/// @param wls is the list of baseline points based on band labels.
|
798cea97
David Mayerich
added progress ba...
|
268
|
bool baseline(std::string outname, std::vector<double> base_pts, unsigned char* mask = NULL, bool PROGRESS = false){
|
1a224b6a
David Mayerich
fixed linux compa...
|
269
|
|
88c3e636
heziqi
Ziqi added big.h
|
270
|
std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file
|
1a224b6a
David Mayerich
fixed linux compa...
|
271
|
|
63fc1df8
David Mayerich
added an inverse ...
|
272
273
274
275
276
277
278
279
|
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 ...
|
280
281
|
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...
|
282
283
284
285
|
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...
|
286
|
|
798cea97
David Mayerich
added progress ba...
|
287
|
pixel(s, n); //retrieve the spectrum s
|
1a224b6a
David Mayerich
fixed linux compa...
|
288
|
base_vals = hsi<T>::interp_spectrum(s, base_pts); //get the values at each baseline point
|
798cea97
David Mayerich
added progress ba...
|
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
|
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...
|
310
|
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 ...
|
311
|
}
|
63fc1df8
David Mayerich
added an inverse ...
|
312
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
313
|
|
63fc1df8
David Mayerich
added an inverse ...
|
314
|
if(PROGRESS) progress = (double)(n+1) / N * 100; //set the current progress
|
6708cc25
heziqi
Ziqi added envi c...
|
315
|
|
63fc1df8
David Mayerich
added an inverse ...
|
316
317
|
target.write((char*)sbc, sizeof(T) * B); //write the corrected data into destination
} //end for each pixel
|
1a224b6a
David Mayerich
fixed linux compa...
|
318
|
|
63fc1df8
David Mayerich
added an inverse ...
|
319
320
|
free(s);
free(sbc);
|
88c3e636
heziqi
Ziqi added big.h
|
321
|
target.close();
|
1a224b6a
David Mayerich
fixed linux compa...
|
322
323
324
|
return true;
|
88c3e636
heziqi
Ziqi added big.h
|
325
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
326
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
327
328
329
330
331
|
/// 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...
|
332
|
bool ratio(std::string outname, double w, unsigned char* mask = NULL, bool PROGRESS = false)
|
88c3e636
heziqi
Ziqi added big.h
|
333
|
{
|
88c3e636
heziqi
Ziqi added big.h
|
334
|
std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file
|
1a224b6a
David Mayerich
fixed linux compa...
|
335
336
|
std::string headername = outname + ".hdr"; //the header file name
|
63fc1df8
David Mayerich
added an inverse ...
|
337
338
339
340
341
|
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...
|
342
|
if(mask != NULL && !mask[n]) //if the normalization band is below threshold
|
63fc1df8
David Mayerich
added an inverse ...
|
343
344
|
memset(s, 0, sizeof(T) * B); //set the output to zero
else{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
345
|
pixel(s, n); //retrieve the spectrum s
|
1a224b6a
David Mayerich
fixed linux compa...
|
346
347
|
nv = hsi<T>::interp_spectrum(s, w); //find the value of the normalization band
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
348
|
for(unsigned long long b = 0; b < B; b++) //for each band in the spectrum
|
63fc1df8
David Mayerich
added an inverse ...
|
349
|
s[b] /= nv; //divide by the normalization value
|
88c3e636
heziqi
Ziqi added big.h
|
350
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
351
|
|
63fc1df8
David Mayerich
added an inverse ...
|
352
|
if(PROGRESS) progress = (double)(n+1) / N * 100; //set the current progress
|
6708cc25
heziqi
Ziqi added envi c...
|
353
|
|
63fc1df8
David Mayerich
added an inverse ...
|
354
355
|
target.write((char*)s, sizeof(T) * B); //write the corrected data into destination
} //end for each pixel
|
1a224b6a
David Mayerich
fixed linux compa...
|
356
|
|
63fc1df8
David Mayerich
added an inverse ...
|
357
|
free(s);
|
88c3e636
heziqi
Ziqi added big.h
|
358
|
target.close();
|
88c3e636
heziqi
Ziqi added big.h
|
359
360
|
return true;
}
|
c30dd3e6
David Mayerich
added vector norm...
|
361
|
|
904614c3
David Mayerich
added normalizati...
|
362
363
364
365
366
367
368
369
370
371
372
373
374
375
|
void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){
std::ofstream target(outfile.c_str(), std::ios::binary); //open the target binary file
file.seekg(0, std::ios::beg); //move the pointer to the current file to the beginning
size_t B = Z(); //number of spectral components
size_t XY = X() * Y(); //calculate the number of pixels
size_t Bb = B * sizeof(T); //number of bytes in a spectrum
T* spec = (T*) malloc(Bb); //allocate space for the spectrum
T len;
for(size_t xy = 0; xy < XY; xy++){ //for each pixel
memset(spec, 0, Bb); //set the spectrum to zero
if(mask == NULL || mask[xy]){ //if the pixel is masked
|
334547d8
sam
changes for mnf: ...
|
376
|
len = 0; //initialize the
|
904614c3
David Mayerich
added normalizati...
|
377
378
379
380
381
382
383
384
385
386
387
|
file.read((char*)spec, Bb); //read a spectrum
for(size_t b = 0; b < B; b++) //for each band
len += spec[b]*spec[b]; //add the square of the spectral band
len = sqrt(len); //calculate the square of the sum of squared components
for(size_t b = 0; b < B; b++) //for each band
spec[b] /= len; //divide by the length
}
else
file.seekg(Bb, std::ios::cur); //otherwise skip a spectrum
target.write((char*)spec, Bb); //output the normalized spectrum
if(PROGRESS) progress = (double)(xy + 1) / (double)XY * 100; //update the progress
|
334547d8
sam
changes for mnf: ...
|
388
|
}
|
904614c3
David Mayerich
added normalizati...
|
389
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
390
|
|
f6169dea
heziqi
Ziqi completed bi...
|
391
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
392
393
394
|
/// 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 ...
|
395
|
bool bil(std::string outname, bool PROGRESS = false)
|
11f177d5
heziqi
Ziqi added functi...
|
396
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
397
|
unsigned long long S = X() * Z() * sizeof(T); //calculate the number of bytes in a ZX slice
|
1a224b6a
David Mayerich
fixed linux compa...
|
398
|
|
11f177d5
heziqi
Ziqi added functi...
|
399
|
std::ofstream target(outname.c_str(), std::ios::binary);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
400
|
//std::string headername = outname + ".hdr";
|
1a224b6a
David Mayerich
fixed linux compa...
|
401
|
|
f6169dea
heziqi
Ziqi completed bi...
|
402
|
T * p; //pointer to the current ZX slice for bip file
|
11f177d5
heziqi
Ziqi added functi...
|
403
|
p = (T*)malloc(S);
|
f6169dea
heziqi
Ziqi completed bi...
|
404
405
|
T * q; //pointer to the current XZ slice for bil file
q = (T*)malloc(S);
|
1a224b6a
David Mayerich
fixed linux compa...
|
406
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
407
|
for ( unsigned long long i = 0; i < Y(); i++)
|
1a224b6a
David Mayerich
fixed linux compa...
|
408
|
{
|
724ec347
David Mayerich
simplified the EN...
|
409
|
read_plane_y(p, i);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
410
|
for ( unsigned long long k = 0; k < Z(); k++)
|
f6169dea
heziqi
Ziqi completed bi...
|
411
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
412
413
|
unsigned long long ks = k * X();
for ( unsigned long long j = 0; j < X(); j++)
|
724ec347
David Mayerich
simplified the EN...
|
414
|
q[ks + j] = p[k + j * Z()];
|
cbce396e
David Mayerich
added support for...
|
415
|
|
63fc1df8
David Mayerich
added an inverse ...
|
416
|
if(PROGRESS) progress = (double)(i * Z() + k+1) / (Y() * Z()) * 100;
|
f6169dea
heziqi
Ziqi completed bi...
|
417
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
418
|
target.write(reinterpret_cast<const char*>(q), S); //write a band data into target file
|
11f177d5
heziqi
Ziqi added functi...
|
419
|
}
|
f6169dea
heziqi
Ziqi completed bi...
|
420
|
|
11f177d5
heziqi
Ziqi added functi...
|
421
|
free(p);
|
f6169dea
heziqi
Ziqi completed bi...
|
422
|
free(q);
|
11f177d5
heziqi
Ziqi added functi...
|
423
424
425
|
target.close();
return true;
}
|
88c3e636
heziqi
Ziqi added big.h
|
426
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
427
428
429
430
431
432
433
434
|
/// 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...
|
435
436
|
bool baseline_band(double lb, double rb, T* lp, T* rp, double wavelength, T* result){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
437
|
unsigned long long XY = X() * Y();
|
20c212c0
heziqi
Ziqi added functi...
|
438
439
440
441
|
band(result, wavelength); //get band
//perform the baseline correction
double r = (double) (wavelength - lb) / (double) (rb - lb);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
442
|
for(unsigned long long i=0; i < XY; i++){
|
20c212c0
heziqi
Ziqi added functi...
|
443
444
445
446
|
result[i] =(T) (result[i] - (rp[i] - lp[i]) * r - lp[i] );
}
return true;
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
447
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
448
449
450
451
452
453
|
/// 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...
|
454
|
bool height(double lb, double rb, double bandwavelength, T* result){
|
20c212c0
heziqi
Ziqi added functi...
|
455
456
457
|
T* lp;
T* rp;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
458
459
|
unsigned long long XY = X() * Y();
unsigned long long S = XY * sizeof(T);
|
70407ea9
heziqi
Ziqi completed he...
|
460
461
|
lp = (T*) malloc(S); //memory allocation
rp = (T*) malloc(S);
|
20c212c0
heziqi
Ziqi added functi...
|
462
|
|
70407ea9
heziqi
Ziqi completed he...
|
463
|
band(lp, lb);
|
1a224b6a
David Mayerich
fixed linux compa...
|
464
|
band(rp, rb);
|
20c212c0
heziqi
Ziqi added functi...
|
465
466
467
|
baseline_band(lb, rb, lp, rp, bandwavelength, result);
|
517876d6
heziqi
metrics finished ...
|
468
469
|
free(lp);
free(rp);
|
20c212c0
heziqi
Ziqi added functi...
|
470
471
472
|
return true;
}
|
c359422d
heziqi
Ziqi added functi...
|
473
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
474
475
476
477
478
479
480
|
/// 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...
|
481
|
bool area(double lb, double rb, double lab, double rab, T* result){
|
20c212c0
heziqi
Ziqi added functi...
|
482
483
484
485
486
|
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...
|
487
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
488
489
|
unsigned long long XY = X() * Y();
unsigned long long S = XY * sizeof(T);
|
20c212c0
heziqi
Ziqi added functi...
|
490
491
492
493
494
|
lp = (T*) malloc(S); //memory allocation
rp = (T*) malloc(S);
cur = (T*) malloc(S);
cur2 = (T*) malloc(S);
|
20c212c0
heziqi
Ziqi added functi...
|
495
496
497
498
|
memset(result, (char)0, S);
//find the wavelenght position in the whole band
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
499
500
501
|
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...
|
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
|
//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...
|
517
|
while (lab >= w[ai]){
|
20c212c0
heziqi
Ziqi added functi...
|
518
519
|
ai++;
}
|
70407ea9
heziqi
Ziqi completed he...
|
520
|
while (rab <= w[bi]){
|
20c212c0
heziqi
Ziqi added functi...
|
521
522
523
524
525
526
527
|
bi--;
}
band(lp, lb);
band(rp, rb);
//calculate the beginning and the ending part
|
70407ea9
heziqi
Ziqi completed he...
|
528
529
|
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...
|
530
531
|
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...
|
532
|
}
|
70407ea9
heziqi
Ziqi completed he...
|
533
534
|
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...
|
535
|
for(unsigned long long j = 0; j < XY; j++){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
536
|
result[j] += (T)((w[ai] - lab) * ((double)cur[j] + (double)cur2[j]) / 2.0);
|
20c212c0
heziqi
Ziqi added functi...
|
537
|
}
|
517876d6
heziqi
metrics finished ...
|
538
|
|
20c212c0
heziqi
Ziqi added functi...
|
539
540
|
//calculate the area
ai++;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
541
|
for(unsigned long long i = ai; i <= bi ;i++)
|
20c212c0
heziqi
Ziqi added functi...
|
542
543
|
{
baseline_band(lb, rb, lp, rp, w[ai], cur2);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
544
|
for(unsigned long long j = 0; j < XY; j++)
|
20c212c0
heziqi
Ziqi added functi...
|
545
|
{
|
1a224b6a
David Mayerich
fixed linux compa...
|
546
|
result[j] += (T)((w[ai] - w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
|
20c212c0
heziqi
Ziqi added functi...
|
547
548
549
|
}
std::swap(cur,cur2); //swap the band pointers
}
|
70407ea9
heziqi
Ziqi completed he...
|
550
|
|
517876d6
heziqi
metrics finished ...
|
551
552
553
554
|
free(lp);
free(rp);
free(cur);
free(cur2);
|
20c212c0
heziqi
Ziqi added functi...
|
555
556
557
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
558
559
560
561
562
563
564
565
566
|
/// 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...
|
567
|
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...
|
568
|
|
1a224b6a
David Mayerich
fixed linux compa...
|
569
570
571
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
70407ea9
heziqi
Ziqi completed he...
|
572
573
574
575
|
//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...
|
576
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
70407ea9
heziqi
Ziqi completed he...
|
577
578
579
580
581
|
if(p1[i] == 0 && p2[i] ==0)
result[i] = 1;
else
result[i] = p1[i] / p2[i];
}
|
517876d6
heziqi
metrics finished ...
|
582
583
584
|
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
585
|
return true;
|
70407ea9
heziqi
Ziqi completed he...
|
586
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
587
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
588
589
590
591
592
593
594
595
596
|
/// 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...
|
597
598
|
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...
|
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);
height(lb2, rb2, pos, 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
618
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
619
620
621
622
623
624
625
626
627
|
/// 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...
|
628
|
/// @param rab2 is the label value for the right bound (end of the integration) of the second peak (denominator)
|
a23c4132
David Mayerich
Doxygen comments ...
|
629
|
/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
|
ba51ae6a
David Mayerich
fixed metric calc...
|
630
631
|
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...
|
632
633
634
635
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
70407ea9
heziqi
Ziqi completed he...
|
636
637
638
639
|
//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...
|
640
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
70407ea9
heziqi
Ziqi completed he...
|
641
642
643
644
645
|
if(p1[i] == 0 && p2[i] ==0)
result[i] = 1;
else
result[i] = p1[i] / p2[i];
}
|
517876d6
heziqi
metrics finished ...
|
646
647
648
|
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
649
650
|
return true;
}
|
517876d6
heziqi
metrics finished ...
|
651
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
652
653
654
655
656
657
658
|
/// 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 ...
|
659
660
661
662
663
664
|
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...
|
665
666
|
unsigned long long XY = X() * Y();
unsigned long long S = XY * sizeof(T);
|
517876d6
heziqi
metrics finished ...
|
667
668
669
670
671
672
673
674
675
|
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...
|
676
677
678
|
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 ...
|
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
|
//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...
|
705
706
|
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 ...
|
707
708
709
|
}
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...
|
710
|
for(unsigned long long j = 0; j < XY; j++){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
711
|
result[j] += (T)((w[ai] - lab) * (w[ai] + lab) * ((double)cur[j] + (double)cur2[j]) / 4.0);
|
517876d6
heziqi
metrics finished ...
|
712
713
714
715
|
}
//calculate f(x) times x
ai++;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
716
|
for(unsigned long long i = ai; i <= bi ;i++)
|
517876d6
heziqi
metrics finished ...
|
717
718
|
{
baseline_band(lb, rb, lp, rp, w[ai], cur2);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
719
|
for(unsigned long long j = 0; j < XY; j++)
|
517876d6
heziqi
metrics finished ...
|
720
|
{
|
1a224b6a
David Mayerich
fixed linux compa...
|
721
|
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 ...
|
722
723
724
725
726
727
728
729
730
731
732
|
}
std::swap(cur,cur2); //swap the band pointers
}
free(lp);
free(rp);
free(cur);
free(cur2);
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
733
734
735
736
737
738
739
|
/// 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...
|
740
|
bool centroid(T* result, double lb, double rb, double lab, double rab, unsigned char* mask = NULL){
|
1a224b6a
David Mayerich
fixed linux compa...
|
741
742
743
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
517876d6
heziqi
metrics finished ...
|
744
745
746
747
|
//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...
|
748
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
ba51ae6a
David Mayerich
fixed metric calc...
|
749
|
if(mask == NULL || mask[i])
|
517876d6
heziqi
metrics finished ...
|
750
751
752
753
754
|
result[i] = p1[i] / p2[i];
}
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
755
|
return true;
|
517876d6
heziqi
metrics finished ...
|
756
757
|
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
758
759
760
761
762
763
764
|
/// Create a mask based on a given band and threshold value.
/// All pixels in the
/// specified band greater than the threshold are true and all pixels less than the threshold are false.
/// @param mask_band is the band used to specify the mask
/// @param threshold is the threshold used to determine if the mask value is true or false
/// @param p is a pointer to a pre-allocated array at least X * Y in size
|
e37ce7bf
David Mayerich
added the ability...
|
765
|
bool build_mask(unsigned char* out_mask, double mask_band, double lower, double upper, unsigned char* mask = NULL, bool PROGRESS = false){
|
4a6f666c
heziqi
Added mask method
|
766
|
|
724ec347
David Mayerich
simplified the EN...
|
767
|
T* temp = (T*)malloc(X() * Y() * sizeof(T)); //allocate memory for the certain band
|
798cea97
David Mayerich
added progress ba...
|
768
|
band(temp, mask_band, PROGRESS);
|
4a6f666c
heziqi
Added mask method
|
769
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
770
|
for (unsigned long long i = 0; i < X() * Y();i++) {
|
e37ce7bf
David Mayerich
added the ability...
|
771
772
773
774
775
776
777
|
if(mask == NULL || mask[i] != 0){
if(temp[i] > lower && temp[i] < upper){
out_mask[i] = 255;
}
else
out_mask[i] = 0;
}
|
4a6f666c
heziqi
Added mask method
|
778
779
|
}
|
517876d6
heziqi
metrics finished ...
|
780
|
free(temp);
|
4a6f666c
heziqi
Added mask method
|
781
782
783
784
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
785
786
787
788
|
/// 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...
|
789
|
bool apply_mask(std::string outfile, unsigned char* p, bool PROGRESS = false){
|
740f8cd2
heziqi
added apply_mask
|
790
791
792
|
std::ofstream target(outfile.c_str(), std::ios::binary);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
793
794
|
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
|
795
|
|
1a55a328
David Mayerich
Added comments fo...
|
796
|
T * temp = (T*)malloc(L); //allocate space for that page
|
740f8cd2
heziqi
added apply_mask
|
797
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
798
|
for (unsigned long long i = 0; i < Y(); i++) //for each page (Y in BIP)
|
740f8cd2
heziqi
added apply_mask
|
799
|
{
|
724ec347
David Mayerich
simplified the EN...
|
800
|
read_plane_y(temp, i); //load that page (it's pointed to by temp)
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
801
|
for ( unsigned long long j = 0; j < X(); j++) //for each X value
|
740f8cd2
heziqi
added apply_mask
|
802
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
803
|
for (unsigned long long k = 0; k < Z(); k++) //for each B value (band)
|
740f8cd2
heziqi
added apply_mask
|
804
|
{
|
724ec347
David Mayerich
simplified the EN...
|
805
806
|
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...
|
807
|
else //otherwise just continue
|
740f8cd2
heziqi
added apply_mask
|
808
809
810
|
continue;
}
}
|
1a55a328
David Mayerich
Added comments fo...
|
811
|
target.write(reinterpret_cast<const char*>(temp), L); //write the edited band data into target file
|
6a46c8ff
David Mayerich
fixed bugs in app...
|
812
|
if(PROGRESS) progress = (double)(i+1) / (double)Y() * 100;
|
740f8cd2
heziqi
added apply_mask
|
813
|
}
|
1a55a328
David Mayerich
Added comments fo...
|
814
815
816
|
target.close(); //close the target file
free(temp); //free allocated memory
return true; //return true
|
740f8cd2
heziqi
added apply_mask
|
817
818
|
}
|
ba51ae6a
David Mayerich
fixed metric calc...
|
819
820
821
822
823
|
/// Copies all spectra corresponding to nonzero values of a mask into a pre-allocated matrix of size (B x P)
/// where P is the number of masked pixels and B is the number of bands. The allocated memory can be accessed
/// using the following indexing: i = p*B + b
/// @param matrix is the destination for the pixel data
/// @param mask is the mask
|
814eb271
David Mayerich
fixed problems wi...
|
824
|
bool sift(T* matrix, unsigned char* mask = NULL, bool PROGRESS = false){
|
ba51ae6a
David Mayerich
fixed metric calc...
|
825
826
827
828
829
|
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...
|
830
|
|
ba51ae6a
David Mayerich
fixed metric calc...
|
831
832
|
size_t p = 0; //create counter variables
for(size_t xy = 0; xy < XY; xy++){ //for each pixel
|
d4aa4cf7
David Mayerich
fixed the in-memo...
|
833
|
if(mask == NULL || mask[xy]){ //if the current pixel is masked
|
ba51ae6a
David Mayerich
fixed metric calc...
|
834
835
836
837
838
839
840
841
|
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...
|
842
|
file.seekg(Bbytes, std::ios::cur); //otherwise skip this band
|
814eb271
David Mayerich
fixed problems wi...
|
843
|
if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
|
ba51ae6a
David Mayerich
fixed metric calc...
|
844
845
846
|
}
return true;
}
|
e843658b
Brad Deutsch
Previous push did...
|
847
|
|
2fcab4f0
David Mayerich
added unsifting f...
|
848
|
/// Saves to disk only those spectra corresponding to mask values != 0
|
c1078e19
David Mayerich
HSIproc bug fixes...
|
849
|
bool sift(std::string outfile, unsigned char* mask, bool PROGRESS = false){
|
8550e94f
David Mayerich
added sifting and...
|
850
851
|
//reset the file pointer to the beginning of the file
|
1a224b6a
David Mayerich
fixed linux compa...
|
852
|
file.seekg(0, std::ios::beg);
|
8550e94f
David Mayerich
added sifting and...
|
853
854
|
// open an output stream
|
e843658b
Brad Deutsch
Previous push did...
|
855
856
|
std::ofstream target(outfile.c_str(), std::ios::binary);
|
8550e94f
David Mayerich
added sifting and...
|
857
|
//allocate space for a single spectrum
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
858
|
unsigned long long B = Z();
|
8550e94f
David Mayerich
added sifting and...
|
859
|
T* spectrum = (T*) malloc(B * sizeof(T));
|
e843658b
Brad Deutsch
Previous push did...
|
860
|
|
8550e94f
David Mayerich
added sifting and...
|
861
|
//calculate the number of pixels in a band
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
862
|
unsigned long long XY = X() * Y();
|
e843658b
Brad Deutsch
Previous push did...
|
863
|
|
8550e94f
David Mayerich
added sifting and...
|
864
|
//for each pixel
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
865
866
|
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...
|
867
868
869
870
871
|
//if the current pixel isn't masked
if( mask[x] == 0){
//increment the number of skipped pixels
skip++;
|
e843658b
Brad Deutsch
Previous push did...
|
872
|
}
|
8550e94f
David Mayerich
added sifting and...
|
873
874
875
876
877
878
879
880
881
882
883
884
885
886
|
//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...
|
887
|
}
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
888
|
if(PROGRESS) progress = (double) (x+1) / XY * 100;
|
8550e94f
David Mayerich
added sifting and...
|
889
|
|
e843658b
Brad Deutsch
Previous push did...
|
890
|
}
|
8550e94f
David Mayerich
added sifting and...
|
891
892
|
//close the output file
|
e843658b
Brad Deutsch
Previous push did...
|
893
|
target.close();
|
8550e94f
David Mayerich
added sifting and...
|
894
|
free(spectrum);
|
cbce396e
David Mayerich
added support for...
|
895
|
|
cbce396e
David Mayerich
added support for...
|
896
|
return true;
|
8550e94f
David Mayerich
added sifting and...
|
897
898
|
}
|
a2bf1d08
David Mayerich
general bug fixes...
|
899
|
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...
|
900
901
902
903
904
|
// 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...
|
905
|
file.seekg(0, std::ios::beg);
|
8550e94f
David Mayerich
added sifting and...
|
906
907
|
//allocate space for a single spectrum
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
908
|
unsigned long long B = Z();
|
8550e94f
David Mayerich
added sifting and...
|
909
910
911
912
913
914
915
|
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...
|
916
|
unsigned long long XY = samples * lines;
|
8550e94f
David Mayerich
added sifting and...
|
917
918
|
//for each pixel
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
919
920
|
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...
|
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
|
//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...
|
938
|
if(PROGRESS) progress = (double)(x + 1) / XY * 100;
|
3150c537
David Mayerich
added further pro...
|
939
|
|
8550e94f
David Mayerich
added sifting and...
|
940
941
942
943
944
945
|
}
//close the output file
target.close();
free(spectrum);
|
a2bf1d08
David Mayerich
general bug fixes...
|
946
|
//progress = 100;
|
3150c537
David Mayerich
added further pro...
|
947
|
|
e843658b
Brad Deutsch
Previous push did...
|
948
|
return true;
|
8550e94f
David Mayerich
added sifting and...
|
949
950
|
|
e843658b
Brad Deutsch
Previous push did...
|
951
952
|
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
953
954
955
956
|
/// 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 ...
|
957
|
bool avg_band(double* p, unsigned char* mask = NULL, bool PROGRESS = false){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
958
959
960
961
962
963
964
965
966
967
968
969
970
971
|
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...
|
972
973
|
}
}
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
974
|
if(PROGRESS) progress = (double)(i+1) / XY * 100; //increment the progress
|
a43c4fe1
heziqi
Added crop in env...
|
975
|
}
|
2aa68315
David Mayerich
major bug fixes f...
|
976
|
|
a43c4fe1
heziqi
Added crop in env...
|
977
|
free(temp);
|
a43c4fe1
heziqi
Added crop in env...
|
978
979
|
return true;
}
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
980
981
|
#ifdef CUDA_FOUND
/// Calculate the covariance matrix for masked pixels using cuBLAS
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
982
|
/// Note that cuBLAS only supports integer-sized arrays, so there may be issues with large spectra
|
63fc1df8
David Mayerich
added an inverse ...
|
983
|
bool co_matrix_cublas(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
984
|
|
1a224b6a
David Mayerich
fixed linux compa...
|
985
|
cudaError_t cudaStat;
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
986
987
|
cublasStatus_t stat;
cublasHandle_t handle;
|
1a224b6a
David Mayerich
fixed linux compa...
|
988
|
|
63fc1df8
David Mayerich
added an inverse ...
|
989
|
progress = 0; //initialize the progress to zero (0)
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
990
991
992
993
994
995
996
997
998
999
1000
|
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...
|
1001
|
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...
|
1002
|
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
1003
1004
1005
1006
1007
1008
1009
1010
1011
|
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...
|
1012
1013
|
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...
|
1014
1015
1016
|
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...
|
1017
|
}
|
63fc1df8
David Mayerich
added an inverse ...
|
1018
|
if(PROGRESS) progress = (double)(xy+1) / XY * 100; //record the current progress
|
5a5359da
David Mayerich
fixed - allows th...
|
1019
|
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
1020
1021
|
}
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1022
|
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...
|
1023
1024
1025
1026
1027
|
cudaFree(A_dev); //clean up allocated device memory
cudaFree(s_dev);
cudaFree(avg_dev);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1028
1029
|
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...
|
1030
1031
1032
1033
|
co[B * i + j] = co[B * j + i];
}
}
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
1034
1035
1036
|
return true;
}
#endif
|
1a224b6a
David Mayerich
fixed linux compa...
|
1037
|
|
2aa68315
David Mayerich
major bug fixes f...
|
1038
|
/// Calculate the covariance matrix for all masked pixels in the image with 64-bit floating point precision.
|
a23c4132
David Mayerich
Doxygen comments ...
|
1039
1040
1041
1042
|
/// @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 ...
|
1043
|
bool co_matrix(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
|
c6f526c2
David Mayerich
added cuBLAS supp...
|
1044
1045
|
#ifdef CUDA_FOUND
|
b37665a9
David Mayerich
cuBLAS won't run ...
|
1046
1047
1048
1049
1050
1051
|
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...
|
1052
|
#endif
|
63fc1df8
David Mayerich
added an inverse ...
|
1053
|
progress = 0;
|
a43c4fe1
heziqi
Added crop in env...
|
1054
|
//memory allocation
|
2aa68315
David Mayerich
major bug fixes f...
|
1055
1056
|
unsigned long long XY = X() * Y();
unsigned long long B = Z();
|
a43c4fe1
heziqi
Added crop in env...
|
1057
|
T* temp = (T*)malloc(sizeof(T) * B);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1058
1059
|
unsigned long long count = nnz(mask); //count the number of masked pixels
|
63fc1df8
David Mayerich
added an inverse ...
|
1060
|
|
2aa68315
David Mayerich
major bug fixes f...
|
1061
1062
1063
1064
|
//initialize covariance matrix
memset(co, 0, B * B * sizeof(double));
//calculate covariance matrix
|
2aa68315
David Mayerich
major bug fixes f...
|
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
|
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...
|
1078
1079
|
}
}
|
63fc1df8
David Mayerich
added an inverse ...
|
1080
|
if(PROGRESS) progress = (double)(xy+1) / XY * 100;
|
a43c4fe1
heziqi
Added crop in env...
|
1081
|
}
|
2aa68315
David Mayerich
major bug fixes f...
|
1082
1083
1084
1085
|
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...
|
1086
1087
1088
|
}
}
|
a43c4fe1
heziqi
Added crop in env...
|
1089
|
free(temp);
|
2aa68315
David Mayerich
major bug fixes f...
|
1090
|
free(temp_precise);
|
a43c4fe1
heziqi
Added crop in env...
|
1091
1092
1093
|
return true;
}
|
334547d8
sam
changes for mnf: ...
|
1094
|
|
58e55e6b
David Mayerich
fixed colormap.h ...
|
1095
|
#ifdef CUDA_FOUND
|
334547d8
sam
changes for mnf: ...
|
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
|
/// Calculate the covariance matrix of Noise for masked pixels using cuBLAS
/// Note that cuBLAS only supports integer-sized arrays, so there may be issues with large spectra
bool coNoise_matrix_cublas(double* coN, double* avg, unsigned char *mask, bool PROGRESS = false){
cudaError_t cudaStat;
cublasStatus_t stat;
cublasHandle_t handle;
progress = 0; //initialize the progress to zero (0)
unsigned long long XY = X() * Y(); //calculate the number of elements in a band image
unsigned long long B = Z(); //calculate the number of spectral elements
double* s = (double*)malloc(sizeof(double) * B); //allocate space for the spectrum that will be pulled from the file
double* s_dev; //declare a device pointer that will store the spectrum on the GPU
double* s2_dev; // device pointer on the GPU
cudaStat = cudaMalloc(&s2_dev, B * sizeof(double)); // allocate space on the CUDA device
cudaStat = cudaMemset(s2_dev, 0, B * sizeof(double)); // initialize s2_dev to zero (0)
double* A_dev; //declare a device pointer that will store the covariance matrix on the GPU
double* avg_dev; //declare a device pointer that will store the average spectrum
cudaStat = cudaMalloc(&s_dev, B * sizeof(double)); //allocate space on the CUDA device for the spectrum
cudaStat = cudaMalloc(&A_dev, B * B * sizeof(double)); //allocate space on the CUDA device for the covariance matrix
cudaStat = cudaMemset(A_dev, 0, B * B * sizeof(double)); //initialize the covariance matrix to zero (0)
cudaStat = cudaMalloc(&avg_dev, B * sizeof(double)); //allocate space on the CUDA device for the average spectrum
stat = cublasSetVector((int)B, sizeof(double), avg, 1, avg_dev, 1); //copy the average spectrum to the CUDA device
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
if (mask == NULL || mask[xy] != 0){
pixeld(s, xy); //retreive the spectrum at the current xy pixel location
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
cudaMemcpy(s2_dev, s_dev + 1 , (B-1) * sizeof(double), cudaMemcpyDeviceToDevice); //copy B-1 elements from shifted source data (s_dev) to device pointer (s2_dev )
stat = cublasDaxpy(handle, (int)B, &axpy_alpha, s2_dev, 1, s_dev, 1); //Minimum/Maximum Autocorrelation Factors (MAF) method : subtranct each pixel from adjacent pixel (z direction is choosed to do so , which is almost the same as x or y direction or even average of them )
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)
}
if(PROGRESS) progress = (double)(xy+1) / XY * 100; //record the current progress
}
cublasGetMatrix((int)B, (int)B, sizeof(double), A_dev, (int)B, coN, (int)B); //copy the result from the GPU to the CPU
cudaFree(A_dev); //clean up allocated device memory
cudaFree(s_dev);
cudaFree(s2_dev);
cudaFree(avg_dev);
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++){
coN[B * i + j] = coN[B * j + i];
}
}
return true;
}
#endif
/// Calculate the covariance of noise matrix for all masked pixels in the image with 64-bit floating point precision.
/// @param coN is a pointer to pre-allocated memory of size [B * B] that stores the resulting covariance matrix
/// @param avg is a pointer to memory of size B that stores the average spectrum
/// @param mask is a pointer to memory of size [X * Y] that stores the mask value at each pixel location
bool coNoise_matrix(double* coN, double* avg, unsigned char *mask, bool PROGRESS = false){
#ifdef CUDA_FOUND
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 coNoise_matrix_cublas(coN, avg, mask, PROGRESS); //use cuBLAS to calculate the covariance matrix
#endif
progress = 0;
//memory allocation
unsigned long long XY = X() * Y();
unsigned long long B = Z();
T* temp = (T*)malloc(sizeof(T) * B);
unsigned long long count = nnz(mask); //count the number of masked pixels
//initialize covariance matrix of noise
memset(coN, 0, B * B * sizeof(double));
//calculate covariance matrix
double* coN_half = (double*) malloc(B * B * sizeof(double)); //allocate space for a higher-precision intermediate matrix
double* temp_precise = (double*) malloc(B * sizeof(double));
memset(coN_half, 0, B * B * sizeof(double)); //initialize the high-precision matrix with zeros
unsigned long long idx; //stores i*B to speed indexing
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];
for(unsigned long long b2 = 0; b2 < B-1; b2++) //Minimum/Maximum Autocorrelation Factors (MAF) method : subtranct each pixel from adjacent pixel (z direction is choosed to do so , which is almost the same as x or y direction or even average of them )
temp_precise[b2] -= temp_precise[b2+1];
idx = 0;
for (unsigned long long b0 = 0; b0 < B; b0++){ //for each band
for (unsigned long long b1 = b0; b1 < B; b1++)
coN_half[idx++] += temp_precise[b0] * temp_precise[b1];
}
}
if(PROGRESS) progress = (double)(xy+1) / XY * 100;
}
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++){
coN[j * B + i] = coN[i * B + j] = coN_half[idx++] / (double) count;
}
}
free(temp);
free(temp_precise);
return true;
}
#ifdef CUDA_FOUND
/// Project the spectra onto a set of basis functions
/// @param outfile is the name of the new binary output file that will be created
/// @param center is a spectrum about which the data set will be rotated (ex. when performing mean centering)
/// @param basis a set of basis vectors that the data set will be projected onto (after centering)
/// @param M is the number of basis vectors
/// @param mask is a character mask used to limit processing to valid pixels
bool project_cublas(std::string outfile, double* center, double* basis, unsigned long long M, unsigned char* mask = NULL, bool PROGRESS = false){
cudaError_t cudaStat;
cublasStatus_t stat;
cublasHandle_t handle;
std::ofstream target(outfile.c_str(), std::ios::binary); //open the target binary file
|
334547d8
sam
changes for mnf: ...
|
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
|
progress = 0; //initialize the progress to zero (0)
unsigned long long XY = X() * Y(); //calculate the number of elements in a band image
unsigned long long B = Z(); //calculate the number of spectral elements
double* s = (double*)malloc(sizeof(double) * B); //allocate space for the spectrum that will be pulled from the file
double* s_dev; //declare a device pointer that will store the spectrum on the GPU
cudaStat = cudaMalloc(&s_dev, B * sizeof(double)); //allocate space on the CUDA device for the spectrum
double* basis_dev; // device pointer on the GPU
cudaStat = cudaMalloc(&basis_dev, M * B * sizeof(double)); // allocate space on the CUDA device
cudaStat = cudaMemset(basis_dev, 0, M * B * sizeof(double)); // initialize basis_dev to zero (0)
/// transposing basis matrix (because cuBLAS is column-major)
double *basis_Transposed = (double*)malloc(M * B * sizeof(double));
memset(basis_Transposed, 0, M * B * sizeof(double));
for (int i = 0; i<M; i++)
for (int j = 0; j<B; j++)
basis_Transposed[i+j*M] = basis[i*B+j];
stat = cublasSetMatrix((int)M, (int)B, sizeof(double),basis_Transposed, (int)M, basis_dev, (int)M); //copy the basis_Transposed matrix to the CUDA device (both matrices are stored in column-major format)
double* center_dev; //declare a device pointer that will store the center (average)
cudaStat = cudaMalloc(¢er_dev, B * sizeof(double)); //allocate space on the CUDA device for the center (average)
stat = cublasSetVector((int)B, sizeof(double), center, 1, center_dev, 1); //copy the center vector (average) to the CUDA device (from host to device)
double* A = (double*)malloc(sizeof(double) * M); //allocate space for the projected pixel on the host
double* A_dev; //declare a device pointer that will store the projected pixel on the GPU
cudaStat = cudaMalloc(&A_dev,M * sizeof(double)); //allocate space on the CUDA device for the projected pixel
cudaStat = cudaMemset(A_dev, 0,M * sizeof(double)); //initialize the projected pixel to zero (0)
double axpy_alpha = -1; //multiplication factor for the center (in order to perform a subtraction)
double axpy_alpha2 = 1; //multiplication factor for the matrix-vector multiplication
double axpy_beta = 0; //multiplication factor for the matrix-vector multiplication (there is no second scalor)
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;
}
|
08d1c3b9
David Mayerich
fixed casting war...
|
1286
1287
1288
|
T* temp = (T*)malloc(sizeof(T) * M); //allocate space for the projected pixel to be written on the disc
size_t i;
for (unsigned long long xy = 0; xy < XY; xy++){ //for each pixel
|
334547d8
sam
changes for mnf: ...
|
1289
|
if (mask == NULL || mask[xy] != 0){
|
08d1c3b9
David Mayerich
fixed casting war...
|
1290
|
pixeld(s, xy); //retreive the spectrum at the current xy pixel location
|
334547d8
sam
changes for mnf: ...
|
1291
1292
1293
1294
|
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, center_dev, 1, s_dev, 1); //subtract the center (average)
stat = cublasDgemv(handle,CUBLAS_OP_N,(int)M,(int)B,&axpy_alpha2,basis_dev,(int)M,s_dev,1,&axpy_beta,A_dev,1); //performs the matrix-vector multiplication
|
3a880d1c
David Mayerich
fixed PCA and BIP...
|
1295
1296
1297
|
stat = cublasGetVector((int)M, sizeof(double), A_dev, 1, A, 1); //copy the projected pixel to the host (from GPU to CPU)
//stat = cublasGetVector((int)B, sizeof(double), s_dev, 1, A, 1); //copy the projected pixel to the host (from GPU to CPU)
|
334547d8
sam
changes for mnf: ...
|
1298
|
|
08d1c3b9
David Mayerich
fixed casting war...
|
1299
1300
|
//std::copy<double*, T*>(A, A + M, temp);
for(i = 0; i < M; i++) temp[i] = (T)A[i]; //casting projected pixel from double to whatever T is
|
334547d8
sam
changes for mnf: ...
|
1301
|
}
|
3a880d1c
David Mayerich
fixed PCA and BIP...
|
1302
1303
|
else
memset(temp, 0, sizeof(T) * M);
|
334547d8
sam
changes for mnf: ...
|
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
|
target.write(reinterpret_cast<const char*>(temp), sizeof(T) * M); //write the projected vector
if(PROGRESS) progress = (double)(xy+1) / XY * 100; //record the current progress
}
//clean up allocated device memory
cudaFree(A_dev);
cudaFree(s_dev);
cudaFree(basis_dev);
cudaFree(center_dev);
free(A);
free(s);
free(temp);
target.close(); //close the output file
|
08d1c3b9
David Mayerich
fixed casting war...
|
1319
|
|
334547d8
sam
changes for mnf: ...
|
1320
1321
1322
1323
|
return true;
}
#endif
|
812a2741
David Mayerich
added documentati...
|
1324
1325
1326
1327
1328
1329
|
/// Project the spectra onto a set of basis functions
/// @param outfile is the name of the new binary output file that will be created
/// @param center is a spectrum about which the data set will be rotated (ex. when performing mean centering)
/// @param basis a set of basis vectors that the data set will be projected onto (after centering)
/// @param M is the number of basis vectors
/// @param mask is a character mask used to limit processing to valid pixels
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1330
|
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 ...
|
1331
|
|
334547d8
sam
changes for mnf: ...
|
1332
1333
1334
1335
1336
1337
1338
1339
|
#ifdef CUDA_FOUND
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 project_cublas(outfile,center,basis,M,mask,PROGRESS); //use cuBLAS to calculate the covariance matrix
#endif
|
63fc1df8
David Mayerich
added an inverse ...
|
1340
|
std::ofstream target(outfile.c_str(), std::ios::binary); //open the target binary file
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1341
|
//std::string headername = outfile + ".hdr"; //the header file name
|
63fc1df8
David Mayerich
added an inverse ...
|
1342
1343
1344
1345
1346
1347
1348
1349
1350
|
//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...
|
1351
1352
1353
1354
1355
1356
1357
1358
|
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 ...
|
1359
1360
|
}
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
1361
|
|
63fc1df8
David Mayerich
added an inverse ...
|
1362
1363
1364
1365
1366
1367
1368
|
target.write(reinterpret_cast<const char*>(rs), sizeof(T) * M); //write the projected vector
if(PROGRESS) progress = (double)(xy+1) / XY * 100;
}
free(s); //free temporary storage arrays
free(rs);
target.close(); //close the output file
|
08d1c3b9
David Mayerich
fixed casting war...
|
1369
|
|
63fc1df8
David Mayerich
added an inverse ...
|
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
|
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...
|
1391
|
for(unsigned long long b = 0; b < B; b++){ //for each component of the basis vector
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1392
|
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 ...
|
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
|
}
}
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
|
1407
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
1408
1409
1410
1411
1412
1413
1414
|
/// 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...
|
1415
1416
1417
1418
1419
|
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 ...
|
1420
1421
|
unsigned long long b1,
bool PROGRESS = false){
|
a43c4fe1
heziqi
Added crop in env...
|
1422
|
|
3033b886
David Mayerich
implemented spect...
|
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
|
//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...
|
1434
|
T* temp = (T*)malloc(L);
|
3033b886
David Mayerich
implemented spect...
|
1435
1436
|
//open an output file for binary writing
|
a43c4fe1
heziqi
Added crop in env...
|
1437
|
std::ofstream out(outfile.c_str(), std::ios::binary);
|
3033b886
David Mayerich
implemented spect...
|
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
|
//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...
|
1453
|
{
|
3033b886
David Mayerich
implemented spect...
|
1454
|
for (unsigned x = 0; x < samples; x++)
|
a43c4fe1
heziqi
Added crop in env...
|
1455
|
{
|
3033b886
David Mayerich
implemented spect...
|
1456
1457
1458
|
//read the cropped spectral region
file.read( (char*) temp, L );
//pixel(temp, sp + x + y * X());
|
1a224b6a
David Mayerich
fixed linux compa...
|
1459
|
out.write(reinterpret_cast<const char*>(temp), L); //write slice data into target file
|
3150c537
David Mayerich
added further pro...
|
1460
|
|
3033b886
David Mayerich
implemented spect...
|
1461
1462
|
file.seekg(jump_sample, std::ios::cur);
|
63fc1df8
David Mayerich
added an inverse ...
|
1463
|
if(PROGRESS) progress = (double)((y+1) * samples + x + 1) / (lines * samples) * 100;
|
a43c4fe1
heziqi
Added crop in env...
|
1464
|
}
|
3033b886
David Mayerich
implemented spect...
|
1465
1466
|
file.seekg(jump_line, std::ios::cur);
|
a43c4fe1
heziqi
Added crop in env...
|
1467
1468
|
}
free(temp);
|
3150c537
David Mayerich
added further pro...
|
1469
|
|
a43c4fe1
heziqi
Added crop in env...
|
1470
1471
1472
|
return true;
}
|
dc8eb8aa
David Mayerich
added band trimmi...
|
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
|
/// 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...
|
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
|
/// 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...
|
1540
|
if(PROGRESS) progress = (double)( (y+1) * S[0] + 1) / (double) (S[0] * S[1]) * 100;
|
2ce6954b
David Mayerich
added the ability...
|
1541
1542
1543
|
}
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1544
1545
1546
1547
1548
1549
1550
1551
|
/// 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...
|
1552
|
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...
|
1553
|
std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing
|
1a224b6a
David Mayerich
fixed linux compa...
|
1554
|
|
9b2a5d71
David Mayerich
implemented finit...
|
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
|
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...
|
1570
1571
1572
1573
1574
|
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...
|
1575
1576
1577
1578
1579
1580
1581
|
}
}
out.write((char*)outspec, Nb); //output the band
if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
}
}
|
9d77bbd9
David Mayerich
updates for HSIvi...
|
1582
|
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...
|
1583
|
std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing
|
1a224b6a
David Mayerich
fixed linux compa...
|
1584
1585
|
|
9b2a5d71
David Mayerich
implemented finit...
|
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
|
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...
|
1607
1608
1609
|
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...
|
1610
1611
|
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...
|
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
|
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...
|
1632
|
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1633
|
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1634
1635
1636
1637
1638
|
out.write((char*)outspec, Bb); //output the band
if(PROGRESS) progress = (double)(xy+1) / (double)XY * 100;
}
}
|
278c622b
David Mayerich
finished Windows ...
|
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
|
bool multiply(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){
std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file
std::string headername = outname + ".hdr"; //the header file name
unsigned long long N = X() * Y(); //calculate the total number of pixels to be processed
unsigned long long B = Z(); //get the number of bands
T* s = (T*)malloc(sizeof(T) * B); //allocate memory to store a pixel
for(unsigned long long n = 0; n < N; n++){ //for each pixel in the image
if(mask == NULL || mask[n]){ //if the pixel is masked
for(size_t b = 0; b < B; b++) //for each band in the spectrum
s[b] *= (T)v; //multiply
}
if(PROGRESS) progress = (double)(n+1) / N * 100; //set the current progress
target.write((char*)s, sizeof(T) * B); //write the corrected data into destination
} //end for each pixel
free(s); //free the spectrum
target.close(); //close the output file
return true;
}
bool add(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){
std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file
std::string headername = outname + ".hdr"; //the header file name
unsigned long long N = X() * Y(); //calculate the total number of pixels to be processed
unsigned long long B = Z(); //get the number of bands
T* s = (T*)malloc(sizeof(T) * B); //allocate memory to store a pixel
for(unsigned long long n = 0; n < N; n++){ //for each pixel in the image
if(mask == NULL || mask[n]){ //if the pixel is masked
for(size_t b = 0; b < B; b++) //for each band in the spectrum
s[b] += (T)v; //multiply
}
if(PROGRESS) progress = (double)(n+1) / N * 100; //set the current progress
target.write((char*)s, sizeof(T) * B); //write the corrected data into destination
} //end for each pixel
free(s); //free the spectrum
target.close(); //close the output file
return true;
}
|
dc8eb8aa
David Mayerich
added band trimmi...
|
1685
|
|
0df38ff3
heziqi
fixed head detached
|
1686
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
1687
|
/// Close the file.
|
f6169dea
heziqi
Ziqi completed bi...
|
1688
1689
1690
1691
1692
|
bool close(){
file.close();
return true;
}
|
88c3e636
heziqi
Ziqi added big.h
|
1693
|
};
|
6aa04ba2
David Mayerich
interleave types
|
1694
|
}
|
e8eb202f
David Mayerich
added a new ENVI ...
|
1695
1696
|
#endif
|