#ifndef STIM_BSQ_H #define STIM_BSQ_H #include #include #include #include #include #include #include #include #include #include namespace stim{ /** The BIP class represents a 3-dimensional binary file stored using band sequential (BSQ) image encoding. The binary file is stored such that X-Y "frames" are stored sequentially to form an image stack along the z-axis. When accessing the data sequentially on disk, the dimensions read, from fastest to slowest, are X, Y, Z. 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 bsq: public hsi { protected: //std::vector w; //band wavelengths unsigned long long offset; using binary::R; 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::read_line_2; using binary::read_plane_2; bsq(){ hsi::init_bsq(); } /// 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 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, stim::iotype io = stim::io_in){ //copy the wavelengths to the BSQ file structure w = wavelengths; //copy the wavelengths to the structure offset = header_offset; return open(filename, vec(X, Y, B), header_offset, io); } /// 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){ return read_plane_2(p, page); //call the binary read_plane function (don't let it update the progress) } /// 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. bool band( T * p, double wavelength, bool PROGRESS = false){ if(PROGRESS) progress = 0; //if there are no wavelengths in the BSQ file if(w.size() == 0){ band_index(p, (unsigned long long)wavelength); if(PROGRESS) progress = 100; return true; } unsigned long long XY = X() * Y(); //calculate the number of pixels in a band unsigned long long page = 0; //get the two neighboring bands (above and below 'wavelength') //if wavelength is smaller than the first one in header file if ( w[page] > wavelength ){ band_index(p, page); if(PROGRESS) progress = 100; return true; } while( w[page] < wavelength ) { page++; //if wavelength is larger than the last wavelength in the header file // (the wavelength is out of bounds) if (page == Z()) { band_index(p, Z()-1); //return the last band if(PROGRESS) progress = 100; return true; } } //when the page counter points to the first band above 'wavelength' if ( wavelength < w[page] ){ //do the interpolation T * p1; T * p2; p1=(T*)malloc( XY * sizeof(T)); //memory allocation p2=(T*)malloc( XY * sizeof(T)); band_index(p1, page - 1); band_index(p2, page ); 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); } //if the wavelength is equal to a wavelength in header file else{ band_index(p, page); //return the band } if(PROGRESS) progress = 100; return true; } /// Retrieve a single spectrum (Z-axis line) at a given (x, y) location and stores it in pre-allocated memory. /// @param p is a pointer to pre-allocated memory at least B * sizeof(T) in size. /// @param x is the x-coordinate (dimension 1) of the spectrum. /// @param y is the y-coordinate (dimension 2) of the spectrum. void spectrum(T* p, size_t n, bool PROGRESS){ read_line_2(p, n, PROGRESS); } void spectrum(T * p, unsigned long long x, unsigned long long y, bool PROGRESS = false){ read_line_2(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){ unsigned long long bandnum = X() * Y(); //calculate numbers in one band if ( n >= bandnum){ //make sure the pixel number is right std::cout<<"ERROR: sample or line out of range"<::header, std::ios::beg); //point to the certain pixel for (unsigned long long i = 0; i < Z(); i++) { file.read((char *)(p + i), sizeof(T)); file.seekg((bandnum - 1) * sizeof(T), std::ios::cur); //go to the next band } return true; } /// 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 ) { size_t 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 B = Z(); //calculate the number of bands 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 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=0; T * a; //pointers to the high and low band images T * b; T * c; //pointer to the current image a = (T*)malloc( S ); //memory allocation b = (T*)malloc( S ); c = (T*)malloc( S ); if (a == NULL || b == NULL || c == NULL){ std::cout<<"ERROR: error allocating memory"; exit(1); } //initialize lownum, highnum, low, high ai=w[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, S); } //else get the low band else{ control += 1; band(a, ai); bi = wls[control]; } //get the high band band(b, bi); //correct every band 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]; band(b, 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, S); //clear the high band ai = bi; bi = w[B - 1]; } } //get the current band band_index(c, cii); ci = w[cii]; //perform the baseline correction for(unsigned long long i=0; i < XY; i++){ if(mask != NULL && !mask[i]) //if the pixel is excluded by a mask c[i] = 0; //set the value to zero else{ double r = (double) (ci - ai) / (double) (bi - ai); c[i] =(T) ( c[i] - (b[i] - a[i]) * r - a[i] ); } } target.write(reinterpret_cast(c), S); //write the corrected data into destination if(PROGRESS)progress = (double)(cii+1) / B * 100; } //header.save(headername); //save the new header file 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 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 std::ofstream target(outname.c_str(), std::ios::binary); //open the target binary file std::string headername = outname + ".hdr"; //the header file name T * b; //pointers to the certain wavelength band T * c; //pointer to the current image b = (T*)malloc( S ); //memory allocation c = (T*)malloc( S ); band(b, w); //get the certain band into memory for(unsigned long long j = 0; j < B; j++) { band_index(c, j); //get the current band into memory for(unsigned long long i = 0; i < XY; i++) { if(mask != NULL && !mask[i]) c[i] = (T)0.0; else c[i] = c[i] / b[i]; } target.write(reinterpret_cast(c), S); //write normalized data into destination if(PROGRESS) progress = (double)(j+1) / B * 100; } //header.save(headername); //save the new header file free(b); free(c); target.close(); return true; } void normalize(std::string outfile, unsigned char* mask = NULL, bool PROGRESS = false){ size_t B = Z(); //calculate the number of bands size_t XY = X() * Y(); //calculate the number of pixels in a band size_t XYb = XY * sizeof(T); //calculate the size of a band in bytes std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file file.seekg(0, std::ios::beg); //move the pointer to the current file to the beginning T* len = (T*)malloc(XYb); //allocate space to store the vector length memset(len, 0, XYb); //initialize the vector length to zero (0) T* band = (T*) malloc(XYb); //allocate space to store a band image for(size_t b = 0; b < B; b++){ file.read((char*)band, XYb); for(size_t xy = 0; xy < XY; xy++){ if(mask == NULL || mask[xy]){ len[xy] += pow(band[xy], 2); //sum the squared value for each pixel value in the band } } if(PROGRESS) progress = (double) (b+1) / (double)B * 50; } for(size_t xy = 0; xy < XY; xy++){ //for each pixel, calculate the square root if(mask == NULL || mask[xy]){ len[xy] += pow(band[xy], 2); //sum the squared value for each pixel value in the band } } file.seekg(0, std::ios::beg); //move the pointer to the current file to the beginning for(size_t b = 0; b < B; b++){ file.read((char*)band, XYb); for(size_t xy = 0; xy < XY; xy++){ if(mask == NULL || mask[xy]){ band[xy] /= len[xy]; //divide the band by the vector length } } out.write((char*)band, XYb); if(PROGRESS) progress = (double) (b+1) / (double)B * 50 + 50; } } bool select(std::string outfile, std::vector bandlist, unsigned char* mask = NULL, bool PROGRESS = false) { std::ofstream out(outfile.c_str(), std::ios::binary); //open the output file if (!out) { std::cout << "ERROR opening output file: " << outfile << std::endl; return false; } file.seekg(0, std::ios::beg); //move the pointer to the current file to the beginning size_t B = Z(); //calculate the number of bands size_t XY = X() * Y(); //calculate the number of pixels in a band size_t in_bytes = XY * sizeof(T); //calculate the size of a band in bytes T* in = (T*)malloc(in_bytes); //allocate space for the band image size_t nb = bandlist.size(); //get the number of bands in the output image for (size_t b = 0; b < nb; b++) { band(in, bandlist[b]); //get the band associated with the given wavelength out.write((char*)in, in_bytes); //write the band to the output file if (PROGRESS) progress = (double)(b + 1) / (double)bandlist.size() * 100; } out.close(); free(in); return true; } size_t readlines(T* dest, size_t start, size_t n){ return hsi::read(dest, 0, start, 0, X(), n, Z()); } size_t writeblock(std::ofstream* f, T* src, size_t n){ auto t0 = std::chrono::high_resolution_clock::now(); f->write((char*)src, n); auto t1 = std::chrono::high_resolution_clock::now(); return std::chrono::duration_cast(t1-t0).count(); } /// Convert this BSQ file to a BIL bool bil(std::string outname, bool PROGRESS = false, bool VERBOSE = false, bool OPTIMIZATION = true){ const size_t buffers = 4; //number of buffers required for this algorithm size_t mem_per_batch = binary::buffer_size / buffers; //calculate the maximum memory available for a batch size_t slice_bytes = X() * Z() * sizeof(T); //number of bytes in an input batch slice (Y-slice in this case) size_t max_slices_per_batch = mem_per_batch / slice_bytes; //maximum number of slices we can process in one batch given memory constraints std::cout<<"maximum memory available for processing: "<<(double)binary::buffer_size/(double)1000000<<" MB"<(Y(), max_slices_per_batch); //start with the maximum number of slices that can be stored (may be the entire data set) std::ofstream target(outname.c_str(), std::ios::binary); //open an output file for writing //initialize with buffer 0 (used for double buffering) size_t y_load = 0; size_t y_proc = 0; std::future rthread; std::future wthread; //create asynchronous threads for reading and writing std::chrono::high_resolution_clock::time_point t_start, pt_start; //high-resolution timers std::chrono::high_resolution_clock::time_point t_end, pt_end; size_t t_batch; //number of milliseconds to process a batch size_t t_total = 0; //total time for operation size_t pt_total = 0; //total time spent processing data size_t rt_total = 0; //total time spent reading data size_t wt_total = 0; size_t dr = 0; rt_total += readlines(src[0], 0, N[0]); //read the first batch into the 0 source buffer y_load += N[0]; //increment the loaded slice counter int b = 1; //initialize the double buffer to 0 while(y_proc < Y()){ //while there are still slices to be processed t_start = std::chrono::high_resolution_clock::now(); //start the timer for this batch if(y_load < Y()){ //if there are still slices to be loaded, load them //if(y_proc > 0){ //} if(y_load + N[b] > Y()) N[b] = Y() - y_load; //if the next batch would process more than the total slices, adjust the batch size rthread = std::async(std::launch::async, &stim::bsq::readlines, this, src[b], y_load, N[b]); //rt_total += rthread.get(); y_load += N[b]; //increment the number of loaded slices } b = !b; //swap the double-buffer pt_total += binary::permute(dst[b], src[b], X(), N[b], Z(), 0, 2, 1); //permute the batch to a BIL file wt_total += writeblock(&target, dst[b], N[b] * slice_bytes); //write the permuted data to the output file y_proc += N[b]; //increment the counter of processed pixels if(PROGRESS) progress = (double)( y_proc + 1 ) / Y() * 100; //increment the progress counter based on the number of processed pixels if(y_proc < Y()) rt_total += rthread.get(); //if a new batch was set to load, make sure it loads after calculations t_end = std::chrono::high_resolution_clock::now(); t_batch = std::chrono::duration_cast(t_end-t_start).count(); t_total += t_batch; if(OPTIMIZATION) N[b] = O.update(N[!b] * slice_bytes, t_batch, binary::data_rate, VERBOSE); //set the batch size based on optimization //binary::data_rate = dr; //std::cout<<"New N = "<::buffer_size / buffers; //calculate the maximum memory available for a batch size_t slice_bytes = X() * Z() * sizeof(T); //number of bytes in an input batch slice (Y-slice in this case) size_t max_slices_per_batch = mem_per_batch / slice_bytes; //maximum number of slices we can process in one batch given memory constraints std::cout<<"maximum memory available for processing: "<<(double)binary::buffer_size/(double)1000000<<" MB"<(Y(), max_slices_per_batch); //start with the maximum number of slices that can be stored (may be the entire data set) std::ofstream target(outname.c_str(), std::ios::binary); //open an output file for writing //initialize with buffer 0 (used for double buffering) size_t y_load = 0; size_t y_proc = 0; std::future rthread; std::future wthread; //create asynchronous threads for reading and writing std::chrono::high_resolution_clock::time_point t_start, pt_start; //high-resolution timers std::chrono::high_resolution_clock::time_point t_end, pt_end; size_t t_batch; //number of milliseconds to process a batch size_t t_total = 0; //total time for operation size_t pt_total = 0; //total time spent processing data size_t rt_total = 0; //total time spent reading data size_t wt_total = 0; size_t dr = 0; rt_total += readlines(src[0], 0, N[0]); //read the first batch into the 0 source buffer y_load += N[0]; //increment the loaded slice counter int b = 1; //initialize the double buffer to 0 while(y_proc < Y()){ //while there are still slices to be processed t_start = std::chrono::high_resolution_clock::now(); //start the timer for this batch if(y_load < Y()){ //if there are still slices to be loaded, load them //if(y_proc > 0){ //} if(y_load + N[b] > Y()) N[b] = Y() - y_load; //if the next batch would process more than the total slices, adjust the batch size rthread = std::async(std::launch::async, &stim::bsq::readlines, this, src[b], y_load, N[b]); //rt_total += rthread.get(); y_load += N[b]; //increment the number of loaded slices } b = !b; //swap the double-buffer pt_total += binary::permute(dst[b], src[b], X(), N[b], Z(), 2, 0, 1); //permute the batch to a BIL file wt_total += writeblock(&target, dst[b], N[b] * slice_bytes); //write the permuted data to the output file y_proc += N[b]; //increment the counter of processed pixels if(PROGRESS) progress = (double)( y_proc + 1 ) / Y() * 100; //increment the progress counter based on the number of processed pixels if(y_proc < Y()) rt_total += rthread.get(); //if a new batch was set to load, make sure it loads after calculations t_end = std::chrono::high_resolution_clock::now(); t_batch = std::chrono::duration_cast(t_end-t_start).count(); t_total += t_batch; if(OPTIMIZATION) N[b] = O.update(N[!b] * slice_bytes, t_batch, binary::data_rate, VERBOSE); //set the batch size based on optimization //binary::data_rate = dr; //std::cout<<"New N = "< w[n-1] || rb >w[n-1]){ if (lb < w[0]) { std::cout << "bsq::area ERROR - left bound " << lb << " is below the minimum available wavelength " << w[0] << std::endl; } if (rb < w[0]) { std::cout << "bsq::area ERROR - right bound " << rb << " is below the minimum available wavelength " << w[0] << std::endl; } if (lb > w[n - 1]) { std::cout << "bsq::area ERROR - left bound " << lb << " is above the maximum available wavelength " << w[n - 1] << std::endl; } if (rb > w[n - 1]) { std::cout << "bsq::area ERROR - right bound " << rb << " is above the maximum available wavelength " << w[0] << std::endl; } return false; } //to make sure right bound is bigger than left bound else if(lb > rb){ std::cout << "bsq::area ERROR - right bound " << rb << " should be larger than left bound " << lb << std::endl; return false; } //find the indices of the left and right baseline points while (lab >= w[ai]){ ai++; } while (rab <= w[bi]){ bi--; } band(lp, lb); //get the band images for the left and right baseline points band(rp, rb); // calculate the average value of the indexed region memset(result, 0, S); //initialize the integral to zero (0) //integrate the region between the specified bands and the closest indexed band // this integrates the "tails" of the spectrum that lie outside the main indexed region baseline_band(lb, rb, lp, rp, rab, cur2); //calculate the image for the right-most band in the integral baseline_band(lb, rb, lp, rp, w[bi], cur); //calculate the image for the right-most indexed band 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); } //integrate the main indexed region ai++; for(unsigned long long i = ai; i <= bi ;i++) //for each band in the integral { baseline_band(lb, rb, lp, rp, w[ai], cur2); //get the baselined band for(unsigned long long j = 0; j < XY; j++){ //for each pixel in that band 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){ size_t XYbytes = X() * Y() * sizeof(T); //calculate the size of the band image (in bytes) T* p1 = (T*)malloc(XYbytes); //allocate space for both bands in the ratio T* p2 = (T*)malloc(XYbytes); memset(result, 0, XYbytes); //initialize the ratio to zero //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(mask == NULL || mask[i]){ 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){ size_t bytes = X() * Y() * sizeof(T); T* p1 = (T*)malloc(bytes); //allocate space for both ratio components T* p2 = (T*)malloc(bytes); memset(result, 0, bytes); //initialize the ratio to zero //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(mask == NULL || mask[i]) 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){ size_t bytes = X() * Y() * sizeof(T); T* p1 = (T*)malloc(bytes); //allocate space for each of the operands T* p2 = (T*)malloc(bytes); memset(result, 0, bytes); //initialize the ratio result to zero (0) //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(mask == NULL || mask[i]) //if the pixel is masked result[i] = p1[i] / p2[i]; //calculate the ratio } free(p1); free(p2); return true; } /// Compute the definite integral of a baseline corrected peak weighted by the corresponding wavelength /// @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); //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); memset(result, (char)0, S); //initialize the integral to zero (0) //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++){ T v = (T)((w[ai] - w[ai-1]) * (w[ai] + w[ai-1]) * ((double)cur[j] + (double)cur2[j]) / 4.0); result[j] += v; } 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. /// Note that the values for the centroid can be outside of [lab, rab] if the spectrum goes negative. /// @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){ size_t bytes = X() * Y() * sizeof(T); //calculate the number of bytes in a band image T* p1 = (T*)malloc(X() * Y() * sizeof(T)); //allocate space for both operands T* p2 = (T*)malloc(X() * Y() * sizeof(T)); memset(result, 0, bytes); //initialize the ratio result to zero (0) //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* out_mask, double mask_band, double lower, double upper, unsigned char* mask = NULL, bool PROGRESS = false){ memset(out_mask, 0, X() * Y()); //initialize the mask to zero T* temp = (T*)malloc(X() * Y() * sizeof(T)); //allocate memory for the certain band band(temp, mask_band); for (unsigned long long i = 0; i < X() * Y(); i++) { if(mask == NULL || mask[i] != 0){ if(temp[i] > lower && temp[i] < upper){ out_mask[i] = 255; } else out_mask[i] = 0; } if(PROGRESS) progress = (double) (i+1) / (X() * Y()) * 100; } 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); unsigned long long XY = X() * Y(); //calculate number of a band unsigned long long L = XY * sizeof(T); T * temp = (T*)malloc(L); for (unsigned long long i = 0; i < Z(); i++) //for each spectral bin { band_index(temp, i); //get the specified band (by index) for ( unsigned long long j = 0; j < XY; j++) // for each pixel { if(p[j] == 0){ //if the mask is 0 at that pixel temp[j] = 0; //set temp to zero } else{ continue; } } target.write(reinterpret_cast(temp), L); //write the XY slice at that band to disk if(PROGRESS) progress = (double)(i + 1) / (double)Z() * 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){ unsigned long long XY = X() * Y(); //Number of XY pixels unsigned long long L = XY * sizeof(T); //size of XY plane (in bytes) //calculate the number of pixels in the mask //unsigned long long P = nnz(mask); T* band_image = (T*) malloc( XY * sizeof(T)); //allocate space for a single band unsigned long long i; //pixel index into the sifted array for(unsigned long long b = 0; b < Z(); b++){ //for each band in the data set band_index(band_image, b); //retrieve an image of that band i = 0; for(unsigned long long xy = 0; xy < XY; xy++){ if(mask == NULL || mask[xy] != 0){ //if the pixel is valid matrix[i*Z() + b] = band_image[xy]; //copy it to the appropriate point in the values[] array i++; } if(PROGRESS) progress = (double)(b * XY + xy+1) / (double)(XY * Z()) * 100; } } return true; } /// Saves to disk only those spectra corresponding to mask values != 0 /// @param outfile is the name of the sifted ENVI file to be written to disk /// @param unsigned char* p is the mask file used for sifting bool sift(std::string outfile, unsigned char* p, bool PROGRESS = false){ std::ofstream target(outfile.c_str(), std::ios::binary); // open a band (XY plane) unsigned long long XY = X() * Y(); //Number of XY pixels unsigned long long L = XY * sizeof(T); //size of XY pixels unsigned long long B = Z(); T * temp = (T*)malloc(L); //allocate memory for a band T * temp_vox = (T*)malloc(sizeof(T)); //allocate memory for one voxel for (unsigned long long i = 0; i < B; i++) //for each spectral bin { band_index(temp, i); //get the specified band (XY sheet by index) for (unsigned long long j = 0; j < XY; j++) // for each pixel { if (p[j] != 0){ //if the mask is != 0 at that pixel temp_vox[0] = temp[j]; target.write(reinterpret_cast(temp_vox), sizeof(T)); //write the XY slice at that band to disk } else{ continue; } } if(PROGRESS) progress = (double)(i+1)/ B * 100; } target.close(); free(temp); progress = 100; return true; } /// Generates a spectral image from a matrix of spectral values in lexicographic order and a mask bool unsift(std::string outfile, unsigned char* p, unsigned long long samples, unsigned long long lines, bool PROGRESS = false){ //create a binary output stream std::ofstream target(outfile.c_str(), std::ios::binary); //make sure that there's only one line if(Y() != 1){ std::cout<<"ERROR in stim::bsq::sift() - number of lines does not equal 1"<(unsifted), sizeof(T) * XY); } //std::cout<<"unsifted"< band_values(count); //create an STD vector of band values //this loops goes through each band in B (Z()) // masked (or valid) pixels from that band are averaged and the average is stored in p size_t k; for (size_t i = 0; i < Z(); i++){ //for each band band_index(temp, i); //get the band image and store it in temp k = 0; //initialize the band_value index to zero for (size_t j = 0; j < XY; j++){ //loop through temp, averaging valid pixels if (mask == NULL || mask[j] != 0){ band_values[k] = temp[j]; //store the value in the band_values array k++; //increment the band_values index } } std::sort(band_values.begin(), band_values.end()); //sort all of the values in the band m[i] = band_values[ count/2 ]; //store the center value in the array if(PROGRESS) progress = (double)(i+1) / Z() * 100; //update the progress counter } 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 number of samples, lines, and bands unsigned long long samples = x1 - x0 + 1; unsigned long long lines = y1 - y0 + 1; unsigned long long bands = b1 - b0 + 1; //calculate the size of a single band unsigned long long L = samples * lines * sizeof(T); //allocate space for a single band T* temp = (T*)malloc(L); //create an output stream to store the output file std::ofstream out(outfile.c_str(), std::ios::binary); //calculate the distance required to jump from the end of one band to the beginning of another unsigned long long jumpb = X() * (Y() - lines) * sizeof(T); //calculate the distance required to jump from the end of one line to the beginning of another unsigned long long jumpl = (X() - samples) * sizeof(T); //seek to the start of the cropped region in the input file file.seekg( (b0 * X() * Y() + y0 * X() + x0) * sizeof(T), std::ios::beg); //for each band for (unsigned long long z = b0; z <= b1; z++) { //std::cout<(temp), L); //write slice data into target file file.seekg(jumpb, std::ios::cur); } free(temp); return true; } ///Crop out several subimages and assemble a new image from these concatenated subimages /// @param outfile is the file name for the output image /// @param sx is the width of each subimage /// @param sy is the height of each subimage /// @mask is the mask used to define subimage positions extracted from the input file void subimages(std::string outfile, size_t sx, size_t sy, unsigned char* mask, bool PROGRESS = false){ size_t N = nnz(mask); //get the number of subimages T* dst = (T*) malloc(N * sx * sy * sizeof(T)); //allocate space for a single band of the output image memset(dst, 0, N*sx*sy*sizeof(T)); //initialize the band image to zero std::ofstream out(outfile, std::ios::binary); //open a file for writing T* src = (T*) malloc(X() * Y() * sizeof(T)); for(size_t b = 0; b < Z(); b++){ //for each band band_index(src, b); //load the band image size_t i = 0; //create an image index and initialize it to zero size_t n = 0; while(n < N){ //for each subimage if(mask[i]){ //if the pixel is masked, copy the surrounding pixels into the destination band size_t yi = i / X(); //determine the y position of the current pixel size_t xi = i - yi * X(); //determine the x position of the current pixel if( xi > sx/2 && xi < X() - sx/2 && //if the subimage is completely within the bounds of the original image yi > sy/2 && yi < Y() - sy/2){ size_t cx = xi - sx/2; //calculate the corner position for the subimage size_t cy = yi - sy/2; for(size_t syi = 0; syi < sy; syi++){ //for each line in the subimage size_t src_i = (cy + syi) * X() + cx; //size_t dst_i = syi * (N * sx) + n * sx; size_t dst_i = (n * sy + syi) * sx; memcpy(&dst[dst_i], &src[src_i], sx * sizeof(T)); //copy one line from the subimage to the destination image } n++; } } i++; if(PROGRESS) progress = (double)( (n+1) * (b+1) ) / (N * Z()) * 100; }//end while n out.write((const char*)dst, N * sx * sy * sizeof(T)); //write the band to memory } free(dst); //free memory free(src); } /// 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 XY = X() * Y(); //calculate the number of elements in a band size_t XYb = XY * sizeof(T); //calculate the number of bytes in a band T* temp = (T*)malloc(XYb); //allocate space to store a band size_t i = 0; //store the first index into the band array for(size_t b = 0; b < Z(); b++){ //for each band if(b != band_array[i]){ //if this band is not trimmed file.read((char*)temp, XYb); //read the band out.write((char*)temp, XYb); //output the band } else{ file.seekg(XYb, std::ios::cur); //otherwise, skip the band i++; } if(PROGRESS) progress = (double)(b+1) / (double) Z() * 100; } free(temp); //free the scratch space for the band } /// 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, bsq* 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(X() * Y() * sizeof(T)); //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(C->X() * C->Y() * sizeof(T)); //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(band_dst_bytes); //allocate space for a band of the destination image memset(dst, 0, band_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 file.read((char*)cur, band_bytes); //read a band from the current image C->file.read((char*)src, band_src_bytes); //read a band from the source image for(size_t y = 0; y < Y(); y++) memcpy( &dst[ (p0[1]+y) * S[0] + p0[0] ], &cur[ y * X() ], line_bytes); //copy the line from the current to the destination image //memset( &dst[ (p0[1]+y) * S[0] + p0[0] ], 0, line_dst_bytes); for(size_t y = 0; y < C->Y(); y++) memcpy( &dst[ (p1[1]+y) * S[0] + p1[0] ], &src[ y * C->X() ], line_src_bytes); //copy the line from the source to the destination image out.write((char*)dst, band_dst_bytes); //write the combined image to an output file if(PROGRESS) progress = (double)(b + 1) / (double) Z() * 100; } out.close(); } /// Append an image to this one along the band dimension void append(std::string outfile, bsq* C, 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); if (PROGRESS) progress = 0; out << file.rdbuf(); //copy the data from this ENVI file if (PROGRESS) progress = (double)(Z() + 1) / (double)(Z() + C->Z()) * 100; out << C->file.rdbuf(); //copy the data from the appending file if (PROGRESS) progress = 100; out.close(); } /// 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::ofstream& out, std::vector C, size_t start, size_t end, unsigned char* mask = NULL, bool PROGRESS = false){ size_t nbands = end - start + 1; size_t XY = X() * Y(); //calculate the number of values in a band size_t XYb = XY * sizeof(T); //calculate the size of a band (frame) in bytes file.seekg(XYb * start, std::ios::beg); //move to the beginning of the 'start' band size_t nframes = C.size(); //get the number of bands that the kernel spans std::deque frame(nframes, NULL); //create an array to store pointers to each frame for(size_t f = 0; f < nframes; f++){ //for each frame frame[f] = (T*)malloc(XYb); //allocate space for the frame file.read((char*)frame[f], XYb); //load the frame } T* outband = (T*)malloc(XYb); //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 b = 0; b < nbands; b++){ //for each band memset(outband, 0, XYb); //set the output band to zero (0) size_t c, xy; double coeff; for(c = 0; c < nframes; c++){ //for each frame (corresponding to each coefficient) coeff = C[c]; for(xy = 0; xy < XY; xy++){ //for each pixel if(mask == NULL || mask[xy]){ outband[xy] += (T)(coeff * frame[c][xy]); //calculate the contribution of the current frame (scaled by the corresponding coefficient) } } } out.write((char*)outband, XYb); //output the band file.read((char*)frame[0], XYb); //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)(b+1) / (double)nbands * 100; } } /// Performs a single convolution and saves it to an output file /// @param outfile is the convolved file to be output /// @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 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 convolve(out, C, start, end, mask, PROGRESS); //start the convolution out.close(); } /// Performs a set of convolutions and chains the results together in a single file /// @param outfile is the convolved file to be output /// @param C is an array containing an array of coefficients for each kernel /// @param start is the list of start bands for each kernel /// @param end is the list of end bands for each kernel void convolve(std::string outfile, std::vector< std::vector > C, std::vector start, std::vector 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 K = C.size(); //get the number of kernels for(size_t k = 0; k < K; k++){ size_t b0 = start[k]; //calculate the range of the convolution size_t b1 = end[k]; convolve(out, C[k], b0, b1, mask, PROGRESS); //perform the convolution with the current kernel in the given range } out.close(); } /// 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 XY = X() * Y(); //calculate the number of values in a band size_t XYb = XY * sizeof(T); //calculate the size of a band (frame) in bytes size_t B = Z(); file.seekg(0, std::ios::beg); //move 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(XYb); //allocate space for the frame file.read((char*)frame[f], XYb); //load the frame } T* outband = (T*)malloc(XYb); //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 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], XYb); //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(outband, 0, XYb); //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 xy = 0; xy < XY; xy++){ //for each pixel if(mask == NULL || mask[xy]){ outband[xy] += (T)(C[c] * frame[c][xy]); //calculate the contribution of the current frame (scaled by the corresponding coefficient) } } } out.write((char*)outband, XYb); //output the band if(PROGRESS) progress = (double)(b+1) / (double)B * 100; } } //end deriv bool multiply(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){ unsigned long long B = Z(); //calculate the number of bands 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 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 image c = (T*)malloc( S ); //allocate memory for the band image for(unsigned long long j = 0; j < B; j++){ //for each band band_index(c, j); //load the current band for(unsigned long long i = 0; i < XY; i++){ //for each pixel if(mask == NULL || mask[i]) //if the pixel is masked c[i] *= (T)v; //perform the multiplication } target.write(reinterpret_cast(c), S); //write normalized data into destination if(PROGRESS) progress = (double)(j+1) / B * 100; //update the progress } free(c); //free the band target.close(); //close the output file return true; } bool add(std::string outname, double v, unsigned char* mask = NULL, bool PROGRESS = false){ unsigned long long B = Z(); //calculate the number of bands 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 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 image c = (T*)malloc( S ); //allocate memory for the band image for(unsigned long long j = 0; j < B; j++){ //for each band band_index(c, j); //load the current band for(unsigned long long i = 0; i < XY; i++){ //for each pixel if(mask == NULL || mask[i]) //if the pixel is masked c[i] += (T)v; //perform the multiplication } target.write(reinterpret_cast(c), S); //write normalized data into destination if(PROGRESS) progress = (double)(j+1) / B * 100; //update the progress } free(c); //free the band target.close(); //close the output file return true; } /// Close the file. bool close(){ file.close(); return true; } }; } #endif