PRACH tested in Matlab with fading and CFO

This commit is contained in:
Ismael Gomez 2016-05-03 10:24:45 +02:00
parent 85c819dff4
commit c3ebc3bdfc
6 changed files with 171 additions and 75 deletions

View File

@ -1,16 +1,22 @@
%% PRACH Detection Conformance Test
clear
detect_factor=5;
d=18;%linspace(4,14,6);
pDetection2 = zeros(2,length(d));
for dd=1:length(d)
detect_factor=d(dd);
numSubframes = 100; % Number of subframes frames to simulate at each SNR
SNRdB = linspace(-12,-6,8); % SNR points to simulate
foffset = 270.0; % Frequency offset in Hertz
numSubframes = 75; % Number of subframes frames to simulate at each SNR
SNRdB = linspace(-14,10,8); % SNR points to simulate
foffset = 300.0; % Frequency offset in Hertz
add_fading=true;
addpath('../../build/srslte/lib/phch/test')
%% UE Configuration
% User Equipment (UE) settings are specified in the structure |ue|.
ue.NULRB = 6; % 6 Resource Blocks
ue.NULRB = 100; % 6 Resource Blocks
ue.DuplexMode = 'FDD'; % Frequency Division Duplexing (FDD)
ue.CyclicPrefixUL = 'Normal'; % Normal cyclic prefix length
ue.NTxAnts = 1; % Number of transmission antennas
@ -18,11 +24,8 @@ ue.NTxAnts = 1; % Number of transmission antennas
%% PRACH Configuration
prach.Format = 0; % PRACH format: TS36.104, Table 8.4.2.1-1
prach.SeqIdx = 22; % Logical sequence index: TS36.141, Table A.6-1
prach.CyclicShiftIdx = 1; % Cyclic shift index: TS36.141, Table A.6-1
prach.HighSpeed = 0; % Normal mode: TS36.104, Table 8.4.2.1-1
prach.FreqOffset = 0; % Default frequency location
prach.PreambleIdx = 32; % Preamble index: TS36.141, Table A.6-1
info = ltePRACHInfo(ue, prach); % PRACH information
@ -31,7 +34,7 @@ info = ltePRACHInfo(ue, prach); % PRACH information
% TS36.104, Table 8.4.2.1-1 [ <#9 1> ].
chcfg.NRxAnts = 1; % Number of receive antenna
chcfg.DelayProfile = 'EPA'; % Delay profile
chcfg.DelayProfile = 'ETU'; % Delay profile
chcfg.DopplerFreq = 70.0; % Doppler frequency
chcfg.MIMOCorrelation = 'Low'; % MIMO correlation
chcfg.Seed = 1; % Channel seed
@ -64,29 +67,39 @@ for nSNR = 1:length(SNRdB)
% Loop for each subframe
for nsf = 1:numSubframes
prach.SeqIdx = randi(838,1,1)-1; % Logical sequence index: TS36.141, Table A.6-1
prach.CyclicShiftIdx = randi(16,1,1)-1; % Cyclic shift index: TS36.141, Table A.6-1
prach.PreambleIdx = randi(64,1,1)-1; % Preamble index: TS36.141, Table A.6-1
info = ltePRACHInfo(ue, prach); % PRACH information
% PRACH transmission
ue.NSubframe = mod(nsf-1, 10);
ue.NFrame = fix((nsf-1)/10);
% Set PRACH timing offset in us as per TS36.141, Figure 8.4.1.4.2-2
prach.TimingOffset = info.BaseOffset + ue.NSubframe/10.0;
prach.TimingOffset = 0;
% prach.TimingOffset = 0;
% Generate transmit wave
txwave = ltePRACH(ue, prach);
[txwave,prachinfo] = ltePRACH(ue, prach);
% Channel modeling
chcfg.InitTime = (nsf-1)/1000;
[rxwave, fadinginfo] = lteFadingChannel(chcfg, ...
if (add_fading)
chcfg.InitTime = (nsf-1)/1000;
[rxwave, fadinginfo] = lteFadingChannel(chcfg, ...
[txwave; zeros(25, 1)]);
else
rxwave = txwave;
end
% Add noise
noise = N*complex(randn(size(rxwave)), randn(size(rxwave)));
rxwave = rxwave + noise;
% Remove the implementation delay of the channel modeling
rxwave = rxwave((fadinginfo.ChannelFilterDelay + 1):1920, :);
if (add_fading)
rxwave = rxwave((fadinginfo.ChannelFilterDelay + 1):end, :);
end
% Apply frequency offset
t = ((0:size(rxwave, 1)-1)/chcfg.SamplingRate).';
rxwave = rxwave .* repmat(exp(1i*2*pi*foffset*t), ...
@ -95,28 +108,26 @@ for nSNR = 1:length(SNRdB)
% PRACH detection for all cell preamble indices
[detected, offsets] = ltePRACHDetect(ue, prach, rxwave, (0:63).');
[detected_srs] = srslte_prach_detect(ue, prach, rxwave, detect_factor);
[detected_srs, offsets_srs] = srslte_prach_detect(ue, prach, rxwave, detect_factor);
% Test for preamble detection
if (length(detected)==1)
% Test for correct preamble detection
if (detected==prach.PreambleIdx)
detectedCount = detectedCount + 1; % Detected preamble
% Calculate timing estimation error. The true offset is
% PRACH offset plus channel delay
% trueOffset = prach.TimingOffset/1e6 + 310e-9;
% measuredOffset = offsets(1)/chcfg.SamplingRate;
% timingerror = abs(measuredOffset-trueOffset);
%
% % Test for acceptable timing error
% if (timingerror<=2.08e-6)
% detectedCount = detectedCount + 1; % Detected preamble
% else
% disp('Timing error');
% end
trueOffset = prach.TimingOffset/1e6 + 310e-9;
measuredOffset = offsets(1)/chcfg.SamplingRate;
timingerror = abs(measuredOffset-trueOffset);
% Test for acceptable timing error
if (timingerror<=2.08e-6)
detectedCount = detectedCount + 1; % Detected preamble
else
disp('Timing error');
end
else
disp('Detected incorrect preamble');
end
@ -124,8 +135,30 @@ for nSNR = 1:length(SNRdB)
disp('Detected multiple or zero preambles');
end
if (length(detected_srs)==1 && detected_srs==prach.PreambleIdx)
detectedCount_srs = detectedCount_srs + 1;
% Test for preamble detection
if (length(detected_srs)==1)
% Test for correct preamble detection
if (detected_srs==prach.PreambleIdx)
% Calculate timing estimation error. The true offset is
% PRACH offset plus channel delay
trueOffset = prach.TimingOffset/1e6 + 310e-9;
measuredOffset = offsets_srs(1)/1e6;
timingerror = abs(measuredOffset-trueOffset);
% Test for acceptable timing error
if (timingerror<=2.08e-6)
detectedCount_srs = detectedCount_srs + 1; % Detected preamble
else
disp('SRS: Timing error');
end
else
disp('SRS: Detected incorrect preamble');
end
else
fprintf('SRS: Detected %d preambles. D=%.1f, Seq=%3d, NCS=%2d, Idx=%2d\n', ...
length(detected_srs),detect_factor, prach.SeqIdx, prach.CyclicShiftIdx, prach.PreambleIdx);
end
@ -137,12 +170,22 @@ for nSNR = 1:length(SNRdB)
end % of SNR loop
pDetection2(1,dd)=pDetection(1,1);
pDetection2(2,dd)=pDetection(2,1);
end
%% Analysis
plot(SNRdB, pDetection)
legend('Matlab','srsLTE')
grid on
xlabel('SNR (dB)')
ylabel('Pdet')
if (length(SNRdB)>1)
plot(SNRdB, pDetection)
legend('Matlab','srsLTE')
grid on
xlabel('SNR (dB)')
ylabel('Pdet')
else
plot(d,pDetection2)
legend('Matlab','srsLTE')
grid on
xlabel('d')
ylabel('Pdet')
fprintf('Pdet=%.4f%%, Pdet_srs=%.4f%%\n',pDetection(1,nSNR),pDetection(2,nSNR))
end

View File

@ -70,7 +70,7 @@ typedef struct SRSLTE_API {
cf_t dft_seqs[64][839]; // DFT-precoded seqs
uint32_t root_seqs_idx[64]; // Indices of root seqs in seqs table
uint32_t N_roots; // Number of root sequences used in this configuration
// Containers
cf_t *ifft_in;
cf_t *ifft_out;
@ -88,7 +88,11 @@ typedef struct SRSLTE_API {
cf_t *signal_fft;
float detect_factor;
uint32_t deadzone;
float peak_values[65];
uint32_t peak_offsets[65];
} srslte_prach_t;
typedef struct SRSLTE_API {
@ -128,9 +132,17 @@ SRSLTE_API int srslte_prach_detect(srslte_prach_t *p,
uint32_t freq_offset,
cf_t *signal,
uint32_t sig_len,
uint32_t *indices,
uint32_t *indices,
uint32_t *ind_len);
SRSLTE_API int srslte_prach_detect_offset(srslte_prach_t *p,
uint32_t freq_offset,
cf_t *signal,
uint32_t sig_len,
uint32_t *indices,
uint32_t *offsets,
uint32_t *ind_len);
SRSLTE_API void srslte_prach_set_detect_factor(srslte_prach_t *p,
float factor);

View File

@ -92,7 +92,9 @@ endif(RF_FOUND)
if(VOLK_FOUND)
target_link_libraries(srslte ${VOLK_LIBRARIES})
target_link_libraries(srslte_static ${VOLK_LIBRARIES})
if(NOT DisableMEX)
target_link_libraries(srslte_static ${VOLK_LIBRARIES})
endif(NOT DisableMEX)
endif(VOLK_FOUND)
INSTALL(TARGETS srslte DESTINATION ${LIBRARY_DIR})

View File

@ -34,7 +34,8 @@
//PRACH detection threshold is PRACH_DETECT_FACTOR*average
#define PRACH_DETECT_FACTOR 10
#define PRACH_DETECT_FACTOR 18
#define CFO_REPLICA_FACTOR 0.3
#define N_SEQS 64 // Number of prach sequences available
#define N_RB_SC 12 // Number of subcarriers per resource block
@ -339,6 +340,7 @@ int srslte_prach_init(srslte_prach_t *p,
p->zczc = zero_corr_zone_config;
p->detect_factor = PRACH_DETECT_FACTOR;
// Determine N_zc and N_cs
if(4 == preamble_format){
p->N_zc = 139;
@ -370,7 +372,7 @@ int srslte_prach_init(srslte_prach_t *p,
return SRSLTE_ERROR;
}
srslte_dft_plan_set_mirror(p->zc_ifft, false);
srslte_dft_plan_set_norm(p->zc_ifft, true);
srslte_dft_plan_set_norm(p->zc_ifft, false);
// Generate our 64 sequences
p->N_roots = 0;
@ -380,7 +382,7 @@ int srslte_prach_init(srslte_prach_t *p,
for(int i=0;i<N_SEQS;i++){
srslte_dft_run(p->zc_fft, p->seqs[i], p->dft_seqs[i]);
}
// Create our FFT objects and buffers
p->N_ifft_ul = N_ifft_ul;
if(4 == preamble_format){
@ -388,6 +390,16 @@ int srslte_prach_init(srslte_prach_t *p,
}else{
p->N_ifft_prach = p->N_ifft_ul * DELTA_F/DELTA_F_RA;
}
/* The deadzone specifies the number of samples at the end of the correlation window
* that will be considered as belonging to the next preamble
*/
p->deadzone = 0;
/*
if(p->N_cs != 0) {
float samp_rate=15000*p->N_ifft_ul;
p->deadzone = (uint32_t) ceil((float) samp_rate/((float) p->N_zc*subcarrier_spacing));
}*/
p->ifft_in = (cf_t*)srslte_vec_malloc(p->N_ifft_prach*sizeof(cf_t));
p->ifft_out = (cf_t*)srslte_vec_malloc(p->N_ifft_prach*sizeof(cf_t));
@ -412,7 +424,7 @@ int srslte_prach_init(srslte_prach_t *p,
}
srslte_dft_plan_set_mirror(p->fft, true);
srslte_dft_plan_set_norm(p->fft, true);
srslte_dft_plan_set_norm(p->fft, false);
p->N_seq = prach_Tseq[p->f]*p->N_ifft_ul/2048;
p->N_cp = prach_Tcp[p->f]*p->N_ifft_ul/2048;
@ -468,11 +480,22 @@ void srslte_prach_set_detect_factor(srslte_prach_t *p, float ratio) {
}
int srslte_prach_detect(srslte_prach_t *p,
uint32_t freq_offset,
cf_t *signal,
uint32_t sig_len,
uint32_t *indices,
uint32_t *n_indices)
uint32_t freq_offset,
cf_t *signal,
uint32_t sig_len,
uint32_t *indices,
uint32_t *n_indices)
{
return srslte_prach_detect_offset(p, freq_offset, signal, sig_len, indices, NULL, n_indices);
}
int srslte_prach_detect_offset(srslte_prach_t *p,
uint32_t freq_offset,
cf_t *signal,
uint32_t sig_len,
uint32_t *indices,
uint32_t *offsets,
uint32_t *n_indices)
{
int ret = SRSLTE_ERROR;
if(p != NULL &&
@ -484,11 +507,10 @@ int srslte_prach_detect(srslte_prach_t *p,
fprintf(stderr, "srslte_prach_detect: Signal is not of length %d", p->N_ifft_prach);
return SRSLTE_ERROR_INVALID_INPUTS;
}
// FFT incoming signal
srslte_dft_run(p->fft, signal, p->signal_fft);
memset(p->prach_bins, 0, sizeof(cf_t)*p->N_zc);
*n_indices = 0;
// Extract bins of interest
@ -497,28 +519,18 @@ int srslte_prach_detect(srslte_prach_t *p,
uint32_t K = DELTA_F/DELTA_F_RA;
uint32_t begin = PHI + (K*k_0) + (K/2);
for(int i=0;i<p->N_zc;i++){
p->prach_bins[i] = p->signal_fft[begin+i];
}
memcpy(p->prach_bins, &p->signal_fft[begin], p->N_zc*sizeof(cf_t));
for(int i=0;i<p->N_roots;i++){
memset(p->corr_spec, 0, sizeof(cf_t)*p->N_zc);
memset(p->corr, 0, sizeof(float)*p->N_zc);
float corr_max = 0;
float corr_ave = 0;
cf_t *root_spec = p->dft_seqs[p->root_seqs_idx[i]];
srslte_vec_prod_conj_ccc(p->prach_bins, root_spec, p->corr_spec, p->N_zc);
srslte_dft_run(p->zc_ifft, p->corr_spec, p->corr_spec);
srslte_vec_abs_cf(p->corr_spec, p->corr, p->N_zc);
float norm = sqrtf(p->N_zc);
srslte_vec_sc_prod_fff(p->corr, 1.0/norm, p->corr, p->N_zc);
corr_ave = srslte_vec_acc_ff(p->corr, p->N_zc)/p->N_zc;
srslte_vec_abs_square_cf(p->corr_spec, p->corr, p->N_zc);
float corr_ave = srslte_vec_acc_ff(p->corr, p->N_zc)/p->N_zc;
uint32_t winsize = 0;
if(p->N_cs != 0){
@ -527,18 +539,37 @@ int srslte_prach_detect(srslte_prach_t *p,
winsize = p->N_zc;
}
uint32_t n_wins = p->N_zc/winsize;
for(int j=0;j<n_wins;j++){
float max_peak = 0;
for(int j=0;j<n_wins;j++) {
uint32_t start = (p->N_zc-(j*p->N_cs))%p->N_zc;
uint32_t end = start+winsize;
corr_max = 0;
if (end>p->deadzone) {
end-=p->deadzone;
}
start += p->deadzone;
p->peak_values[j] = 0;
for(int k=start;k<end;k++){
if(p->corr[k] > corr_max){
corr_max = p->corr[k];
if(p->corr[k] > p->peak_values[j]) {
p->peak_values[j] = p->corr[k];
p->peak_offsets[j] = k-start;
if (p->peak_values[j] > max_peak) {
max_peak = p->peak_values[j];
}
}
}
if(corr_max > p->detect_factor*corr_ave){
indices[*n_indices] = (i*n_wins)+j;
(*n_indices)++;
}
if (max_peak > p->detect_factor*corr_ave) {
for (int j=0;j<n_wins;j++) {
if(p->peak_values[j] > p->detect_factor*corr_ave &&
p->peak_values[j] >= CFO_REPLICA_FACTOR*max_peak)
{
indices[*n_indices] = (i*n_wins)+j;
if (offsets) {
offsets[*n_indices] = p->peak_offsets[j];
}
(*n_indices)++;
}
}
}
}

View File

@ -51,6 +51,8 @@ void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
help();
return;
}
srslte_use_standard_symbol_size(true);
uint32_t n_ul_rb = 0;
if (mexutils_read_uint32_struct(UECFG, "NULRB", &n_ul_rb)) {
@ -92,6 +94,7 @@ void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
int nof_samples = mexutils_read_cf(INPUT, &input_signal);
uint32_t preambles[64];
uint32_t offsets[64];
uint32_t nof_detected = 0;
if (nrhs > NOF_INPUTS) {
@ -99,7 +102,7 @@ void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
srslte_prach_set_detect_factor(&prach, factor);
}
if (srslte_prach_detect(&prach, frequency_offset, &input_signal[prach.N_cp], nof_samples, preambles, &nof_detected)) {
if (srslte_prach_detect_offset(&prach, frequency_offset, &input_signal[prach.N_cp], nof_samples, preambles, offsets, &nof_detected)) {
mexErrMsgTxt("Error detecting PRACH\n");
return;
}
@ -107,6 +110,9 @@ void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
if (nlhs >= 1) {
mexutils_write_int((int*) preambles, &plhs[0], nof_detected, 1);
}
if (nlhs >= 2) {
mexutils_write_int((int*) offsets, &plhs[1], nof_detected, 1);
}
free(input_signal);
srslte_prach_free(&prach);

View File

@ -105,6 +105,8 @@ int main(int argc, char **argv) {
for(int i=0;i<64;i++)
indices[i] = 0;
srslte_prach_set_detect_factor(p, 10);
for(seq_index=0;seq_index<n_seqs;seq_index++)
{
srslte_prach_gen(p,