Commit source
This commit is contained in:
9
TransmissionDetection/CMakeLists.txt
Normal file
9
TransmissionDetection/CMakeLists.txt
Normal file
@@ -0,0 +1,9 @@
|
||||
project(TranDetection)
|
||||
find_package (OpenMP REQUIRED)
|
||||
file(GLOB_RECURSE cpp_files ./src/*.cpp)
|
||||
add_library(TranDetection SHARED ${cpp_files} )
|
||||
target_include_directories(TranDetection PRIVATE ./src)
|
||||
target_link_libraries(TranDetection PRIVATE OpenMP::OpenMP_CXX fftw3f)
|
||||
target_compile_options(TranDetection PRIVATE ${OpenMP_CXX_FLAGS} -march=native)
|
||||
set_target_properties(TranDetection PROPERTIES PUBLIC_HEADER
|
||||
${CMAKE_CURRENT_LIST_DIR}/src/calculateBankDetectAndHilbertTransformation.hpp)
|
||||
@@ -0,0 +1,191 @@
|
||||
#include <stdio.h>
|
||||
#include <stdlib.h>
|
||||
#include <string.h>
|
||||
#include <complex.h>
|
||||
#include "fftw3.h"
|
||||
#include <omp.h>
|
||||
#include "mathMethods.hpp"
|
||||
|
||||
|
||||
#define REAL 0
|
||||
#define IMAG 1
|
||||
#define FFTW_WISDOM_TYPE FFTW_MEASURE
|
||||
|
||||
/**
|
||||
* Calculation of hilbert transformation and ...
|
||||
*
|
||||
* @param aScan_r, array containing aScans
|
||||
* @param aScansRef_r, array containing corresponding reference aScans
|
||||
* @param numberScans, defining the number of scans in the ascn arrays
|
||||
* @param numberSamples, defining the number of samples of the sacns in the aScan array
|
||||
* @param RESAMPLE_FACTOR, used resample factor
|
||||
* @param nthreads, number of threads used for fft and omp for loop parallelisation
|
||||
* @param[out] resDetect, result of maximum detection
|
||||
* @param[out] resEnvelope, result of envelope of aScan
|
||||
* @param[out] resEnvelopeRef, result of envelope of reference aScan
|
||||
*
|
||||
**/
|
||||
int calculateBankDetectAndHilbertTransformation(float* aScans_r, float* aScansRef_r,int numberScans, int numberSamples, int RESAMPLE_FACTOR, int nthreads, float* resDetect, float* resEnvelope, float* resEnvelopeRef) {
|
||||
|
||||
|
||||
// resampling infos
|
||||
int nresample_c; // for complex hermetian symmetry for upsample 2 -> stays the same (!)
|
||||
int nresample_r; // for real from hermetian symmetry for upsample 2 -> stays the same (!)
|
||||
|
||||
float scale;
|
||||
|
||||
bool even;
|
||||
|
||||
// indices
|
||||
int index;
|
||||
int index_outer;
|
||||
int firstHalf;
|
||||
int endIndex;
|
||||
|
||||
// interim results
|
||||
static fftwf_complex* resultCrossCor_c = NULL;
|
||||
static fftwf_complex* aScans_c_res = NULL;
|
||||
static fftwf_complex* aScansRef_c_res = NULL;
|
||||
static fftwf_complex* aScans_c = NULL;
|
||||
static fftwf_complex* aScansRef_c = NULL;
|
||||
float* resultCrossCor_r;
|
||||
float* aScans;
|
||||
|
||||
// fftw plans
|
||||
static fftwf_plan plan_fftAScans_rc = NULL;
|
||||
static fftwf_plan plan_ifftAScans_cr = NULL;
|
||||
static fftwf_plan plan_ifftAScans_cc = NULL;
|
||||
|
||||
// fftw wisdom
|
||||
char filenameFftwWisdom[200] = "";
|
||||
|
||||
// precalculations
|
||||
nresample_r = numberSamples * RESAMPLE_FACTOR;
|
||||
nresample_c = numberSamples * RESAMPLE_FACTOR / 2;
|
||||
|
||||
scale = float((1.0 / nresample_r)* RESAMPLE_FACTOR);
|
||||
|
||||
even = (nresample_r / 2.0);
|
||||
|
||||
// load wisdom
|
||||
sprintf(filenameFftwWisdom, "fftwf_wisdom_detection_%d.wis", FFTW_WISDOM_TYPE);
|
||||
int loadedWisdomUsed = 0;
|
||||
if(fftwf_import_wisdom_from_filename(filenameFftwWisdom) == 0) {
|
||||
// printf("wisdom not loaded.\n");
|
||||
} else {
|
||||
loadedWisdomUsed = 1;
|
||||
// printf("wisdom loaded.\n");
|
||||
}
|
||||
|
||||
// mem alloc
|
||||
resultCrossCor_c = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numberScans * nresample_r);
|
||||
aScans_c_res = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numberScans * nresample_r);
|
||||
aScansRef_c_res = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numberScans * nresample_r);
|
||||
resultCrossCor_r = (float*)malloc(numberScans * nresample_r * sizeof(float));
|
||||
aScans_c = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numberScans * numberSamples);
|
||||
aScansRef_c = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numberScans * numberSamples);
|
||||
aScans = (float*)malloc(numberScans * numberSamples * sizeof(float));
|
||||
|
||||
/* fftw initializations */
|
||||
// thread ini
|
||||
if ((nthreads > 0)) {
|
||||
if (fftwf_init_threads() == 0){
|
||||
printf("Data input are not the same size. Exiting.");
|
||||
return 1;
|
||||
}
|
||||
fftwf_plan_with_nthreads(nthreads);
|
||||
}
|
||||
|
||||
|
||||
// plan creations
|
||||
plan_fftAScans_rc = fftwf_plan_many_dft_r2c(1, &numberSamples, numberScans, aScans, NULL, 1, numberSamples, aScans_c_res, &nresample_r, 1, nresample_r, FFTW_WISDOM_TYPE);
|
||||
plan_ifftAScans_cr = fftwf_plan_many_dft_c2r(1, &nresample_r, numberScans, resultCrossCor_c, NULL, 1, nresample_r, resultCrossCor_r, &nresample_r, 1, nresample_r, FFTW_WISDOM_TYPE);
|
||||
plan_ifftAScans_cc = fftwf_plan_many_dft(1, &numberSamples, numberScans, aScans_c_res, NULL, 1, nresample_r, aScans_c, &numberSamples, 1, numberSamples, 1, FFTW_WISDOM_TYPE);
|
||||
|
||||
// DFT of input signals
|
||||
fftwf_execute_dft_r2c(plan_fftAScans_rc, aScans_r, aScans_c_res);
|
||||
fftwf_execute_dft_r2c(plan_fftAScans_rc, aScansRef_r, aScansRef_c_res);
|
||||
|
||||
/* Calculus of fft(tab1)* conj(fft(tab2)) (first part) */
|
||||
/* and calculations for hilbert transform */
|
||||
/* even: 0, n/2 -> times 1; 1 to n/2 - 1 -> times2, other zero */
|
||||
/* odd: 0 -> times 1; 1 to (n + 1)/2 - 1 -> times2, other zero */
|
||||
|
||||
#pragma omp parallel for num_threads(nthreads) default(shared) private(index, index_outer, firstHalf, endIndex, even)
|
||||
for (int i = 0; i < numberScans; i++) {
|
||||
|
||||
index_outer = i * nresample_r;
|
||||
|
||||
resultCrossCor_c[index_outer][REAL] = (aScans_c_res[index_outer][REAL] * aScansRef_c_res[index_outer][REAL] + aScans_c_res[index_outer][IMAG] * aScansRef_c_res[index_outer][IMAG]);// / Normalization;// pow(Normalization[REAL], powerFactor);
|
||||
resultCrossCor_c[index_outer][IMAG] = (aScans_c_res[index_outer][IMAG] * aScansRef_c_res[index_outer][REAL] - aScans_c_res[index_outer][REAL] * aScansRef_c_res[index_outer][IMAG]);// / Normalization;// pow(Normalization[REAL], powerFactor);
|
||||
|
||||
firstHalf = index_outer + (nresample_r / 2.0);
|
||||
for (index = index_outer + 1; index < firstHalf; index++) // from second element onwards
|
||||
{
|
||||
|
||||
resultCrossCor_c[index][REAL] = (aScans_c_res[index][REAL] * aScansRef_c_res[index][REAL] + aScans_c_res[index][IMAG] * aScansRef_c_res[index][IMAG]);// / Normalization;// pow(Normalization[REAL], powerFactor);
|
||||
resultCrossCor_c[index][IMAG] = (aScans_c_res[index][IMAG] * aScansRef_c_res[index][REAL] - aScans_c_res[index][REAL] * aScansRef_c_res[index][IMAG]);// / Normalization;// pow(Normalization[REAL], powerFactor);
|
||||
|
||||
aScans_c_res[index][REAL] *= 2.0;
|
||||
aScans_c_res[index][IMAG] *= 2.0;
|
||||
aScansRef_c_res[index][REAL] *= 2.0;
|
||||
aScansRef_c_res[index][IMAG] *= 2.0;
|
||||
}
|
||||
|
||||
if (even){
|
||||
firstHalf++;
|
||||
}else{
|
||||
aScans_c_res[firstHalf][REAL] *= 2.0;
|
||||
aScans_c_res[firstHalf][IMAG] *= 2.0;
|
||||
aScansRef_c_res[firstHalf][REAL] *= 2.0;
|
||||
aScansRef_c_res[firstHalf][IMAG] *= 2.0;
|
||||
|
||||
}
|
||||
|
||||
endIndex = index_outer + nresample_r;
|
||||
for(index = firstHalf; index<endIndex; index++) {
|
||||
aScans_c_res[index][REAL] = 0.0;
|
||||
aScans_c_res[index][IMAG] = 0.0;
|
||||
aScansRef_c_res[index][REAL] = 0.0;
|
||||
aScansRef_c_res[index][IMAG] = 0.0;
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
|
||||
// Execute a IDFT plan
|
||||
fftwf_execute(plan_ifftAScans_cr);
|
||||
fftwf_execute_dft(plan_ifftAScans_cc,aScans_c_res,aScans_c);
|
||||
fftwf_execute_dft(plan_ifftAScans_cc,aScansRef_c_res,aScansRef_c);
|
||||
|
||||
// maximum detection
|
||||
maximumDetection(resultCrossCor_r, numberScans, nresample_r, resDetect);
|
||||
|
||||
// abs calculation, adapt scaling
|
||||
for (int i = 0; i < numberSamples * numberScans; i++) {
|
||||
resEnvelope[i] = sqrt(aScans_c[i][REAL] * aScans_c[i][REAL] + aScans_c[i][IMAG] * aScans_c[i][IMAG]) * scale;
|
||||
resEnvelopeRef[i] = sqrt(aScansRef_c[i][REAL] * aScansRef_c[i][REAL] + aScansRef_c[i][IMAG] * aScansRef_c[i][IMAG]) * scale;
|
||||
}
|
||||
|
||||
// Store Wisdom
|
||||
fftwf_export_wisdom_to_filename(filenameFftwWisdom);
|
||||
|
||||
// clean
|
||||
fftwf_destroy_plan(plan_fftAScans_rc);
|
||||
fftwf_destroy_plan(plan_ifftAScans_cr);
|
||||
fftwf_destroy_plan(plan_ifftAScans_cc);
|
||||
|
||||
fftwf_free(resultCrossCor_c);
|
||||
fftwf_free(aScans_c_res);
|
||||
fftwf_free(aScansRef_c_res);
|
||||
fftwf_free(aScans_c);
|
||||
fftwf_free(aScansRef_c);
|
||||
free(resultCrossCor_r);
|
||||
free(aScans);
|
||||
|
||||
// fftwf_cleanup_threads();
|
||||
// fftwf_cleanup();
|
||||
fftwf_forget_wisdom();
|
||||
|
||||
return 0;
|
||||
}
|
||||
@@ -0,0 +1 @@
|
||||
extern int calculateBankDetectAndHilbertTransformation(float* aScans_r, float* aScansRef_r,int numberScans, int numberSamples, int resampleFactor, int nthreads, float* resDetect, float* resEnvelope, float* resEnvelopeRef);
|
||||
28
TransmissionDetection/src/mathMethods.cpp
Normal file
28
TransmissionDetection/src/mathMethods.cpp
Normal file
@@ -0,0 +1,28 @@
|
||||
#include <stdio.h>
|
||||
|
||||
/**
|
||||
* Returns index of maximum value, for each row (scan)
|
||||
*
|
||||
* @param inArray array containing block of data (assumed to be 2D)
|
||||
* @param n number data (scans)
|
||||
* @param m number data points (samples)
|
||||
* @param[out] outVector pointer to output array, calculated idx for each 0:n-1
|
||||
*
|
||||
**/
|
||||
void maximumDetection(float* inArray, int n, int m, float* outVector) {
|
||||
|
||||
float maxVal;
|
||||
for (int j = 0; j < n; j++) {
|
||||
outVector[j] = 0;
|
||||
maxVal = inArray[j*m];
|
||||
for (int i = 1; i < m; i++) {
|
||||
if (inArray[j * m + i] > maxVal) {
|
||||
maxVal = inArray[j * m + i];
|
||||
outVector[j] = i;
|
||||
}
|
||||
}
|
||||
// printf("\n %i \n",outVector[j]);
|
||||
|
||||
}
|
||||
|
||||
}
|
||||
1
TransmissionDetection/src/mathMethods.hpp
Normal file
1
TransmissionDetection/src/mathMethods.hpp
Normal file
@@ -0,0 +1 @@
|
||||
extern void maximumDetection(float* inArray, int n, int m, float* outVector);
|
||||
Reference in New Issue
Block a user