33 #define INFINITY std::numeric_limits<double>::infinity()
38 void SISO::gen_rsctrellis(
void)
42 int mem_len = gen.cols()-1;
47 rsctrellis.numStates = (1<<mem_len);
48 rsctrellis.prevStates =
new int[2*rsctrellis.numStates];
49 rsctrellis.nextStates =
new int[2*rsctrellis.numStates];
50 rsctrellis.PARout =
new double[2*rsctrellis.numStates];
51 rsctrellis.fm =
new itpp::bin[rsctrellis.numStates];
53 itpp::bvec cases(mem_len);
56 #pragma omp parallel for private(k, j, cases, feedback, out, buffer)
57 for (k=0; k<rsctrellis.numStates; k++)
62 for (j=1; j<(mem_len+1); j++)
64 feedback ^= (gen(0,j)*cases[j-1]);
67 out = feedback*gen(1,0);
68 for (j=1; j<(mem_len+1); j++)
70 out ^= (gen(1,j)*cases[j-1]);
72 rsctrellis.PARout[k+n*rsctrellis.numStates] = (out?1.0:0.0);
75 for (j=mem_len-1; j>0; j--)
77 cases[j] = cases[j-1];
82 rsctrellis.nextStates[k+n*rsctrellis.numStates] = buffer;
83 rsctrellis.prevStates[buffer+n*rsctrellis.numStates] = k;
88 void SISO::rsc_logMAP(itpp::vec &extrinsic_coded, itpp::vec &extrinsic_data,
89 const itpp::vec &intrinsic_coded,
const itpp::vec &apriori_data)
101 int N = apriori_data.length();
113 double* Lc1I =
new double[N];
114 double* Lc2I =
new double[N];
115 #pragma omp parallel for private(n)
118 Lc1I[n] = intrinsic_coded[2*n];
119 Lc2I[n] = intrinsic_coded[2*n+1];
121 double* A0 =
new double[rsctrellis.numStates*N];
122 double* A1 =
new double[rsctrellis.numStates*N];
123 double* A_mid =
new double[N];
124 double* B0 =
new double[rsctrellis.numStates*N];
125 double* B1 =
new double[rsctrellis.numStates*N];
126 buffer = (tail?-INFINITY:0);
127 #pragma omp parallel for private(n,k)
130 for (k=0; k<rsctrellis.numStates; k++)
132 A0[k+n*rsctrellis.numStates] = -INFINITY;
133 A1[k+n*rsctrellis.numStates] = -INFINITY;
134 B0[k+n*rsctrellis.numStates] = buffer;
135 B1[k+n*rsctrellis.numStates] = buffer;
141 A0[0] = Lc2I[0]*rsctrellis.PARout[0];
142 A1[0] = Lc1I[0]+apriori_data[0]+Lc2I[0]*rsctrellis.PARout[rsctrellis.numStates];
147 for (k=0; k<rsctrellis.numStates; k++)
149 buffer =
itpp::log_add(A0[rsctrellis.prevStates[k]+(n-1)*rsctrellis.numStates],
150 A1[rsctrellis.prevStates[k+rsctrellis.numStates]+(n-1)*rsctrellis.numStates]);
151 A0[k+rsctrellis.numStates*n] = Lc2I[n]*rsctrellis.PARout[k]+buffer;
152 A1[k+rsctrellis.numStates*n] = Lc1I[n]+apriori_data[n]+
153 Lc2I[n]*rsctrellis.PARout[k+rsctrellis.numStates]+buffer;
155 A_min =
std::min(A_min, A0[k+rsctrellis.numStates*n]);
157 A_max =
std::max(A_max, A0[k+rsctrellis.numStates*n]);
160 A_mid[n] = (A_min+A_max)/2;
161 if (std::isinf(A_mid[n]))
165 for (k=0; k<rsctrellis.numStates; k++)
167 A0[k+rsctrellis.numStates*n] -= A_mid[n];
168 A1[k+rsctrellis.numStates*n] -= A_mid[n];
173 B0[rsctrellis.prevStates[0]+(N-1)*rsctrellis.numStates] = 0;
174 B1[rsctrellis.prevStates[rsctrellis.numStates]+(N-1)*rsctrellis.numStates] = 0;
175 for (n=N-2; n>=0; n--)
177 for (k=0; k<rsctrellis.numStates; k++)
179 index = rsctrellis.nextStates[k];
180 B0[k+rsctrellis.numStates*n] =
itpp::log_add(B0[index+(n+1)*rsctrellis.numStates]+
181 Lc2I[n+1]*rsctrellis.PARout[index],
182 B1[index+(n+1)*rsctrellis.numStates]+
183 Lc1I[n+1]+apriori_data[n+1]+Lc2I[n+1]*rsctrellis.PARout[index+rsctrellis.numStates]);
184 index = rsctrellis.nextStates[k+rsctrellis.numStates];
185 B1[k+rsctrellis.numStates*n] =
itpp::log_add(B0[index+(n+1)*rsctrellis.numStates]+
186 Lc2I[n+1]*rsctrellis.PARout[index],
187 B1[index+(n+1)*rsctrellis.numStates]+
188 Lc1I[n+1]+apriori_data[n+1]+Lc2I[n+1]*rsctrellis.PARout[index+rsctrellis.numStates]);
191 if (std::isinf(A_mid[n+1]))
195 for (k=0; k<rsctrellis.numStates; k++)
197 B0[k+rsctrellis.numStates*n] -= A_mid[n+1];
198 B1[k+rsctrellis.numStates*n] -= A_mid[n+1];
203 extrinsic_data.set_size(N);
204 extrinsic_coded.set_size(2*N);
205 #pragma omp parallel for private(n, k, sum0, sum1)
210 for (k=0; k<rsctrellis.numStates; k++)
212 sum1 +=
std::exp(A1[k+n*rsctrellis.numStates]+B1[k+n*rsctrellis.numStates]);
213 sum0 +=
std::exp(A0[k+n*rsctrellis.numStates]+B0[k+n*rsctrellis.numStates]);
215 extrinsic_data[n] =
std::log(sum1/sum0)-apriori_data[n];
216 extrinsic_coded[2*n] =
std::log(sum1/sum0)-Lc1I[n];
220 #pragma omp parallel for private(n, k, sum0, sum1)
225 for (k=0; k<rsctrellis.numStates; k++)
227 if (rsctrellis.fm[k])
229 sum1 +=
std::exp(A1[k+n*rsctrellis.numStates]+B1[k+n*rsctrellis.numStates]);
230 sum0 +=
std::exp(A0[k+n*rsctrellis.numStates]+B0[k+n*rsctrellis.numStates]);
234 sum0 +=
std::exp(A1[k+n*rsctrellis.numStates]+B1[k+n*rsctrellis.numStates]);
235 sum1 +=
std::exp(A0[k+n*rsctrellis.numStates]+B0[k+n*rsctrellis.numStates]);
238 extrinsic_coded[2*n+1] =
std::log(sum0/sum1)-Lc2I[n];
242 delete[] rsctrellis.prevStates;
243 delete[] rsctrellis.nextStates;
244 delete[] rsctrellis.PARout;
245 delete[] rsctrellis.fm;
256 void SISO::rsc_maxlogMAP(itpp::vec &extrinsic_coded, itpp::vec &extrinsic_data,
257 const itpp::vec &intrinsic_coded,
const itpp::vec &apriori_data)
269 int N = apriori_data.length();
281 double* Lc1I =
new double[N];
282 double* Lc2I =
new double[N];
283 #pragma omp parallel for private(n)
286 Lc1I[n] = intrinsic_coded[2*n];
287 Lc2I[n] = intrinsic_coded[2*n+1];
289 double* A0 =
new double[rsctrellis.numStates*N];
290 double* A1 =
new double[rsctrellis.numStates*N];
291 double* A_mid =
new double[N];
292 double* B0 =
new double[rsctrellis.numStates*N];
293 double* B1 =
new double[rsctrellis.numStates*N];
294 buffer = (tail?-INFINITY:0);
295 #pragma omp parallel for private(n,k)
298 for (k=0; k<rsctrellis.numStates; k++)
300 A0[k+n*rsctrellis.numStates] = -INFINITY;
301 A1[k+n*rsctrellis.numStates] = -INFINITY;
302 B0[k+n*rsctrellis.numStates] = buffer;
303 B1[k+n*rsctrellis.numStates] = buffer;
309 A0[0] = Lc2I[0]*rsctrellis.PARout[0];
310 A1[0] = Lc1I[0]+apriori_data[0]+Lc2I[0]*rsctrellis.PARout[rsctrellis.numStates];
315 for (k=0; k<rsctrellis.numStates; k++)
317 buffer =
std::max(A0[rsctrellis.prevStates[k]+(n-1)*rsctrellis.numStates],
318 A1[rsctrellis.prevStates[k+rsctrellis.numStates]+(n-1)*rsctrellis.numStates]);
319 A0[k+rsctrellis.numStates*n] = Lc2I[n]*rsctrellis.PARout[k]+buffer;
320 A1[k+rsctrellis.numStates*n] = Lc1I[n]+apriori_data[n]+
321 Lc2I[n]*rsctrellis.PARout[k+rsctrellis.numStates]+buffer;
323 A_min =
std::min(A_min, A0[k+rsctrellis.numStates*n]);
325 A_max =
std::max(A_max, A0[k+rsctrellis.numStates*n]);
328 A_mid[n] = (A_min+A_max)/2;
329 if (std::isinf(A_mid[n]))
331 for (k=0; k<rsctrellis.numStates; k++)
333 A0[k+rsctrellis.numStates*n] -= A_mid[n];
334 A1[k+rsctrellis.numStates*n] -= A_mid[n];
338 B0[rsctrellis.prevStates[0]+(N-1)*rsctrellis.numStates] = 0;
339 B1[rsctrellis.prevStates[rsctrellis.numStates]+(N-1)*rsctrellis.numStates] = 0;
340 for (n=N-2; n>=0; n--)
342 for (k=0; k<rsctrellis.numStates; k++)
344 index = rsctrellis.nextStates[k];
345 B0[k+rsctrellis.numStates*n] =
std::max(B0[index+(n+1)*rsctrellis.numStates]+
346 Lc2I[n+1]*rsctrellis.PARout[index],
347 B1[index+(n+1)*rsctrellis.numStates]+
348 Lc1I[n+1]+apriori_data[n+1]+Lc2I[n+1]*rsctrellis.PARout[index+rsctrellis.numStates]);
349 index = rsctrellis.nextStates[k+rsctrellis.numStates];
350 B1[k+rsctrellis.numStates*n] =
std::max(B0[index+(n+1)*rsctrellis.numStates]+
351 Lc2I[n+1]*rsctrellis.PARout[index],
352 B1[index+(n+1)*rsctrellis.numStates]+
353 Lc1I[n+1]+apriori_data[n+1]+Lc2I[n+1]*rsctrellis.PARout[index+rsctrellis.numStates]);
356 if (std::isinf(A_mid[n+1]))
358 for (k=0; k<rsctrellis.numStates; k++)
360 B0[k+rsctrellis.numStates*n] -= A_mid[n+1];
361 B1[k+rsctrellis.numStates*n] -= A_mid[n+1];
366 extrinsic_data.set_size(N);
367 extrinsic_coded.set_size(2*N);
368 #pragma omp parallel for private(n, k, sum0, sum1)
373 for (k=0; k<rsctrellis.numStates; k++)
375 sum1 =
std::max(sum1, A1[k+n*rsctrellis.numStates]+B1[k+n*rsctrellis.numStates]);
376 sum0 =
std::max(sum0, A0[k+n*rsctrellis.numStates]+B0[k+n*rsctrellis.numStates]);
378 extrinsic_data[n] = (sum1-sum0)-apriori_data[n];
379 extrinsic_coded[2*n] = (sum1-sum0)-Lc1I[n];
383 #pragma omp parallel for private(n, k, sum0, sum1)
388 for (k=0; k<rsctrellis.numStates; k++)
390 if (rsctrellis.fm[k])
392 sum1 =
std::max(sum1, A1[k+n*rsctrellis.numStates]+B1[k+n*rsctrellis.numStates]);
393 sum0 =
std::max(sum0, A0[k+n*rsctrellis.numStates]+B0[k+n*rsctrellis.numStates]);
397 sum0 =
std::max(sum0, A1[k+n*rsctrellis.numStates]+B1[k+n*rsctrellis.numStates]);
398 sum1 =
std::max(sum1, A0[k+n*rsctrellis.numStates]+B0[k+n*rsctrellis.numStates]);
401 extrinsic_coded[2*n+1] = (sum0-sum1)-Lc2I[n];
404 delete[] rsctrellis.prevStates;
405 delete[] rsctrellis.nextStates;
406 delete[] rsctrellis.PARout;
407 delete[] rsctrellis.fm;
418 void SISO::rsc_sova(itpp::vec &extrinsic_data,
const itpp::vec &intrinsic_coded,
419 const itpp::vec &apriori_data,
const int &win_len)
436 int nb_outputs = gen.rows();
441 int nb_states = rsctrellis.numStates;
443 itpp::imat next_state(nb_states,2);
446 bin_out(i).set_size(nb_states, nb_outputs);
447 for (j=0; j<nb_states; j++)
449 bin_out(i)(j,0) =
double(i);
450 bin_out(i)(j,1) = rsctrellis.PARout[j+i*nb_states];
451 next_state(j,i) = rsctrellis.nextStates[j+i*nb_states];
454 itpp::vec bin_inp(
"0 1");
456 int len = apriori_data.length();
459 itpp::mat metr(nb_states,win_len+1);
463 itpp::mat surv(nb_states,win_len+1);
465 itpp::mat inpt(nb_states,win_len+1);
467 itpp::mat diff(nb_states,win_len+1);
469 itpp::mat comp(nb_states,win_len+1);
471 itpp::mat inpc(nb_states,win_len+1);
479 int Cur,Nxt,nxt,sur,b,tmp,idx;
480 itpp::vec buf(nb_outputs);
481 double llb,mtr,dif,cmp,inc,srv,inp;
482 itpp::vec bin(nb_outputs);
483 itpp::ivec cyclic_buffer(win_len);
484 extrinsic_data.set_size(len);
485 for (i = 0; i < len; i++)
489 Nxt = (i+1)%(win_len+1);
490 buf = intrinsic_coded(i*nb_outputs,(i+1)*nb_outputs-1);
491 llb = apriori_data(i);
492 metr.set_col(Nxt, -INFINITY*
itpp::ones(nb_states));
495 for (s = 0; s<nb_states; s++)
497 for (j = 0; j<2; j++)
499 nxt = next_state(s,j);
500 bin = bin_out(j).get_row(s);
501 mtr = bin*buf+metr(s,Cur);
502 mtr += bin_inp(j)*llb;
504 if (metr(nxt,Nxt) < mtr)
506 diff(nxt,Nxt) = mtr-metr(nxt,Nxt);
507 comp(nxt,Nxt) = surv(nxt,Nxt);
508 inpc(nxt,Nxt) = inpt(nxt,Nxt);
516 dif = metr(nxt,Nxt)-mtr;
517 if (dif <= diff(nxt,Nxt))
534 for (j=0; j<win_len; j++)
536 cyclic_buffer(j) = (Nxt-j)%(win_len+1);
537 cyclic_buffer(j) = (cyclic_buffer(j)<0)?(cyclic_buffer(j)+win_len+1):cyclic_buffer(j);
540 for (j=0; j<win_len; j++)
542 inp = inpt(sur,cyclic_buffer(j));
543 extrinsic_data(b) = inp;
545 tmp = cyclic_buffer(j);
546 cmp = comp(sur, tmp);
547 inc = inpc(sur, tmp);
548 dif = diff(sur, tmp);
549 srv = surv(sur, tmp);
551 for (s=j+1; s<=win_len; s++)
566 tmp = cyclic_buffer(s);
567 inp = inpt(
int(srv), tmp);
568 inc = inpt(
int(cmp), tmp);
569 srv = surv(
int(srv), tmp);
570 cmp = surv(
int(cmp), tmp);
572 sur = int(surv(sur, cyclic_buffer(j)));
580 itpp::elem_mult((2.0*extrinsic_data-1.0), SOVA_scaling_factor*sft)-apriori_data;
583 delete[] rsctrellis.prevStates;
584 delete[] rsctrellis.nextStates;
585 delete[] rsctrellis.PARout;
586 delete[] rsctrellis.fm;
589 void SISO::rsc_viterbi(itpp::vec &extrinsic_coded, itpp::vec &extrinsic_data,
590 const itpp::vec &intrinsic_coded,
const itpp::vec &apriori_data,
const int &win_len)
603 int nb_outputs = gen.rows();
608 int nb_states = rsctrellis.numStates;
610 itpp::imat next_state(nb_states,2);
613 bin_out(i).set_size(nb_states, nb_outputs);
614 for (j=0; j<nb_states; j++)
616 bin_out(i)(j,0) =
double(i);
617 bin_out(i)(j,1) = rsctrellis.PARout[j+i*nb_states];
618 next_state(j,i) = rsctrellis.nextStates[j+i*nb_states];
621 itpp::vec bin_inp(
"0 1");
623 int len = apriori_data.length();
626 itpp::mat metr(nb_states,win_len+1);
630 itpp::mat surv(nb_states,win_len+1);
632 itpp::mat inpt(nb_states,win_len+1);
634 itpp::mat parity(nb_states,win_len+1);
638 int Cur,Nxt,nxt,sur,b;
639 itpp::vec buf(nb_outputs);
641 itpp::vec bin(nb_outputs);
642 itpp::ivec cyclic_buffer(win_len);
643 extrinsic_coded.set_size(2*len);
644 extrinsic_data.set_size(len);
645 for (i = 0; i < len; i++)
649 Nxt = (i+1)%(win_len+1);
650 buf = intrinsic_coded(i*nb_outputs,(i+1)*nb_outputs-1);
651 llb = apriori_data(i);
652 metr.set_col(Nxt, -INFINITY*
itpp::ones(nb_states));
655 for (s = 0; s<nb_states; s++)
657 for (j = 0; j<2; j++)
659 nxt = next_state(s,j);
660 bin = bin_out(j).get_row(s);
661 mtr = bin*buf+metr(s,Cur);
662 mtr += bin_inp(j)*llb;
664 if (metr(nxt,Nxt) < mtr)
668 inpt(nxt,Nxt) = bin(0);
669 parity(nxt,Nxt) = bin(1);
681 for (j=0; j<win_len; j++)
683 cyclic_buffer(j) = (Nxt-j)%(win_len+1);
684 cyclic_buffer(j) = (cyclic_buffer(j)<0)?(cyclic_buffer(j)+win_len+1):cyclic_buffer(j);
687 for (j=0; j<win_len; j++)
689 extrinsic_data(b) = inpt(sur,cyclic_buffer(j));
690 extrinsic_coded(2*b) = extrinsic_data(b);
691 extrinsic_coded(2*b+1) = parity(sur,cyclic_buffer(j));
692 sur = int(surv(sur, cyclic_buffer(j)));
698 delete[] rsctrellis.prevStates;
699 delete[] rsctrellis.nextStates;
700 delete[] rsctrellis.PARout;
701 delete[] rsctrellis.fm;
704 if (Viterbi_hard_output_flag==
true)
712 for (i=0; i<len; i++)
714 tmp(i) = ((apriori_data(i)+intrinsic_coded(2*i))>=0)^(extrinsic_data(i)==1.0);
716 itpp::ivec *ptr_idx_matching =
new itpp::ivec;
717 itpp::ivec *ptr_idx_nonmatching =
new itpp::ivec;
721 int idx_nonmatching_len = ptr_idx_nonmatching->length();
722 double error_rate = double(idx_nonmatching_len)/double(len);
728 }
else if (error_rate==1.0)
733 LLR =
std::log((1.0-error_rate)/error_rate);
735 for (i=0; i<ptr_idx_matching->length(); i++)
737 extrinsic_data(ptr_idx_matching->get(i)) = Viterbi_scaling_factor[0]*
738 (2.0*extrinsic_data(ptr_idx_matching->get(i))-1.0)*LLR;
740 for (i=0; i<idx_nonmatching_len; i++)
742 extrinsic_data(ptr_idx_nonmatching->get(i)) = Viterbi_scaling_factor[1]*
743 (2.0*extrinsic_data(ptr_idx_nonmatching->get(i))-1.0)*LLR;
746 extrinsic_data -= apriori_data;
749 tmp.set_length(2*len);
750 for (i=0; i<(2*len); i++)
752 tmp(i) = (intrinsic_coded(i)>=0)^(extrinsic_coded(i)==1.0);
754 delete ptr_idx_matching;
755 delete ptr_idx_nonmatching;
756 ptr_idx_matching =
new itpp::ivec;
757 ptr_idx_nonmatching =
new itpp::ivec;
761 idx_nonmatching_len = ptr_idx_nonmatching->length();
762 error_rate = double(idx_nonmatching_len)/double(2*len);
767 }
else if (error_rate==1.0)
772 LLR =
std::log((1.0-error_rate)/error_rate);
774 for (i=0; i<ptr_idx_matching->length(); i++)
776 extrinsic_coded(ptr_idx_matching->get(i)) = Viterbi_scaling_factor[0]*
777 (2.0*extrinsic_coded(ptr_idx_matching->get(i))-1.0)*LLR;
779 for (i=0; i<idx_nonmatching_len; i++)
781 extrinsic_coded(ptr_idx_nonmatching->get(i)) = Viterbi_scaling_factor[1]*
782 (2.0*extrinsic_coded(ptr_idx_nonmatching->get(i))-1.0)*LLR;
784 delete ptr_idx_matching;
785 delete ptr_idx_nonmatching;
787 extrinsic_coded -= intrinsic_coded;