#ifndef STIM_BIL_H #define STIM_BIL_H #include "../envi/envi_header.h" #include "../envi/hsi.h" #include "../math/fd_coefficients.h" #include #include #include namespace stim{ /** 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. */ template class bil: public hsi { protected: using hsi::w; //use the wavelength array in stim::hsi using hsi::nnz; using binary::progress; using hsi::X; using hsi::Y; using hsi::Z; public: using binary::open; using binary::file; using binary::R; bil(){ hsi::init_bil(); } /// 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 bool open(std::string filename, unsigned long long X, unsigned long long Y, unsigned long long B, unsigned long long header_offset, std::vector wavelengths){ w = wavelengths; return open(filename, vec(X, B, Y), header_offset); } /// 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. bool band_index( T * p, unsigned long long page, bool PROGRESS = false){ //return binary::read_plane_1(p, page); if(PROGRESS) progress = 0; 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); if (page >= Z()){ //make sure the bank number is right std::cout<<"ERROR: page out of range"< wavelength ){ band_index(p, page, PROGRESS); return true; } while( w[page] < wavelength ) { page++; //if wavelength is larger than the last wavelength in header file if (page == Z()) { band_index(p, Z()-1); return true; } } if ( wavelength < w[page] ) { T * p1; T * p2; p1=(T*)malloc(S); //memory allocation p2=(T*)malloc(S); band_index(p1, page - 1); band_index(p2, page, PROGRESS); 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]); } free(p1); free(p2); } else //if the wavelength is equal to a wavelength in header file { band_index(p, page, PROGRESS); } return true; } /// 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) { unsigned long long B = Z(); unsigned long long L = X() * sizeof(T); unsigned long long page=0; //samples around the wavelength T * p1; T * p2; //get the bands numbers around the wavelength //if wavelength is smaller than the first one in header file if ( w[page] > wavelength ){ memcpy(p, c, L); return true; } while( w[page] < wavelength ) { page++; //if wavelength is larger than the last wavelength in header file if (page == B) { memcpy(p, c + (B - 1) * X(), L); return true; } } if ( wavelength < w[page] ) { p1=(T*)malloc( L ); //memory allocation p2=(T*)malloc( L ); memcpy(p1, c + (page - 1) * X(), L); memcpy(p2, c + page * X(), L); for(unsigned long long i=0; i < X(); i++){ double r = (double) (wavelength - w[page-1]) / (double) (w[page] - w[page-1]); p[i] = (T)(((double)p2[i] - (double)p1[i]) * r + (double)p1[i]); } } else //if the wavelength is equal to a wavelength in header file memcpy(p, c + page * X(), L); return true; } /// 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. bool spectrum(T* p, size_t n, bool PROGRESS = false){ size_t y = n / X(); size_t x = n - y * X(); return binary::read_line_1(p, x, y, PROGRESS); } bool spectrum(T * p, unsigned long long x, unsigned long long y, bool PROGRESS = false){ return binary::read_line_1(p, x, y, PROGRESS); } /// 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. bool pixel(T * p, unsigned long long n){ //calculate the corresponding x, y unsigned long long x = n % X(); unsigned long long y = n / X(); //get the pixel return spectrum(p, x, y); } //given a Y ,return a XZ slice bool read_plane_y(T * p, unsigned long long y){ return binary::read_plane_2(p, y); } /// 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. bool baseline(std::string outname, std::vector wls, unsigned char* mask = NULL, bool PROGRESS = false){ unsigned long long N = wls.size(); //get the number of baseline points std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file std::string headername = outname + ".hdr"; //the header file name //simplify image resolution 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(); T* c; //pointer to the current Y slice c = (T*)malloc(L); //memory allocation T* a; //pointer to the two YZ lines surrounding the current YZ line T* b; a = (T*)malloc(X() * sizeof(T)); b = (T*)malloc(X() * sizeof(T)); double ai, bi; //stores the two baseline points wavelength surrounding the current band double ci; //stores the current band's wavelength unsigned long long control; if (a == NULL || b == NULL || c == NULL){ std::cout<<"ERROR: error allocating memory"; exit(1); } // loop start correct every y slice for (unsigned long long k =0; k < Y(); k++) { //get the current y slice read_plane_y(c, k); //initialize lownum, highnum, low, high ai = w[0]; control=0; //if no baseline point is specified at band 0, //set the baseline point at band 0 to 0 if(wls[0] != w[0]){ bi = wls[control]; memset(a, (char)0, X() * sizeof(T) ); } //else get the low band else{ control++; read_x_from_xz(a, c, ai); bi = wls[control]; } //get the high band read_x_from_xz(b, c, bi); //correct every YZ line for(unsigned long long cii = 0; cii < B; cii++){ //update baseline points, if necessary if( w[cii] >= bi && cii != B - 1) { //if the high band is now on the last BL point if (control != N-1) { control++; //increment the index std::swap(a, b); //swap the baseline band pointers ai = bi; bi = wls[control]; read_x_from_xz(b, c, bi); } //if the last BL point on the last band of the file? else if ( wls[control] < w[B - 1]) { std::swap(a, b); //swap the baseline band pointers memset(b, (char)0, X() * sizeof(T) ); //clear the high band ai = bi; bi = w[B - 1]; } } ci = w[cii]; unsigned long long jump = cii * X(); //perform the baseline correction for(unsigned i=0; i < X(); i++) { double r = (double) (ci - ai) / (double) (bi - ai); c[i + jump] =(T) ( c[i + jump] - (b[i] - a[i]) * r - a[i] ); } }//loop for YZ line end target.write(reinterpret_cast(c), L); //write the corrected data into destination if(PROGRESS) progress = (double)(k+1) / Y() * 100; }//loop for Y slice end free(a); free(b); free(c); target.close(); return true; } /// 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. bool ratio(std::string outname, double w, unsigned char* mask = NULL, bool PROGRESS = false) { 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 unsigned long long S = XY * sizeof(T); //calculate the number of bytes in a band unsigned long long L = ZX * sizeof(T); 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 c = (T*)malloc( L ); band(b, w); //get the certain band into memory for(unsigned long long j = 0; j < Y(); j++) { read_plane_y(c, j); for(unsigned long long i = 0; i < B; i++) { for(unsigned long long m = 0; m < X(); m++) { if( mask != NULL && !mask[m + j * X()] ) c[m + i * X()] = (T)0.0; else c[m + i * X()] = c[m + i * X()] / b[m + j * X()]; } } target.write(reinterpret_cast(c), L); //write normalized data into destination if(PROGRESS) progress = (double)(j+1) / Y() * 100; } free(b); free(c); target.close(); return true; } void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){ std::cout<<"ERROR: algorithm not implemented"<(p), S); //write a band data into target file if(PROGRESS) progress = (double)(i+1) / Z() * 100; //store the progress for the current operation } free(p); target.close(); return true; } /// 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. bool bip(std::string outname, bool PROGRESS = false) { unsigned long long S = X() * Z() * sizeof(T); //calculate the number of bytes in a ZX slice std::ofstream target(outname.c_str(), std::ios::binary); std::string headername = outname + ".hdr"; 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); for ( unsigned long long i = 0; i < Y(); i++) { read_plane_y(p, i); for ( unsigned long long k = 0; k < Z(); k++) { unsigned long long ks = k * X(); for ( unsigned long long j = 0; j < X(); j++) q[k + j * Z()] = p[ks + j]; if(PROGRESS) progress = (double)((i+1) * Z() + k+1) / (Z() * Y()) * 100; //store the progress for the current operation } target.write(reinterpret_cast(q), S); //write a band data into target file } free(p); free(q); target.close(); return true; } /// 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. bool baseline_band(double lb, double rb, T* lp, T* rp, double wavelength, T* result){ unsigned long long XY = X() * Y(); band(result, wavelength); //get band //perform the baseline correction double r = (double) (wavelength - lb) / (double) (rb - lb); for(unsigned long long i=0; i < XY; i++){ result[i] =(T) (result[i] - (rp[i] - lp[i]) * r - lp[i] ); } return true; } /// 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. bool height(double lb, double rb, double bandwavelength, T* result){ T* lp; T* rp; unsigned long long XY = X() * Y(); unsigned long long S = XY * sizeof(T); lp = (T*) malloc(S); //memory allocation rp = (T*) malloc(S); band(lp, lb); band(rp, rb); baseline_band(lb, rb, lp, rp, bandwavelength, result); free(lp); free(rp); return true; } /// 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 bool 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 unsigned long long XY = X() * Y(); unsigned long long S = XY * sizeof(T); 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 unsigned long long n = w.size(); unsigned long long ai = 0; //left bound position unsigned long long bi = n - 1; //right bound position //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"< rb){ std::cout<<"ERROR: right bound should be bigger than left bound"<= 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); for(unsigned long long j = 0; j < XY; j++){ result[j] += (T)((rab - w[bi]) * ((double)cur[j] + (double)cur2[j]) / 2.0); } baseline_band(lb, rb, lp, rp, lab, cur2); //beginnning part baseline_band(lb, rb, lp, rp, w[ai], cur); for(unsigned long long j = 0; j < XY; j++){ result[j] += (T)((w[ai] - lab) * ((double)cur[j] + (double)cur2[j]) / 2.0); } //calculate the area ai++; for(unsigned long long i = ai; i <= bi ;i++) { baseline_band(lb, rb, lp, rp, w[ai], cur2); for(unsigned long long j = 0; j < XY; j++) { result[j] += (T)((w[ai] - w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 2.0); } std::swap(cur,cur2); //swap the band pointers } free(lp); free(rp); free(cur); free(cur2); return true; } /// 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 bool ph_to_ph(T* result, double lb1, double rb1, double pos1, double lb2, double rb2, double pos2, unsigned char* mask = NULL){ T* p1 = (T*)malloc(X() * Y() * sizeof(T)); T* p2 = (T*)malloc(X() * Y() * sizeof(T)); //get the two peak band height(lb1, rb1, pos1, p1); height(lb2, rb2, pos2, p2); //calculate the ratio in result for(unsigned long long i = 0; i < X() * Y(); i++){ if(p1[i] == 0 && p2[i] ==0) result[i] = 1; else result[i] = p1[i] / p2[i]; } free(p1); free(p2); return true; } /// 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 bool pa_to_ph(T* result, double lb1, double rb1, double lab1, double rab1, double lb2, double rb2, double pos, unsigned char* mask = NULL){ T* p1 = (T*)malloc(X() * Y() * sizeof(T)); T* p2 = (T*)malloc(X() * Y() * sizeof(T)); //get the area and the peak band area(lb1, rb1, lab1, rab1, p1); height(lb2, rb2, pos, p2); //calculate the ratio in result for(unsigned long long i = 0; i < X() * Y(); i++){ if(p1[i] == 0 && p2[i] ==0) result[i] = 1; else result[i] = p1[i] / p2[i]; } free(p1); free(p2); return true; } /// Compute the ratio between two peak areas. /// @param lb1 is the label value for the left baseline point for the first peak (numerator) /// @param rb1 is the label value for the right baseline point for the first peak (numerator) /// @param lab1 is the label value for the left bound (start of the integration) of the first peak (numerator) /// @param rab1 is the label value for the right bound (end of the integration) of the first peak (numerator) /// @param lb2 is the label value for the left baseline point for the second peak (denominator) /// @param rb2 is the label value for the right baseline point for the second peak (denominator) /// @param lab2 is the label value for the left bound (start of the integration) of the second peak (denominator) /// @param rab2 is the label value for the right bound (end of the integration) of the second peak (denominator) /// @param result is a pointer to a pre-allocated array at least X * Y * sizeof(T) in size 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){ T* p1 = (T*)malloc(X() * Y() * sizeof(T)); T* p2 = (T*)malloc(X() * Y() * sizeof(T)); //get the area and the peak band area(lb1, rb1, lab1, rab1, p1); area(lb2, rb2, lab2, rab2, p2); //calculate the ratio in result for(unsigned long long i = 0; i < X() * Y(); i++){ if(p1[i] == 0 && p2[i] ==0) result[i] = 1; else result[i] = p1[i] / p2[i]; } free(p1); free(p2); return true; } /// 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 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 unsigned long long XY = X() * Y(); unsigned long long S = XY * sizeof(T); 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 unsigned long long n = w.size(); unsigned long long ai = 0; //left bound position unsigned long long bi = n - 1; //right bound position //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"< rb){ std::cout<<"ERROR: right bound should be bigger than left bound"<= 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); 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); } baseline_band(lb, rb, lp, rp, lab, cur2); //beginnning part baseline_band(lb, rb, lp, rp, w[ai], cur); for(unsigned long long j = 0; j < XY; j++){ result[j] += (T)((w[ai] - lab) * (w[ai] + lab) * ((double)cur[j] + (double)cur2[j]) / 4.0); } //calculate f(x) times x ai++; for(unsigned long long i = ai; i <= bi ;i++) { baseline_band(lb, rb, lp, rp, w[ai], cur2); for(unsigned long long j = 0; j < XY; j++) { result[j] += (T)((w[ai] - w[ai-1]) * (w[ai] + w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 4.0); } std::swap(cur,cur2); //swap the band pointers } free(lp); free(rp); free(cur); free(cur2); return true; } /// 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 bool centroid(T* result, double lb, double rb, double lab, double rab, unsigned char* mask = NULL){ T* p1 = (T*)malloc(X() * Y() * sizeof(T)); T* p2 = (T*)malloc(X() * Y() * sizeof(T)); //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 for(unsigned long long i = 0; i < X() * Y(); i++){ if(mask == NULL || mask[i]) result[i] = p1[i] / p2[i]; } free(p1); free(p2); return true; } /// 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 bool build_mask(unsigned char* mask, double mask_band, double threshold, bool PROGRESS = false){ T* temp = (T*)malloc(X() * Y() * sizeof(T)); //allocate memory for the certain band band(temp, mask_band, PROGRESS); for (unsigned long long i = 0; i < X() * Y(); i++) { if (temp[i] < threshold) mask[i] = 0; else mask[i] = 255; } free(temp); return true; } /// 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. bool apply_mask(std::string outfile, unsigned char* p, bool PROGRESS = false){ std::ofstream target(outfile.c_str(), std::ios::binary); //I THINK THIS IS WRONG 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) T * temp = (T*)malloc(L); //allocate memory for a temporary page for (unsigned long long i = 0; i < Y(); i++) //for each value in Y() (BIP should be X) { read_plane_y(temp, i); //retrieve an ZX slice, stored in temp for ( unsigned long long j = 0; j < Z(); j++) //for each Z() (Y) { for (unsigned long long k = 0; k < X(); k++) //for each band { if(p[i * X() + k] == 0) temp[j * X() + k] = 0; else continue; } } target.write(reinterpret_cast(temp), L); //write a band data into target file if(PROGRESS) progress = (double)(i+1) / (double)Y() * 100; } target.close(); free(temp); return true; } /// 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 = NULL, bool PROGRESS = false){ 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 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 == NULL || 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 } } if(PROGRESS) progress = (double)( (y+1)*Z() + 1) / (double)(Y() * Z()) * 100; } p += pl; //add the line increment to the running pixel index } return true; } /// 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 bool sift(std::string outfile, unsigned char* p, bool PROGRESS = false){ // Assume X() = X, Y() = Y, Z() = Z. std::ofstream target(outfile.c_str(), std::ios::binary); //for loading pages: 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) //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)); //for each slice along the y axis for (unsigned long long y = 0; y < Y(); y++) //Select a page by choosing Y coordinate, Y() { read_plane_y(slice, y); //retrieve an ZX page, store in "slice" //for each sample along X for (unsigned long long x = 0; x < X(); x++) //Select a pixel by choosing X coordinate in the page, X() { //if the mask != 0 at that xy pixel if (p[y * X() + x] != 0) //if the mask != 0 at that XY pixel { //for each band at that pixel for (unsigned long long b = 0; b < B; b++) //Select a voxel by choosing Z coordinate at the pixel { spec[b] = slice[b*X() + x]; //Pass the correct spectral value from XZ page into the spectrum to be saved. } target.write((char*)spec, B * sizeof(T)); //write that spectrum to disk. Size is L2. } if(PROGRESS) progress = (double) ((y+1) * X() + x+1) / (Y() * X()) * 100; } } target.close(); free(slice); free(spec); return true; } /// 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. bool band_avg(T* p){ unsigned long long XZ = X() * Z(); T* temp = (T*)malloc(sizeof(T) * XZ); T* line = (T*)malloc(sizeof(T) * X()); for (unsigned long long i = 0; i < Y(); i++){ getY(temp, i); //initialize x-line for (unsigned long long j = 0; j < X(); j++){ line[j] = 0; } unsigned long long c = 0; for (unsigned long long j = 0; j < Z(); j++){ for (unsigned long long k = 0; k < X(); k++){ line[k] += temp[c] / (T)Z(); c++; } } for (unsigned long long j = 0; j < X(); j++){ p[j + i * X()] = line[j]; } } free(temp); return true; } /// 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 bool avg_band(double* p, unsigned char* mask = NULL, bool PROGRESS = false){ unsigned long long XZ = X() * Z(); unsigned long long XY = X() * Y(); T* temp = (T*)malloc(sizeof(T) * XZ); for (unsigned long long j = 0; j < Z(); j++){ p[j] = 0; } //calculate vaild number in a band unsigned long long count = 0; for (unsigned long long j = 0; j < XY; j++){ if (mask == NULL || mask[j] != 0){ count++; } } for (unsigned long long k = 0; k < Y(); k++){ read_plane_y(temp, k); 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; } } } if(PROGRESS) progress = (double)(k+1) / Y() * 100; } free(temp); return true; } /// 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 bool co_matrix(double* co, double* avg, unsigned char *mask, bool PROGRESS = false){ progress = 0; //memory allocation unsigned long long xy = X() * Y(); unsigned long long B = Z(); T* temp = (T*)malloc(sizeof(T) * B); //count vaild pixels in a band unsigned long long count = 0; for (unsigned long long j = 0; j < xy; j++){ if (mask == NULL || mask[j] != 0){ count++; } } //initialize correlation matrix for (unsigned long long i = 0; i < B; i++){ for (unsigned long long k = 0; k < B; k++){ co[i * B + k] = 0; } } //calculate correlation coefficient matrix for (unsigned long long j = 0; j < xy; j++){ if (mask == NULL || mask[j] != 0){ pixel(temp, j); 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; } } } if(PROGRESS) progress = (double)(j+1) / xy * 100; } //because correlation matrix is symmetric for (unsigned long long i = 0; i < B; i++){ for (unsigned long long k = i + 1; k < B; k++){ co[k * B + i] = co[i * B + k]; } } free(temp); return true; } /// 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 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, unsigned long long b1, bool PROGRESS = false){ //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 std::ofstream out(outfile.c_str(), std::ios::binary); //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); for (unsigned long long x = 0; x < lines; x++) { for (unsigned long long z = b0; z < b1; z++) { file.read((char *)(temp + z * samples), sizeof(T) * samples); file.seekg(jumpb, std::ios::cur); //go to the next band if(PROGRESS) progress = (double)(x * Z() + z+1) / (lines * Z()) * 100; } //write slice data into target file out.write(reinterpret_cast(temp), bands * L); //seek to the beginning of the next X-Z slice file.seekg(longjump, std::ios::cur); } //free the temporary frame free(temp); return true; } /// 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 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); } /// 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* 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::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; } } } /// 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) void convolve(std::string outfile, std::vector C, size_t start, size_t end, unsigned char* mask = NULL, bool PROGRESS = false){ 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 frame(N, NULL); //create an array to store pointers to each frame for(size_t f = 0; f < N; f++) frame[f] = (T*)malloc(Xb); //allocate space for the frame 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 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) } } } 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; } } /// Approximate the spectral derivative of the image void deriv(std::string outfile, size_t d, size_t order, const std::vector w = std::vector(), unsigned char* mask = NULL, bool PROGRESS = false){ std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file for writing size_t Xb = X() * sizeof(T); //calculate the size of a line (frame) in bytes size_t XBb = Xb * Z(); size_t B = Z(); file.seekg(0, std::ios::beg); //seek to the beginning of the file size_t N = order + d; //approximating a derivative requires order + d samples std::deque frame(N, NULL); //create an array to store pointers to each frame for(size_t f = 0; f < N; f++) //for each frame frame[f] = (T*)malloc(Xb); //allocate space for the frame 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 w_pts(w.begin() + iw, w.begin() + iw + N); //get the wavelengths corresponding to each sample std::vector 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 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) } } } out.write((char*)outline, Xb); //output the band } if(PROGRESS) progress = (double)(y+1) / (double)Y() * 100; } } /// Close the file. bool close(){ file.close(); return true; } }; } #endif