6aa04ba2
David Mayerich
interleave types
|
1
2
3
|
#ifndef STIM_BIL_H
#define STIM_BIL_H
|
e8eb202f
David Mayerich
added a new ENVI ...
|
4
|
#include "../envi/envi_header.h"
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
5
|
#include "../envi/hsi.h"
|
1a224b6a
David Mayerich
fixed linux compa...
|
6
|
#include "../math/fd_coefficients.h"
|
c25e7d0d
heziqi
speed of bip.base...
|
7
8
|
#include <cstring>
#include <utility>
|
1a224b6a
David Mayerich
fixed linux compa...
|
9
|
#include <deque>
|
c25e7d0d
heziqi
speed of bip.base...
|
10
|
|
8a86bd56
David Mayerich
changed rts names...
|
11
|
namespace stim{
|
c25e7d0d
heziqi
speed of bip.base...
|
12
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
13
14
15
16
17
18
19
20
21
|
/**
The BIL class represents a 3-dimensional binary file stored using band interleaved by line (BIL) image encoding. The binary file is stored
such that X-Z "frames" are stored sequentially to form an image stack along the y-axis. When accessing the data sequentially on disk,
the dimensions read, from fastest to slowest, are X, Z, Y.
This class is optimized for data streaming, and therefore supports extremely large (terabyte-scale) files. Data is loaded from disk
on request. Functions used to access data are written to support efficient reading.
*/
|
c25e7d0d
heziqi
speed of bip.base...
|
22
23
|
template <typename T>
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
24
|
class bil: public hsi<T> {
|
c25e7d0d
heziqi
speed of bip.base...
|
25
26
|
protected:
|
c25e7d0d
heziqi
speed of bip.base...
|
27
|
|
1a224b6a
David Mayerich
fixed linux compa...
|
28
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
29
30
|
using hsi<T>::w; //use the wavelength array in stim::hsi
using hsi<T>::nnz;
|
63fc1df8
David Mayerich
added an inverse ...
|
31
|
using binary<T>::progress;
|
1a224b6a
David Mayerich
fixed linux compa...
|
32
33
34
|
using hsi<T>::X;
using hsi<T>::Y;
using hsi<T>::Z;
|
cbce396e
David Mayerich
added support for...
|
35
|
|
c25e7d0d
heziqi
speed of bip.base...
|
36
37
38
39
|
public:
using binary<T>::open;
using binary<T>::file;
|
6708cc25
heziqi
Ziqi added envi c...
|
40
|
using binary<T>::R;
|
c25e7d0d
heziqi
speed of bip.base...
|
41
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
42
43
|
bil(){ hsi<T>::init_bil(); }
|
a23c4132
David Mayerich
Doxygen comments ...
|
44
45
46
47
48
49
50
51
|
/// 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...
|
52
53
54
55
56
|
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...
|
57
|
std::vector<double> wavelengths){
|
c25e7d0d
heziqi
speed of bip.base...
|
58
|
|
6708cc25
heziqi
Ziqi added envi c...
|
59
|
w = wavelengths;
|
c25e7d0d
heziqi
speed of bip.base...
|
60
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
61
|
return open(filename, vec<unsigned long long>(X, B, Y), header_offset);
|
6708cc25
heziqi
Ziqi added envi c...
|
62
|
|
c25e7d0d
heziqi
speed of bip.base...
|
63
64
|
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
65
66
67
68
|
/// 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...
|
69
|
bool band_index( T * p, unsigned long long page, bool PROGRESS = false){
|
724ec347
David Mayerich
simplified the EN...
|
70
|
//return binary<T>::read_plane_1(p, page);
|
c25e7d0d
heziqi
speed of bip.base...
|
71
|
|
63fc1df8
David Mayerich
added an inverse ...
|
72
|
if(PROGRESS) progress = 0;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
73
74
|
unsigned long long L = X() * sizeof(T); //caculate the number of bytes in a sample line
unsigned long long jump = X() * (Z() - 1) * sizeof(T);
|
c25e7d0d
heziqi
speed of bip.base...
|
75
|
|
724ec347
David Mayerich
simplified the EN...
|
76
|
if (page >= Z()){ //make sure the bank number is right
|
c25e7d0d
heziqi
speed of bip.base...
|
77
78
79
80
|
std::cout<<"ERROR: page out of range"<<std::endl;
return false;
}
|
724ec347
David Mayerich
simplified the EN...
|
81
|
file.seekg(X() * page * sizeof(T), std::ios::beg);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
82
|
for (unsigned long long i = 0; i < Y(); i++)
|
c25e7d0d
heziqi
speed of bip.base...
|
83
|
{
|
724ec347
David Mayerich
simplified the EN...
|
84
|
file.read((char *)(p + i * X()), L);
|
c25e7d0d
heziqi
speed of bip.base...
|
85
|
file.seekg( jump, std::ios::cur);
|
63fc1df8
David Mayerich
added an inverse ...
|
86
|
if(PROGRESS) progress = (double)(i+1) / Y() * 100;
|
c25e7d0d
heziqi
speed of bip.base...
|
87
88
89
90
91
|
}
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
92
93
94
95
|
/// 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 ...
|
96
|
bool band( T * p, double wavelength, bool PROGRESS = false){
|
c25e7d0d
heziqi
speed of bip.base...
|
97
|
|
ee4dea28
David Mayerich
fixed errors in c...
|
98
99
|
//if there are no wavelengths in the BSQ file
if(w.size() == 0)
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
100
|
return band_index(p, (unsigned long long)wavelength);
|
ee4dea28
David Mayerich
fixed errors in c...
|
101
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
102
103
|
unsigned long long XY = X() * Y(); //calculate the number of pixels in a band
unsigned long long S = XY * sizeof(T); //calculate the number of bytes of a band
|
c25e7d0d
heziqi
speed of bip.base...
|
104
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
105
|
unsigned long long page=0; //bands around the wavelength
|
1a224b6a
David Mayerich
fixed linux compa...
|
106
|
|
c25e7d0d
heziqi
speed of bip.base...
|
107
108
109
110
|
//get the bands numbers around the wavelength
//if wavelength is smaller than the first one in header file
|
6708cc25
heziqi
Ziqi added envi c...
|
111
|
if ( w[page] > wavelength ){
|
63fc1df8
David Mayerich
added an inverse ...
|
112
|
band_index(p, page, PROGRESS);
|
c25e7d0d
heziqi
speed of bip.base...
|
113
114
115
|
return true;
}
|
6708cc25
heziqi
Ziqi added envi c...
|
116
|
while( w[page] < wavelength )
|
c25e7d0d
heziqi
speed of bip.base...
|
117
118
119
|
{
page++;
//if wavelength is larger than the last wavelength in header file
|
724ec347
David Mayerich
simplified the EN...
|
120
121
|
if (page == Z()) {
band_index(p, Z()-1);
|
c25e7d0d
heziqi
speed of bip.base...
|
122
123
124
|
return true;
}
}
|
6708cc25
heziqi
Ziqi added envi c...
|
125
|
if ( wavelength < w[page] )
|
c25e7d0d
heziqi
speed of bip.base...
|
126
|
{
|
517876d6
heziqi
metrics finished ...
|
127
128
|
T * p1;
T * p2;
|
c25e7d0d
heziqi
speed of bip.base...
|
129
130
131
|
p1=(T*)malloc(S); //memory allocation
p2=(T*)malloc(S);
band_index(p1, page - 1);
|
63fc1df8
David Mayerich
added an inverse ...
|
132
|
band_index(p2, page, PROGRESS);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
133
134
135
|
for(unsigned long long i=0; i < XY; i++){
double r = (wavelength - w[page-1]) / (w[page] - w[page-1]);
p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]);
|
c25e7d0d
heziqi
speed of bip.base...
|
136
|
}
|
517876d6
heziqi
metrics finished ...
|
137
138
|
free(p1);
free(p2);
|
c25e7d0d
heziqi
speed of bip.base...
|
139
140
141
|
}
else //if the wavelength is equal to a wavelength in header file
{
|
63fc1df8
David Mayerich
added an inverse ...
|
142
|
band_index(p, page, PROGRESS);
|
c25e7d0d
heziqi
speed of bip.base...
|
143
144
145
146
147
|
}
return true;
}
|
724ec347
David Mayerich
simplified the EN...
|
148
149
150
151
152
153
|
/// Retrieves a band of x values from a given xz plane.
/// @param p is a pointer to pre-allocated memory at least Z * sizeof(T) in size
/// @param c is a pointer to an existing XZ plane (size X*Z*sizeof(T))
/// @param wavelength is the wavelength to retrieve
bool read_x_from_xz(T* p, T* c, double wavelength)
|
c25e7d0d
heziqi
speed of bip.base...
|
154
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
155
156
|
unsigned long long B = Z();
unsigned long long L = X() * sizeof(T);
|
c25e7d0d
heziqi
speed of bip.base...
|
157
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
158
|
unsigned long long page=0; //samples around the wavelength
|
c25e7d0d
heziqi
speed of bip.base...
|
159
160
161
162
163
164
|
T * p1;
T * p2;
//get the bands numbers around the wavelength
//if wavelength is smaller than the first one in header file
|
6708cc25
heziqi
Ziqi added envi c...
|
165
|
if ( w[page] > wavelength ){
|
1a224b6a
David Mayerich
fixed linux compa...
|
166
|
memcpy(p, c, L);
|
c25e7d0d
heziqi
speed of bip.base...
|
167
168
169
|
return true;
}
|
6708cc25
heziqi
Ziqi added envi c...
|
170
|
while( w[page] < wavelength )
|
c25e7d0d
heziqi
speed of bip.base...
|
171
172
173
174
|
{
page++;
//if wavelength is larger than the last wavelength in header file
if (page == B) {
|
724ec347
David Mayerich
simplified the EN...
|
175
|
memcpy(p, c + (B - 1) * X(), L);
|
c25e7d0d
heziqi
speed of bip.base...
|
176
177
178
|
return true;
}
}
|
6708cc25
heziqi
Ziqi added envi c...
|
179
|
if ( wavelength < w[page] )
|
c25e7d0d
heziqi
speed of bip.base...
|
180
181
182
183
|
{
p1=(T*)malloc( L ); //memory allocation
p2=(T*)malloc( L );
|
724ec347
David Mayerich
simplified the EN...
|
184
185
|
memcpy(p1, c + (page - 1) * X(), L);
memcpy(p2, c + page * X(), L);
|
1a224b6a
David Mayerich
fixed linux compa...
|
186
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
187
|
for(unsigned long long i=0; i < X(); i++){
|
6708cc25
heziqi
Ziqi added envi c...
|
188
|
double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
189
|
p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]);
|
c25e7d0d
heziqi
speed of bip.base...
|
190
191
|
}
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
192
|
else //if the wavelength is equal to a wavelength in header file
|
724ec347
David Mayerich
simplified the EN...
|
193
|
memcpy(p, c + page * X(), L);
|
1a224b6a
David Mayerich
fixed linux compa...
|
194
195
|
return true;
|
c25e7d0d
heziqi
speed of bip.base...
|
196
197
|
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
198
199
200
201
202
|
/// Retrieve a single spectrum (B-axis line) at a given (x, y) location and stores it in pre-allocated memory.
/// @param p is a pointer to pre-allocated memory at least B * sizeof(T) in size.
/// @param x is the x-coordinate (dimension 1) of the spectrum.
/// @param y is the y-coordinate (dimension 2) of the spectrum.
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
203
|
bool spectrum(T * p, unsigned long long x, unsigned long long y, bool PROGRESS = false){
|
63fc1df8
David Mayerich
added an inverse ...
|
204
|
return binary<T>::read_line_1(p, x, y, PROGRESS);
|
c25e7d0d
heziqi
speed of bip.base...
|
205
|
}
|
bfe9c6db
heziqi
added feature_mat...
|
206
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
207
208
209
210
|
/// 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...
|
211
|
bool pixel(T * p, unsigned long long n){
|
bfe9c6db
heziqi
added feature_mat...
|
212
213
|
//calculate the corresponding x, y
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
214
215
|
unsigned long long x = n % X();
unsigned long long y = n / X();
|
bfe9c6db
heziqi
added feature_mat...
|
216
217
218
219
|
//get the pixel
return spectrum(p, x, y);
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
220
|
|
c25e7d0d
heziqi
speed of bip.base...
|
221
|
//given a Y ,return a XZ slice
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
222
|
bool read_plane_y(T * p, unsigned long long y){
|
724ec347
David Mayerich
simplified the EN...
|
223
|
return binary<T>::read_plane_2(p, y);
|
c25e7d0d
heziqi
speed of bip.base...
|
224
|
}
|
e843658b
Brad Deutsch
Previous push did...
|
225
|
|
1a224b6a
David Mayerich
fixed linux compa...
|
226
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
227
228
229
230
|
/// Perform baseline correction given a list of baseline points and stores the result in a new BSQ file.
/// @param outname is the name of the output file used to store the resulting baseline-corrected data.
/// @param wls is the list of baseline points based on band labels.
|
798cea97
David Mayerich
added progress ba...
|
231
|
bool baseline(std::string outname, std::vector<double> wls, unsigned char* mask = NULL, bool PROGRESS = false){
|
1a224b6a
David Mayerich
fixed linux compa...
|
232
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
233
|
unsigned long long N = wls.size(); //get the number of baseline points
|
1a224b6a
David Mayerich
fixed linux compa...
|
234
|
|
c25e7d0d
heziqi
speed of bip.base...
|
235
|
std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file
|
1a224b6a
David Mayerich
fixed linux compa...
|
236
237
|
std::string headername = outname + ".hdr"; //the header file name
|
c25e7d0d
heziqi
speed of bip.base...
|
238
|
//simplify image resolution
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
239
240
241
|
unsigned long long ZX = Z() * X(); //calculate the number of points in a Y slice
unsigned long long L = ZX * sizeof(T); //calculate the number of bytes of a Y slice
unsigned long long B = Z();
|
724ec347
David Mayerich
simplified the EN...
|
242
|
|
c25e7d0d
heziqi
speed of bip.base...
|
243
244
|
T* c; //pointer to the current Y slice
c = (T*)malloc(L); //memory allocation
|
1a224b6a
David Mayerich
fixed linux compa...
|
245
|
|
c25e7d0d
heziqi
speed of bip.base...
|
246
247
|
T* a; //pointer to the two YZ lines surrounding the current YZ line
T* b;
|
1a224b6a
David Mayerich
fixed linux compa...
|
248
|
|
724ec347
David Mayerich
simplified the EN...
|
249
250
|
a = (T*)malloc(X() * sizeof(T));
b = (T*)malloc(X() * sizeof(T));
|
c25e7d0d
heziqi
speed of bip.base...
|
251
252
253
254
|
double ai, bi; //stores the two baseline points wavelength surrounding the current band
double ci; //stores the current band's wavelength
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
255
|
unsigned long long control;
|
c25e7d0d
heziqi
speed of bip.base...
|
256
257
258
259
260
261
|
if (a == NULL || b == NULL || c == NULL){
std::cout<<"ERROR: error allocating memory";
exit(1);
}
// loop start correct every y slice
|
1a224b6a
David Mayerich
fixed linux compa...
|
262
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
263
|
for (unsigned long long k =0; k < Y(); k++)
|
c25e7d0d
heziqi
speed of bip.base...
|
264
265
|
{
//get the current y slice
|
724ec347
David Mayerich
simplified the EN...
|
266
|
read_plane_y(c, k);
|
1a224b6a
David Mayerich
fixed linux compa...
|
267
268
|
//initialize lownum, highnum, low, high
|
6708cc25
heziqi
Ziqi added envi c...
|
269
|
ai = w[0];
|
f6169dea
heziqi
Ziqi completed bi...
|
270
|
control=0;
|
c25e7d0d
heziqi
speed of bip.base...
|
271
272
|
//if no baseline point is specified at band 0,
//set the baseline point at band 0 to 0
|
6708cc25
heziqi
Ziqi added envi c...
|
273
|
if(wls[0] != w[0]){
|
1a224b6a
David Mayerich
fixed linux compa...
|
274
|
bi = wls[control];
|
724ec347
David Mayerich
simplified the EN...
|
275
|
memset(a, (char)0, X() * sizeof(T) );
|
c25e7d0d
heziqi
speed of bip.base...
|
276
277
278
279
|
}
//else get the low band
else{
control++;
|
724ec347
David Mayerich
simplified the EN...
|
280
|
read_x_from_xz(a, c, ai);
|
c25e7d0d
heziqi
speed of bip.base...
|
281
282
283
|
bi = wls[control];
}
//get the high band
|
724ec347
David Mayerich
simplified the EN...
|
284
|
read_x_from_xz(b, c, bi);
|
1a224b6a
David Mayerich
fixed linux compa...
|
285
|
|
c25e7d0d
heziqi
speed of bip.base...
|
286
|
//correct every YZ line
|
1a224b6a
David Mayerich
fixed linux compa...
|
287
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
288
|
for(unsigned long long cii = 0; cii < B; cii++){
|
c25e7d0d
heziqi
speed of bip.base...
|
289
290
|
//update baseline points, if necessary
|
6708cc25
heziqi
Ziqi added envi c...
|
291
|
if( w[cii] >= bi && cii != B - 1) {
|
c25e7d0d
heziqi
speed of bip.base...
|
292
293
|
//if the high band is now on the last BL point
if (control != N-1) {
|
1a224b6a
David Mayerich
fixed linux compa...
|
294
|
|
c25e7d0d
heziqi
speed of bip.base...
|
295
|
control++; //increment the index
|
1a224b6a
David Mayerich
fixed linux compa...
|
296
|
|
c25e7d0d
heziqi
speed of bip.base...
|
297
|
std::swap(a, b); //swap the baseline band pointers
|
1a224b6a
David Mayerich
fixed linux compa...
|
298
|
|
c25e7d0d
heziqi
speed of bip.base...
|
299
300
|
ai = bi;
bi = wls[control];
|
724ec347
David Mayerich
simplified the EN...
|
301
|
read_x_from_xz(b, c, bi);
|
1a224b6a
David Mayerich
fixed linux compa...
|
302
|
|
c25e7d0d
heziqi
speed of bip.base...
|
303
304
|
}
//if the last BL point on the last band of the file?
|
6708cc25
heziqi
Ziqi added envi c...
|
305
|
else if ( wls[control] < w[B - 1]) {
|
1a224b6a
David Mayerich
fixed linux compa...
|
306
|
|
c25e7d0d
heziqi
speed of bip.base...
|
307
|
std::swap(a, b); //swap the baseline band pointers
|
1a224b6a
David Mayerich
fixed linux compa...
|
308
|
|
724ec347
David Mayerich
simplified the EN...
|
309
|
memset(b, (char)0, X() * sizeof(T) ); //clear the high band
|
1a224b6a
David Mayerich
fixed linux compa...
|
310
|
|
c25e7d0d
heziqi
speed of bip.base...
|
311
|
ai = bi;
|
6708cc25
heziqi
Ziqi added envi c...
|
312
|
bi = w[B - 1];
|
c25e7d0d
heziqi
speed of bip.base...
|
313
314
|
}
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
315
|
|
6708cc25
heziqi
Ziqi added envi c...
|
316
|
ci = w[cii];
|
1a224b6a
David Mayerich
fixed linux compa...
|
317
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
318
|
unsigned long long jump = cii * X();
|
c25e7d0d
heziqi
speed of bip.base...
|
319
|
//perform the baseline correction
|
724ec347
David Mayerich
simplified the EN...
|
320
|
for(unsigned i=0; i < X(); i++)
|
c25e7d0d
heziqi
speed of bip.base...
|
321
322
|
{
double r = (double) (ci - ai) / (double) (bi - ai);
|
92e4cf05
heziqi
include file fixed
|
323
|
c[i + jump] =(T) ( c[i + jump] - (b[i] - a[i]) * r - a[i] );
|
c25e7d0d
heziqi
speed of bip.base...
|
324
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
325
326
327
|
}//loop for YZ line end
target.write(reinterpret_cast<const char*>(c), L); //write the corrected data into destination
|
cbce396e
David Mayerich
added support for...
|
328
|
|
63fc1df8
David Mayerich
added an inverse ...
|
329
|
if(PROGRESS) progress = (double)(k+1) / Y() * 100;
|
1a224b6a
David Mayerich
fixed linux compa...
|
330
331
|
}//loop for Y slice end
|
c25e7d0d
heziqi
speed of bip.base...
|
332
333
334
335
|
free(a);
free(b);
free(c);
target.close();
|
cbce396e
David Mayerich
added support for...
|
336
|
|
1a224b6a
David Mayerich
fixed linux compa...
|
337
338
|
return true;
|
c25e7d0d
heziqi
speed of bip.base...
|
339
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
340
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
341
342
343
344
345
|
/// 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...
|
346
|
bool ratio(std::string outname, double w, unsigned char* mask = NULL, bool PROGRESS = false)
|
c25e7d0d
heziqi
speed of bip.base...
|
347
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
348
349
350
|
unsigned long long B = Z(); //calculate the number of bands
unsigned long long ZX = Z() * X();
unsigned long long XY = X() * Y(); //calculate the number of pixels in a band
|
1a224b6a
David Mayerich
fixed linux compa...
|
351
|
unsigned long long S = XY * sizeof(T); //calculate the number of bytes in a band
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
352
|
unsigned long long L = ZX * sizeof(T);
|
c25e7d0d
heziqi
speed of bip.base...
|
353
354
355
356
357
358
359
360
|
std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file
std::string headername = outname + ".hdr"; //the header file name
T * c; //pointer to the current ZX slice
T * b; //pointer to the standard band
b = (T*)malloc( S ); //memory allocation
|
1a224b6a
David Mayerich
fixed linux compa...
|
361
|
c = (T*)malloc( L );
|
c25e7d0d
heziqi
speed of bip.base...
|
362
|
|
e933c14e
David Mayerich
cleaned up the code
|
363
|
band(b, w); //get the certain band into memory
|
c25e7d0d
heziqi
speed of bip.base...
|
364
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
365
|
for(unsigned long long j = 0; j < Y(); j++)
|
c25e7d0d
heziqi
speed of bip.base...
|
366
|
{
|
724ec347
David Mayerich
simplified the EN...
|
367
|
read_plane_y(c, j);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
368
|
for(unsigned long long i = 0; i < B; i++)
|
c25e7d0d
heziqi
speed of bip.base...
|
369
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
370
|
for(unsigned long long m = 0; m < X(); m++)
|
c25e7d0d
heziqi
speed of bip.base...
|
371
|
{
|
798cea97
David Mayerich
added progress ba...
|
372
|
if( mask != NULL && !mask[m + j * X()] )
|
724ec347
David Mayerich
simplified the EN...
|
373
|
c[m + i * X()] = (T)0.0;
|
93503ec6
David Mayerich
allowed masking d...
|
374
|
else
|
724ec347
David Mayerich
simplified the EN...
|
375
|
c[m + i * X()] = c[m + i * X()] / b[m + j * X()];
|
1a224b6a
David Mayerich
fixed linux compa...
|
376
|
}
|
c25e7d0d
heziqi
speed of bip.base...
|
377
378
|
}
target.write(reinterpret_cast<const char*>(c), L); //write normalized data into destination
|
cbce396e
David Mayerich
added support for...
|
379
|
|
63fc1df8
David Mayerich
added an inverse ...
|
380
|
if(PROGRESS) progress = (double)(j+1) / Y() * 100;
|
c25e7d0d
heziqi
speed of bip.base...
|
381
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
382
|
|
c25e7d0d
heziqi
speed of bip.base...
|
383
384
385
|
free(b);
free(c);
target.close();
|
c25e7d0d
heziqi
speed of bip.base...
|
386
387
|
return true;
}
|
c30dd3e6
David Mayerich
added vector norm...
|
388
|
|
904614c3
David Mayerich
added normalizati...
|
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
|
void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){
std::cout<<"ERROR: algorithm not implemented"<<std::endl;
exit(1);
//This code is almost done but has to be debugged
/*std::ofstream target(outfile.c_str(), std::ios::binary); //open the target binary file
size_t B = Z(); //calculate the number of bands
size_t L = X(); //get the number of items in a line
size_t Lb = L * sizeof(T); //calculate the number of bytes in a line
size_t XY = X() * Y(); //calculate the number of pixels in an image
T* line = (T*) malloc(Lb); //allocate space for a line
T* len = (T*) malloc(Lb); //allocate space for the lengths in a line
for(size_t y = 0; y < Y(); y++){ //for each line in the image
file.seekg(y * Lb, std::ios::beg); //move the pointer to the current file to the beginning
for(size_t b = 0; b < B; b++){ //for each band in this line
file.read((char*)line, Lb); //read a line of pixels
for(size_t l = 0; l < L; l++) //for each element in the line for this band
len[l] += line[l] * line[l]; //add the square of the spectral component
}
for(size_t l = 0; l < L; l++) //for each length element in the line
len[l] = sqrt(len[l]); //calculate the square root of the sum of squares
file.seekg(y * Lb, std::ios::beg); //move the pointer to the current file to the beginning
for(size_t b = 0; b < B; b++){ //for each band in this line
file.read((char*)line, Lb); //read a line of pixels
for(size_t l = 0; l < L; l++) //for each element in the line for this band
line[l] /= len[l]; //divide each element by the length
}
target.write((char*)line, Lb); //write the normalized line to the target file
if(PROGRESS) progress = (double)(y+1) / (double)Y();
}*/
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
423
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
424
425
426
|
/// Convert the current BIL file to a BSQ file with the specified file name.
/// @param outname is the name of the output BSQ file to be saved to disk.
|
63fc1df8
David Mayerich
added an inverse ...
|
427
|
bool bsq(std::string outname, bool PROGRESS = false)
|
c25e7d0d
heziqi
speed of bip.base...
|
428
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
429
|
unsigned long long S = X() * Y() * sizeof(T); //calculate the number of bytes in a band
|
1a224b6a
David Mayerich
fixed linux compa...
|
430
|
|
c25e7d0d
heziqi
speed of bip.base...
|
431
432
|
std::ofstream target(outname.c_str(), std::ios::binary);
std::string headername = outname + ".hdr";
|
1a224b6a
David Mayerich
fixed linux compa...
|
433
|
|
c25e7d0d
heziqi
speed of bip.base...
|
434
435
|
T * p; //pointer to the current band
p = (T*)malloc(S);
|
1a224b6a
David Mayerich
fixed linux compa...
|
436
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
437
|
for ( unsigned long long i = 0; i < Z(); i++)
|
1a224b6a
David Mayerich
fixed linux compa...
|
438
|
{
|
c25e7d0d
heziqi
speed of bip.base...
|
439
|
band_index(p, i);
|
cbce396e
David Mayerich
added support for...
|
440
441
|
target.write(reinterpret_cast<const char*>(p), S); //write a band data into target file
|
1a224b6a
David Mayerich
fixed linux compa...
|
442
|
if(PROGRESS) progress = (double)(i+1) / Z() * 100; //store the progress for the current operation
|
c25e7d0d
heziqi
speed of bip.base...
|
443
|
}
|
6708cc25
heziqi
Ziqi added envi c...
|
444
|
|
c25e7d0d
heziqi
speed of bip.base...
|
445
446
447
448
449
|
free(p);
target.close();
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
450
451
452
|
/// Convert the current BIL file to a BIP file with the specified file name.
/// @param outname is the name of the output BIP file to be saved to disk.
|
63fc1df8
David Mayerich
added an inverse ...
|
453
|
bool bip(std::string outname, bool PROGRESS = false)
|
f6169dea
heziqi
Ziqi completed bi...
|
454
|
{
|
c1078e19
David Mayerich
HSIproc bug fixes...
|
455
|
unsigned long long S = X() * Z() * sizeof(T); //calculate the number of bytes in a ZX slice
|
1a224b6a
David Mayerich
fixed linux compa...
|
456
|
|
f6169dea
heziqi
Ziqi completed bi...
|
457
458
|
std::ofstream target(outname.c_str(), std::ios::binary);
std::string headername = outname + ".hdr";
|
1a224b6a
David Mayerich
fixed linux compa...
|
459
|
|
f6169dea
heziqi
Ziqi completed bi...
|
460
461
462
463
|
T * p; //pointer to the current XZ slice for bil file
p = (T*)malloc(S);
T * q; //pointer to the current ZX slice for bip file
q = (T*)malloc(S);
|
1a224b6a
David Mayerich
fixed linux compa...
|
464
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
465
|
for ( unsigned long long i = 0; i < Y(); i++)
|
1a224b6a
David Mayerich
fixed linux compa...
|
466
|
{
|
724ec347
David Mayerich
simplified the EN...
|
467
|
read_plane_y(p, i);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
468
|
for ( unsigned long long k = 0; k < Z(); k++)
|
f6169dea
heziqi
Ziqi completed bi...
|
469
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
470
471
|
unsigned long long ks = k * X();
for ( unsigned long long j = 0; j < X(); j++)
|
724ec347
David Mayerich
simplified the EN...
|
472
|
q[k + j * Z()] = p[ks + j];
|
cbce396e
David Mayerich
added support for...
|
473
|
|
1a224b6a
David Mayerich
fixed linux compa...
|
474
|
if(PROGRESS) progress = (double)((i+1) * Z() + k+1) / (Z() * Y()) * 100; //store the progress for the current operation
|
f6169dea
heziqi
Ziqi completed bi...
|
475
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
476
477
|
target.write(reinterpret_cast<const char*>(q), S); //write a band data into target file
|
f6169dea
heziqi
Ziqi completed bi...
|
478
479
|
}
|
f6169dea
heziqi
Ziqi completed bi...
|
480
481
482
483
484
485
|
free(p);
free(q);
target.close();
return true;
}
|
20c212c0
heziqi
Ziqi added functi...
|
486
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
487
488
489
490
491
492
493
494
|
/// 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...
|
495
496
|
bool baseline_band(double lb, double rb, T* lp, T* rp, double wavelength, T* result){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
497
|
unsigned long long XY = X() * Y();
|
20c212c0
heziqi
Ziqi added functi...
|
498
499
500
501
|
band(result, wavelength); //get band
//perform the baseline correction
double r = (double) (wavelength - lb) / (double) (rb - lb);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
502
|
for(unsigned long long i=0; i < XY; i++){
|
20c212c0
heziqi
Ziqi added functi...
|
503
504
505
506
|
result[i] =(T) (result[i] - (rp[i] - lp[i]) * r - lp[i] );
}
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
507
508
509
510
511
512
513
|
/// 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...
|
514
|
bool height(double lb, double rb, double bandwavelength, T* result){
|
20c212c0
heziqi
Ziqi added functi...
|
515
516
517
|
T* lp;
T* rp;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
518
519
|
unsigned long long XY = X() * Y();
unsigned long long S = XY * sizeof(T);
|
70407ea9
heziqi
Ziqi completed he...
|
520
521
|
lp = (T*) malloc(S); //memory allocation
rp = (T*) malloc(S);
|
20c212c0
heziqi
Ziqi added functi...
|
522
|
|
70407ea9
heziqi
Ziqi completed he...
|
523
|
band(lp, lb);
|
1a224b6a
David Mayerich
fixed linux compa...
|
524
|
band(rp, rb);
|
20c212c0
heziqi
Ziqi added functi...
|
525
526
527
|
baseline_band(lb, rb, lp, rp, bandwavelength, result);
|
517876d6
heziqi
metrics finished ...
|
528
529
|
free(lp);
free(rp);
|
20c212c0
heziqi
Ziqi added functi...
|
530
531
532
|
return true;
}
|
70407ea9
heziqi
Ziqi completed he...
|
533
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
534
535
536
537
538
539
540
541
|
/// 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...
|
542
|
bool area(double lb, double rb, double lab, double rab, T* result){
|
20c212c0
heziqi
Ziqi added functi...
|
543
544
545
546
547
|
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...
|
548
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
549
550
|
unsigned long long XY = X() * Y();
unsigned long long S = XY * sizeof(T);
|
20c212c0
heziqi
Ziqi added functi...
|
551
552
553
554
555
|
lp = (T*) malloc(S); //memory allocation
rp = (T*) malloc(S);
cur = (T*) malloc(S);
cur2 = (T*) malloc(S);
|
20c212c0
heziqi
Ziqi added functi...
|
556
557
558
559
|
memset(result, (char)0, S);
//find the wavelenght position in the whole band
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
560
561
562
|
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...
|
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
|
//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...
|
578
|
while (lab >= w[ai]){
|
20c212c0
heziqi
Ziqi added functi...
|
579
580
|
ai++;
}
|
70407ea9
heziqi
Ziqi completed he...
|
581
|
while (rab <= w[bi]){
|
20c212c0
heziqi
Ziqi added functi...
|
582
583
584
585
586
587
588
|
bi--;
}
band(lp, lb);
band(rp, rb);
//calculate the beginning and the ending part
|
70407ea9
heziqi
Ziqi completed he...
|
589
590
|
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...
|
591
592
|
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...
|
593
|
}
|
70407ea9
heziqi
Ziqi completed he...
|
594
595
|
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...
|
596
|
for(unsigned long long j = 0; j < XY; j++){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
597
|
result[j] += (T)((w[ai] - lab) * ((double)cur[j] + (double)cur2[j]) / 2.0);
|
20c212c0
heziqi
Ziqi added functi...
|
598
|
}
|
517876d6
heziqi
metrics finished ...
|
599
|
|
20c212c0
heziqi
Ziqi added functi...
|
600
601
|
//calculate the area
ai++;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
602
|
for(unsigned long long i = ai; i <= bi ;i++)
|
20c212c0
heziqi
Ziqi added functi...
|
603
604
|
{
baseline_band(lb, rb, lp, rp, w[ai], cur2);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
605
|
for(unsigned long long j = 0; j < XY; j++)
|
20c212c0
heziqi
Ziqi added functi...
|
606
|
{
|
1a224b6a
David Mayerich
fixed linux compa...
|
607
|
result[j] += (T)((w[ai] - w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 2.0);
|
20c212c0
heziqi
Ziqi added functi...
|
608
609
610
|
}
std::swap(cur,cur2); //swap the band pointers
}
|
70407ea9
heziqi
Ziqi completed he...
|
611
|
|
517876d6
heziqi
metrics finished ...
|
612
613
614
615
|
free(lp);
free(rp);
free(cur);
free(cur2);
|
20c212c0
heziqi
Ziqi added functi...
|
616
617
618
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
619
620
621
622
623
624
625
626
627
|
/// 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...
|
628
|
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...
|
629
|
|
724ec347
David Mayerich
simplified the EN...
|
630
631
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
1a224b6a
David Mayerich
fixed linux compa...
|
632
|
|
70407ea9
heziqi
Ziqi completed he...
|
633
634
635
636
|
//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...
|
637
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
70407ea9
heziqi
Ziqi completed he...
|
638
639
640
641
642
|
if(p1[i] == 0 && p2[i] ==0)
result[i] = 1;
else
result[i] = p1[i] / p2[i];
}
|
517876d6
heziqi
metrics finished ...
|
643
644
645
|
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
646
|
return true;
|
70407ea9
heziqi
Ziqi completed he...
|
647
|
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
648
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
649
650
651
652
653
654
655
656
657
|
/// 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...
|
658
659
|
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...
|
660
|
|
724ec347
David Mayerich
simplified the EN...
|
661
662
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
1a224b6a
David Mayerich
fixed linux compa...
|
663
|
|
70407ea9
heziqi
Ziqi completed he...
|
664
665
666
667
|
//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...
|
668
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
70407ea9
heziqi
Ziqi completed he...
|
669
670
671
672
673
|
if(p1[i] == 0 && p2[i] ==0)
result[i] = 1;
else
result[i] = p1[i] / p2[i];
}
|
517876d6
heziqi
metrics finished ...
|
674
675
676
|
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
677
678
679
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
680
681
682
683
684
685
686
687
688
|
/// 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...
|
689
|
/// @param rab2 is the label value for the right bound (end of the integration) of the second peak (denominator)
|
a23c4132
David Mayerich
Doxygen comments ...
|
690
|
/// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size
|
ba51ae6a
David Mayerich
fixed metric calc...
|
691
692
|
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...
|
693
|
|
724ec347
David Mayerich
simplified the EN...
|
694
695
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
1a224b6a
David Mayerich
fixed linux compa...
|
696
|
|
70407ea9
heziqi
Ziqi completed he...
|
697
698
699
700
|
//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...
|
701
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
70407ea9
heziqi
Ziqi completed he...
|
702
703
704
705
706
|
if(p1[i] == 0 && p2[i] ==0)
result[i] = 1;
else
result[i] = p1[i] / p2[i];
}
|
517876d6
heziqi
metrics finished ...
|
707
708
709
|
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
710
711
|
return true;
}
|
70407ea9
heziqi
Ziqi completed he...
|
712
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
713
714
715
716
717
718
719
|
/// 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 ...
|
720
721
722
723
724
725
|
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...
|
726
727
|
unsigned long long XY = X() * Y();
unsigned long long S = XY * sizeof(T);
|
517876d6
heziqi
metrics finished ...
|
728
729
730
731
732
733
734
735
736
|
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...
|
737
738
739
|
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 ...
|
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
|
//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...
|
766
767
|
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 ...
|
768
769
770
|
}
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...
|
771
|
for(unsigned long long j = 0; j < XY; j++){
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
772
|
result[j] += (T)((w[ai] - lab) * (w[ai] + lab) * ((double)cur[j] + (double)cur2[j]) / 4.0);
|
517876d6
heziqi
metrics finished ...
|
773
774
775
776
|
}
//calculate f(x) times x
ai++;
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
777
|
for(unsigned long long i = ai; i <= bi ;i++)
|
517876d6
heziqi
metrics finished ...
|
778
779
|
{
baseline_band(lb, rb, lp, rp, w[ai], cur2);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
780
|
for(unsigned long long j = 0; j < XY; j++)
|
517876d6
heziqi
metrics finished ...
|
781
|
{
|
1a224b6a
David Mayerich
fixed linux compa...
|
782
|
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 ...
|
783
784
785
786
787
788
789
790
791
792
793
|
}
std::swap(cur,cur2); //swap the band pointers
}
free(lp);
free(rp);
free(cur);
free(cur2);
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
794
795
796
797
798
799
800
|
/// 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...
|
801
|
bool centroid(T* result, double lb, double rb, double lab, double rab, unsigned char* mask = NULL){
|
724ec347
David Mayerich
simplified the EN...
|
802
803
|
T* p1 = (T*)malloc(X() * Y() * sizeof(T));
T* p2 = (T*)malloc(X() * Y() * sizeof(T));
|
1a224b6a
David Mayerich
fixed linux compa...
|
804
|
|
517876d6
heziqi
metrics finished ...
|
805
806
807
808
|
//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...
|
809
|
for(unsigned long long i = 0; i < X() * Y(); i++){
|
ba51ae6a
David Mayerich
fixed metric calc...
|
810
|
if(mask == NULL || mask[i])
|
517876d6
heziqi
metrics finished ...
|
811
812
813
814
815
|
result[i] = p1[i] / p2[i];
}
free(p1);
free(p2);
|
1a224b6a
David Mayerich
fixed linux compa...
|
816
|
return true;
|
517876d6
heziqi
metrics finished ...
|
817
818
|
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
819
820
821
822
823
824
825
|
/// Create a mask based on a given band and threshold value.
/// All pixels in the
/// specified band greater than the threshold are true and all pixels less than the threshold are false.
/// @param mask_band is the band used to specify the mask
/// @param threshold is the threshold used to determine if the mask value is true or false
/// @param p is a pointer to a pre-allocated array at least X * Y in size
|
ba51ae6a
David Mayerich
fixed metric calc...
|
826
|
bool build_mask(unsigned char* mask, double mask_band, double threshold, bool PROGRESS = false){
|
4a6f666c
heziqi
Added mask method
|
827
|
|
724ec347
David Mayerich
simplified the EN...
|
828
|
T* temp = (T*)malloc(X() * Y() * sizeof(T)); //allocate memory for the certain band
|
798cea97
David Mayerich
added progress ba...
|
829
|
band(temp, mask_band, PROGRESS);
|
4a6f666c
heziqi
Added mask method
|
830
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
831
|
for (unsigned long long i = 0; i < X() * Y(); i++) {
|
798cea97
David Mayerich
added progress ba...
|
832
|
if (temp[i] < threshold)
|
ba51ae6a
David Mayerich
fixed metric calc...
|
833
|
mask[i] = 0;
|
798cea97
David Mayerich
added progress ba...
|
834
|
else
|
ba51ae6a
David Mayerich
fixed metric calc...
|
835
|
mask[i] = 255;
|
4a6f666c
heziqi
Added mask method
|
836
837
|
}
|
517876d6
heziqi
metrics finished ...
|
838
|
free(temp);
|
4a6f666c
heziqi
Added mask method
|
839
840
841
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
842
843
844
845
|
/// 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...
|
846
|
bool apply_mask(std::string outfile, unsigned char* p, bool PROGRESS = false){
|
740f8cd2
heziqi
added apply_mask
|
847
848
849
|
std::ofstream target(outfile.c_str(), std::ios::binary);
|
1a55a328
David Mayerich
Added comments fo...
|
850
|
//I THINK THIS IS WRONG
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
851
852
|
unsigned long long XZ = X() * Z(); //calculate the number of values in a page on disk
unsigned long long L = XZ * sizeof(T); //calculate the size of the page (in bytes)
|
740f8cd2
heziqi
added apply_mask
|
853
|
|
1a55a328
David Mayerich
Added comments fo...
|
854
|
T * temp = (T*)malloc(L); //allocate memory for a temporary page
|
740f8cd2
heziqi
added apply_mask
|
855
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
856
|
for (unsigned long long i = 0; i < Y(); i++) //for each value in Y() (BIP should be X)
|
740f8cd2
heziqi
added apply_mask
|
857
|
{
|
724ec347
David Mayerich
simplified the EN...
|
858
|
read_plane_y(temp, i); //retrieve an ZX slice, stored in temp
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
859
|
for ( unsigned long long j = 0; j < Z(); j++) //for each Z() (Y)
|
740f8cd2
heziqi
added apply_mask
|
860
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
861
|
for (unsigned long long k = 0; k < X(); k++) //for each band
|
740f8cd2
heziqi
added apply_mask
|
862
|
{
|
724ec347
David Mayerich
simplified the EN...
|
863
864
|
if(p[i * X() + k] == 0)
temp[j * X() + k] = 0;
|
740f8cd2
heziqi
added apply_mask
|
865
866
867
868
|
else
continue;
}
}
|
1a224b6a
David Mayerich
fixed linux compa...
|
869
|
target.write(reinterpret_cast<const char*>(temp), L); //write a band data into target file
|
6a46c8ff
David Mayerich
fixed bugs in app...
|
870
|
if(PROGRESS) progress = (double)(i+1) / (double)Y() * 100;
|
740f8cd2
heziqi
added apply_mask
|
871
872
873
874
875
876
|
}
target.close();
free(temp);
return true;
}
|
ba51ae6a
David Mayerich
fixed metric calc...
|
877
878
879
880
881
882
883
884
885
886
|
/// Copies all spectra corresponding to nonzero values of a mask into a pre-allocated matrix of size (B x P)
/// where P is the number of masked pixels and B is the number of bands. The allocated memory can be accessed
/// using the following indexing: i = p*B + b
/// @param matrix is the destination for the pixel data
/// @param mask is the mask
bool sift(T* matrix, unsigned char* mask){
size_t Lbytes = sizeof(T) * X();
T* line = (T*) malloc( Lbytes ); //allocate space for a line
file.seekg(0, std::ios::beg); //seek to the beginning of the file
|
1a224b6a
David Mayerich
fixed linux compa...
|
887
|
|
ba51ae6a
David Mayerich
fixed metric calc...
|
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
|
size_t pl;
size_t p = 0; //create counter variables
for(size_t y = 0; y < Y(); y++){ //for each line in the data set
for(size_t b = 0; b < Z(); b++){ //for each band in the data set
pl = 0; //initialize the pixel offset for the current line to zero (0)
file.read( (char*)line, Lbytes ); //read the current line
for(size_t x = 0; x < X(); x++){
if(mask[y * X() + x]){ //if the current pixel is masked
size_t i = (p + pl) * Z() + b; //calculate the index in the sifted matrix
matrix[i] = line[x]; //store the current value in the line at the correct matrix location
pl++; //increment the pixel pointer
}
}
}
p += pl; //add the line increment to the running pixel index
}
return true;
}
|
f7e99517
David Mayerich
added sifting for...
|
907
908
|
/// Saves to disk only those spectra corresponding to mask values != 0
/// Unlike the BIP and BSQ versions of this function, this version saves a different format: the BIL file is saved as a BIP
|
c1078e19
David Mayerich
HSIproc bug fixes...
|
909
|
bool sift(std::string outfile, unsigned char* p, bool PROGRESS = false){
|
724ec347
David Mayerich
simplified the EN...
|
910
|
// Assume X() = X, Y() = Y, Z() = Z.
|
e843658b
Brad Deutsch
Previous push did...
|
911
912
913
|
std::ofstream target(outfile.c_str(), std::ios::binary);
//for loading pages:
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
914
915
916
|
unsigned long long XZ = X() * Z(); //calculate the number of values in an XZ page on disk
unsigned long long B = Z(); //calculate the number of bands
unsigned long long L = XZ * sizeof(T); //calculate the size of the page (in bytes)
|
e843658b
Brad Deutsch
Previous push did...
|
917
|
|
f7e99517
David Mayerich
added sifting for...
|
918
919
920
921
922
|
//allocate temporary memory for a XZ slice
T* slice = (T*) malloc(L);
//allocate temporary memory for a spectrum
T* spec = (T*) malloc(B * sizeof(T));
|
e843658b
Brad Deutsch
Previous push did...
|
923
|
|
f7e99517
David Mayerich
added sifting for...
|
924
|
//for each slice along the y axis
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
925
|
for (unsigned long long y = 0; y < Y(); y++) //Select a page by choosing Y coordinate, Y()
|
e843658b
Brad Deutsch
Previous push did...
|
926
|
{
|
f7e99517
David Mayerich
added sifting for...
|
927
928
|
read_plane_y(slice, y); //retrieve an ZX page, store in "slice"
|
1a224b6a
David Mayerich
fixed linux compa...
|
929
|
//for each sample along X
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
930
|
for (unsigned long long x = 0; x < X(); x++) //Select a pixel by choosing X coordinate in the page, X()
|
e843658b
Brad Deutsch
Previous push did...
|
931
|
{
|
f7e99517
David Mayerich
added sifting for...
|
932
933
|
//if the mask != 0 at that xy pixel
if (p[y * X() + x] != 0) //if the mask != 0 at that XY pixel
|
e843658b
Brad Deutsch
Previous push did...
|
934
|
{
|
f7e99517
David Mayerich
added sifting for...
|
935
|
//for each band at that pixel
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
936
|
for (unsigned long long b = 0; b < B; b++) //Select a voxel by choosing Z coordinate at the pixel
|
e843658b
Brad Deutsch
Previous push did...
|
937
|
{
|
f7e99517
David Mayerich
added sifting for...
|
938
|
spec[b] = slice[b*X() + x]; //Pass the correct spectral value from XZ page into the spectrum to be saved.
|
e843658b
Brad Deutsch
Previous push did...
|
939
|
}
|
f7e99517
David Mayerich
added sifting for...
|
940
|
target.write((char*)spec, B * sizeof(T)); //write that spectrum to disk. Size is L2.
|
e843658b
Brad Deutsch
Previous push did...
|
941
|
}
|
3150c537
David Mayerich
added further pro...
|
942
|
|
c1078e19
David Mayerich
HSIproc bug fixes...
|
943
|
if(PROGRESS) progress = (double) ((y+1) * X() + x+1) / (Y() * X()) * 100;
|
e843658b
Brad Deutsch
Previous push did...
|
944
945
946
|
}
}
target.close();
|
f7e99517
David Mayerich
added sifting for...
|
947
948
|
free(slice);
free(spec);
|
3150c537
David Mayerich
added further pro...
|
949
|
|
e843658b
Brad Deutsch
Previous push did...
|
950
951
952
|
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
953
954
955
|
/// Calculate the mean band value (average along B) at each pixel location.
/// @param p is a pointer to memory of size X * Y * sizeof(T) that will store the band averages.
|
a43c4fe1
heziqi
Added crop in env...
|
956
|
bool band_avg(T* p){
|
724ec347
David Mayerich
simplified the EN...
|
957
|
unsigned long long XZ = X() * Z();
|
a43c4fe1
heziqi
Added crop in env...
|
958
|
T* temp = (T*)malloc(sizeof(T) * XZ);
|
724ec347
David Mayerich
simplified the EN...
|
959
|
T* line = (T*)malloc(sizeof(T) * X());
|
a43c4fe1
heziqi
Added crop in env...
|
960
|
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
961
|
for (unsigned long long i = 0; i < Y(); i++){
|
a43c4fe1
heziqi
Added crop in env...
|
962
963
|
getY(temp, i);
//initialize x-line
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
964
|
for (unsigned long long j = 0; j < X(); j++){
|
a43c4fe1
heziqi
Added crop in env...
|
965
966
|
line[j] = 0;
}
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
967
968
969
|
unsigned long long c = 0;
for (unsigned long long j = 0; j < Z(); j++){
for (unsigned long long k = 0; k < X(); k++){
|
724ec347
David Mayerich
simplified the EN...
|
970
|
line[k] += temp[c] / (T)Z();
|
a43c4fe1
heziqi
Added crop in env...
|
971
972
973
|
c++;
}
}
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
974
|
for (unsigned long long j = 0; j < X(); j++){
|
724ec347
David Mayerich
simplified the EN...
|
975
|
p[j + i * X()] = line[j];
|
a43c4fe1
heziqi
Added crop in env...
|
976
977
978
979
980
981
|
}
}
free(temp);
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
982
983
984
985
|
/// 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 ...
|
986
|
bool avg_band(double* p, unsigned char* mask = NULL, bool PROGRESS = false){
|
724ec347
David Mayerich
simplified the EN...
|
987
988
|
unsigned long long XZ = X() * Z();
unsigned long long XY = X() * Y();
|
a43c4fe1
heziqi
Added crop in env...
|
989
|
T* temp = (T*)malloc(sizeof(T) * XZ);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
990
|
for (unsigned long long j = 0; j < Z(); j++){
|
a43c4fe1
heziqi
Added crop in env...
|
991
992
993
|
p[j] = 0;
}
//calculate vaild number in a band
|
2aa68315
David Mayerich
major bug fixes f...
|
994
995
996
|
unsigned long long count = 0;
for (unsigned long long j = 0; j < XY; j++){
if (mask == NULL || mask[j] != 0){
|
a43c4fe1
heziqi
Added crop in env...
|
997
998
999
|
count++;
}
}
|
2aa68315
David Mayerich
major bug fixes f...
|
1000
|
for (unsigned long long k = 0; k < Y(); k++){
|
724ec347
David Mayerich
simplified the EN...
|
1001
|
read_plane_y(temp, k);
|
2aa68315
David Mayerich
major bug fixes f...
|
1002
1003
1004
1005
1006
|
unsigned long long kx = k * X();
for (unsigned long long i = 0; i < X(); i++){
if (mask == NULL || mask[kx + i] != 0){
for (unsigned long long j = 0; j < Z(); j++){
p[j] += temp[j * X() + i] / (double)count;
|
a43c4fe1
heziqi
Added crop in env...
|
1007
1008
1009
|
}
}
}
|
63fc1df8
David Mayerich
added an inverse ...
|
1010
|
if(PROGRESS) progress = (double)(k+1) / Y() * 100;
|
a43c4fe1
heziqi
Added crop in env...
|
1011
1012
1013
1014
1015
|
}
free(temp);
return true;
}
|
a23c4132
David Mayerich
Doxygen comments ...
|
1016
1017
1018
1019
1020
|
/// Calculate the covariance matrix for all masked pixels in the image.
/// @param co is a pointer to pre-allocated memory of size [B * B] that stores the resulting covariance matrix
/// @param avg is a pointer to memory of size B that stores the average spectrum
/// @param mask is a pointer to memory of size [X * Y] that stores the mask value at each pixel location
|
63fc1df8
David Mayerich
added an inverse ...
|
1021
1022
|
bool co_matrix(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){
progress = 0;
|
a43c4fe1
heziqi
Added crop in env...
|
1023
|
//memory allocation
|
724ec347
David Mayerich
simplified the EN...
|
1024
|
unsigned long long xy = X() * Y();
|
2aa68315
David Mayerich
major bug fixes f...
|
1025
|
unsigned long long B = Z();
|
a43c4fe1
heziqi
Added crop in env...
|
1026
1027
|
T* temp = (T*)malloc(sizeof(T) * B);
//count vaild pixels in a band
|
2aa68315
David Mayerich
major bug fixes f...
|
1028
1029
1030
|
unsigned long long count = 0;
for (unsigned long long j = 0; j < xy; j++){
if (mask == NULL || mask[j] != 0){
|
a43c4fe1
heziqi
Added crop in env...
|
1031
1032
1033
1034
|
count++;
}
}
//initialize correlation matrix
|
2aa68315
David Mayerich
major bug fixes f...
|
1035
1036
|
for (unsigned long long i = 0; i < B; i++){
for (unsigned long long k = 0; k < B; k++){
|
a43c4fe1
heziqi
Added crop in env...
|
1037
1038
1039
1040
|
co[i * B + k] = 0;
}
}
//calculate correlation coefficient matrix
|
2aa68315
David Mayerich
major bug fixes f...
|
1041
1042
|
for (unsigned long long j = 0; j < xy; j++){
if (mask == NULL || mask[j] != 0){
|
a43c4fe1
heziqi
Added crop in env...
|
1043
|
pixel(temp, j);
|
2aa68315
David Mayerich
major bug fixes f...
|
1044
1045
1046
|
for (unsigned long long i = 0; i < B; i++){
for (unsigned long long k = i; k < B; k++){
co[i * B + k] += ((double)temp[i] - (double)avg[i]) * ((double)temp[k] - (double)avg[k]) / (double)count;
|
a43c4fe1
heziqi
Added crop in env...
|
1047
1048
1049
|
}
}
}
|
63fc1df8
David Mayerich
added an inverse ...
|
1050
|
if(PROGRESS) progress = (double)(j+1) / xy * 100;
|
a43c4fe1
heziqi
Added crop in env...
|
1051
1052
|
}
//because correlation matrix is symmetric
|
2aa68315
David Mayerich
major bug fixes f...
|
1053
1054
|
for (unsigned long long i = 0; i < B; i++){
for (unsigned long long k = i + 1; k < B; k++){
|
a43c4fe1
heziqi
Added crop in env...
|
1055
1056
1057
1058
1059
|
co[k * B + i] = co[i * B + k];
}
}
free(temp);
|
a43c4fe1
heziqi
Added crop in env...
|
1060
1061
1062
|
return true;
}
|
0df38ff3
heziqi
fixed head detached
|
1063
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
1064
1065
1066
1067
1068
1069
1070
|
/// 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...
|
1071
1072
1073
1074
1075
|
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 ...
|
1076
1077
|
unsigned long long b1,
bool PROGRESS = false){
|
a43c4fe1
heziqi
Added crop in env...
|
1078
|
|
3033b886
David Mayerich
implemented spect...
|
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
|
//calculate the new image parameters
unsigned long long samples = x1 - x0;
unsigned long long lines = y1 - y0;
unsigned long long bands = b1 - b0;
//calculate the size of a line
unsigned long long L = samples * sizeof(T);
//allocate space for a line
T* temp = (T*)malloc(bands * L);
//create an output stream to store the cropped file
|
a43c4fe1
heziqi
Added crop in env...
|
1092
|
std::ofstream out(outfile.c_str(), std::ios::binary);
|
3033b886
David Mayerich
implemented spect...
|
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
|
//calculate the distance between bands
unsigned long long jumpb = (X() - samples) * sizeof(T);
//distance needed to jump from the previous line of the last band to the next line of the first band
unsigned long long longjump = ((Z() - b1) * X() + b0 * X()) * sizeof(T);
//set the start position for the cropped region
file.seekg((y0 * X() * Z() + b0 * X() + x0) * sizeof(T), std::ios::beg);
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1103
|
for (unsigned long long x = 0; x < lines; x++)
|
a43c4fe1
heziqi
Added crop in env...
|
1104
|
{
|
9d3ba0b1
David Mayerich
added stim::hsi a...
|
1105
|
for (unsigned long long z = b0; z < b1; z++)
|
a43c4fe1
heziqi
Added crop in env...
|
1106
|
{
|
3033b886
David Mayerich
implemented spect...
|
1107
|
file.read((char *)(temp + z * samples), sizeof(T) * samples);
|
1a224b6a
David Mayerich
fixed linux compa...
|
1108
|
file.seekg(jumpb, std::ios::cur); //go to the next band
|
3150c537
David Mayerich
added further pro...
|
1109
|
|
63fc1df8
David Mayerich
added an inverse ...
|
1110
|
if(PROGRESS) progress = (double)(x * Z() + z+1) / (lines * Z()) * 100;
|
a43c4fe1
heziqi
Added crop in env...
|
1111
|
}
|
3033b886
David Mayerich
implemented spect...
|
1112
1113
1114
1115
1116
1117
|
//write slice data into target file
out.write(reinterpret_cast<const char*>(temp), bands * L);
//seek to the beginning of the next X-Z slice
file.seekg(longjump, std::ios::cur);
|
a43c4fe1
heziqi
Added crop in env...
|
1118
|
}
|
3033b886
David Mayerich
implemented spect...
|
1119
1120
|
//free the temporary frame
|
a43c4fe1
heziqi
Added crop in env...
|
1121
|
free(temp);
|
3150c537
David Mayerich
added further pro...
|
1122
|
|
a43c4fe1
heziqi
Added crop in env...
|
1123
1124
1125
|
return true;
}
|
dc8eb8aa
David Mayerich
added band trimmi...
|
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
|
/// Remove a list of bands from the ENVI file
/// @param outfile is the file name for the output hyperspectral image (with trimmed bands)
/// @param b is an array of bands to be eliminated
void trim(std::string outfile, std::vector<size_t> band_array, bool PROGRESS = false){
std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing
file.seekg(0, std::ios::beg); //move to the beginning of the input file
size_t Xb = X() * sizeof(T); //calculate the number of bytes in a line
T* line = (T*)malloc(Xb); //allocate space for a line
size_t i; //create an index into the band array
for(size_t y = 0; y < Y(); y++){ //for each Y plane
i = 0;
for(size_t b = 0; b < Z(); b++){ //for each band
if(b != band_array[i]){ //if this band isn't trimmed
file.read((char*)line, Xb); //read a line
out.write((char*)line, Xb); //write the line
}
else{
file.seekg(Xb, std::ios::cur); //if this band is trimmed, skip it
i++;
}
}
if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100;
}
free(line);
}
|
2ce6954b
David Mayerich
added the ability...
|
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
|
/// Combine two BSQ images along the Y axis
/// @param outfile is the combined file to be output
/// @param infile is the input file stream for the image to combine with this one
/// @param Sx is the size of the second image along X
/// @param Sy is the size of the second image along Y
/// @param offset is a shift (negative or positive) in the combined image to the left or right
void combine(std::string outfile, bil<T>* C, long long xp, long long yp, bool PROGRESS = false){
std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing
file.seekg(0, std::ios::beg); //move to the beginning of both files
C->file.seekg(0, std::ios::beg);
size_t S[2]; //size of the output band image
size_t p0[2]; //position of the current image in the output
size_t p1[2]; //position of the source image in the output
hsi<T>::calc_combined_size(xp, yp, C->X(), C->Y(), S[0], S[1], p0[0], p0[1], p1[0], p1[1]); //calculate the image placement parameters
size_t line_bytes = X() * sizeof(T);
size_t band_bytes = X() * Y() * sizeof(T);
T* cur = (T*)malloc(line_bytes); //allocate space for a band of the current image
size_t line_src_bytes = C->X() * sizeof(T);
size_t band_src_bytes = C->X() * C->Y() * sizeof(T);
T* src = (T*)malloc(line_src_bytes); //allocate space for a band of the source image
size_t line_dst_bytes = S[0] * sizeof(T);
size_t band_dst_bytes = S[0] * S[1] * sizeof(T);
T* dst = (T*)malloc(line_dst_bytes); //allocate space for a band of the destination image
for(size_t y = 0; y < S[1]; y++){ //for each line in the destination file
memset(dst, 0, line_dst_bytes); //set all values to zero (0) in the destination image
for(size_t b = 0; b < Z(); b++){ //for each band in both images
if(y >= p0[1] && y < p0[1] + Y()){ //if the current image crosses this line
file.read((char*)cur, line_bytes); //read a line from the current image
memcpy(&dst[p0[0]], cur, line_bytes); //copy the current line to the correct spot in the destination line
}
if(y >= p1[1] && y < p1[1] + C->Y()){ //if the source image crosses this line
C->file.read((char*)src, line_src_bytes); //read a line from the source image
memcpy(&dst[p1[0]], src, line_src_bytes); //copy the source line into the correct spot in the destination line
}
out.write((char*)dst, line_dst_bytes);
if(PROGRESS) progress = (double)((y + 1)*Z() + b + 1) / (double) (Z() * S[1]) * 100;
}
}
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1202
|
/// Convolve the given band range with a kernel specified by a vector of coefficients.
|
2ce6954b
David Mayerich
added the ability...
|
1203
|
|
9b2a5d71
David Mayerich
implemented finit...
|
1204
1205
1206
1207
1208
1209
|
/// @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...
|
1210
|
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...
|
1211
1212
1213
1214
1215
1216
1217
1218
|
std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing
size_t B = end - start + 1;
size_t Xb = X() * sizeof(T); //calculate the size of a band (frame) in bytes
size_t XBb = Xb * Z();
size_t N = C.size(); //get the number of bands that the kernel spans
std::deque<T*> frame(N, NULL); //create an array to store pointers to each frame
for(size_t f = 0; f < N; f++)
|
1a224b6a
David Mayerich
fixed linux compa...
|
1219
|
frame[f] = (T*)malloc(Xb); //allocate space for the frame
|
9b2a5d71
David Mayerich
implemented finit...
|
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
|
T* outline = (T*)malloc(Xb); //allocate space for the output band
//Implementation: In order to minimize reads from secondary storage, each band is only loaded once into the 'frame' deque.
// When a new band is loaded, the last band is popped, a new frame is copied to the pointer, and it is
// re-inserted into the deque.
for(size_t y = 0; y < Y(); y++){
file.seekg(y * XBb + start * Xb, std::ios::beg); //move to the beginning of the 'start' band
for(size_t f = 0; f < N; f++) //for each frame
file.read((char*)frame[f], Xb); //load the frame
for(size_t b = 0; b < B; b++){ //for each band
memset(outline, 0, Xb); //set the output band to zero (0)
for(size_t c = 0; c < N; c++){ //for each frame (corresponding to each coefficient)
for(size_t x = 0; x < X(); x++){ //for each pixel
|
9d77bbd9
David Mayerich
updates for HSIvi...
|
1234
1235
1236
|
if(mask == NULL || mask[y * X() + x]){
outline[x] += (T)(C[c] * frame[c][x]); //calculate the contribution of the current frame (scaled by the corresponding coefficient)
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1237
1238
1239
1240
1241
1242
1243
1244
|
}
}
out.write((char*)outline, Xb); //output the band
file.read((char*)frame[0], Xb); //read the next band
frame.push_back(frame.front()); //put the first element in the back
frame.pop_front(); //pop the first element
}
if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100;
|
1a224b6a
David Mayerich
fixed linux compa...
|
1245
|
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1246
1247
1248
|
}
/// Approximate the spectral derivative of the image
|
9d77bbd9
David Mayerich
updates for HSIvi...
|
1249
|
void deriv(std::string outfile, size_t d, size_t order, const std::vector<double> w = std::vector<double>(), unsigned char* mask = NULL, bool PROGRESS = false){
|
9b2a5d71
David Mayerich
implemented finit...
|
1250
|
std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing
|
1a224b6a
David Mayerich
fixed linux compa...
|
1251
|
|
9b2a5d71
David Mayerich
implemented finit...
|
1252
1253
1254
1255
|
size_t Xb = X() * sizeof(T); //calculate the size of a line (frame) in bytes
size_t XBb = Xb * Z();
size_t B = Z();
|
1a224b6a
David Mayerich
fixed linux compa...
|
1256
|
file.seekg(0, std::ios::beg); //seek to the beginning of the file
|
9b2a5d71
David Mayerich
implemented finit...
|
1257
1258
1259
1260
|
size_t N = order + d; //approximating a derivative requires order + d samples
std::deque<T*> frame(N, NULL); //create an array to store pointers to each frame
for(size_t f = 0; f < N; f++) //for each frame
|
1a224b6a
David Mayerich
fixed linux compa...
|
1261
|
frame[f] = (T*)malloc(Xb); //allocate space for the frame
|
9b2a5d71
David Mayerich
implemented finit...
|
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
|
T* outline = (T*)malloc(Xb); //allocate space for the output band
//Implementation: In order to minimize reads from secondary storage, each band is only loaded once into the 'frame' deque.
// When a new band is loaded, the last band is popped, a new frame is copied to the pointer, and it is
// re-inserted into the deque.
size_t mid = (size_t)(N / 2); //calculate the mid point of the kernel
size_t iw; //index to the first wavelength used to evaluate the derivative at this band
for(size_t y = 0; y < Y(); y++){
//file.seekg(y * XBb, std::ios::beg); //seek to the beginning of the current Y plane
for(size_t f = 0; f < N; f++) //for each frame
file.read((char*)frame[f], Xb); //load a line
for(size_t b = 0; b < B; b++){ //for each band
if(b < mid) //calculate the first wavelength used to evaluate the derivative at this band
iw = 0;
else if(b > B - (N - mid + 1))
iw = B - N;
else{
iw = b - mid;
file.read((char*)frame[0], Xb); //read the next band
frame.push_back(frame.front()); //put the first element in the back
frame.pop_front(); //pop the first element
}
std::vector<double> w_pts(w.begin() + iw, w.begin() + iw + N); //get the wavelengths corresponding to each sample
std::vector<double> C = diff_coefficients(w[b], w_pts, d); //get the optimal sample weights
memset(outline, 0, Xb); //set the output band to zero (0)
for(size_t c = 0; c < N; c++){ //for each frame (corresponding to each coefficient)
for(size_t x = 0; x < X(); x++){ //for each pixel
|
9d77bbd9
David Mayerich
updates for HSIvi...
|
1293
1294
1295
|
if(mask == NULL || mask[y * X() + x]){
outline[x] += (T)(C[c] * frame[c][x]); //calculate the contribution of the current frame (scaled by the corresponding coefficient)
}
|
9b2a5d71
David Mayerich
implemented finit...
|
1296
1297
1298
|
}
}
out.write((char*)outline, Xb); //output the band
|
1a224b6a
David Mayerich
fixed linux compa...
|
1299
|
|
9b2a5d71
David Mayerich
implemented finit...
|
1300
1301
1302
1303
1304
|
}
if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100;
}
}
|
0df38ff3
heziqi
fixed head detached
|
1305
|
|
a23c4132
David Mayerich
Doxygen comments ...
|
1306
|
/// Close the file.
|
f6169dea
heziqi
Ziqi completed bi...
|
1307
1308
1309
1310
1311
|
bool close(){
file.close();
return true;
}
|
c25e7d0d
heziqi
speed of bip.base...
|
1312
|
};
|
cac62fd3
David Mayerich
modified to work ...
|
1313
|
}
|
6aa04ba2
David Mayerich
interleave types
|
1314
1315
|
#endif
|