slight cleaning
This commit is contained in:
parent
8a4cebe729
commit
fc209d3345
|
@ -409,8 +409,6 @@ public:
|
|||
gf2m fld; //we fix q=2^fld.m=fld.n, n=q/2
|
||||
uint T; //the QD's t parameter is 2^T.
|
||||
permutation block_perm; //order of blocks
|
||||
//TODO block_count is (easily) derivable from hperm.
|
||||
uint block_count; //blocks >= block_count are shortened-out
|
||||
permutation hperm; //block permutation of H block used to get G
|
||||
std::vector<uint> block_perms; //dyadic permutations of blocks
|
||||
|
||||
|
|
74
lib/fwht.cpp
74
lib/fwht.cpp
|
@ -1,74 +0,0 @@
|
|||
|
||||
/*
|
||||
* This file is part of Codecrypt.
|
||||
*
|
||||
* Codecrypt is free software: you can redistribute it and/or modify it
|
||||
* under the terms of the GNU Lesser General Public License as published by
|
||||
* the Free Software Foundation, either version 3 of the License, or (at
|
||||
* your option) any later version.
|
||||
*
|
||||
* Codecrypt is distributed in the hope that it will be useful, but WITHOUT
|
||||
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
|
||||
* License for more details.
|
||||
*
|
||||
* You should have received a copy of the GNU Lesser General Public License
|
||||
* along with Codecrypt. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include "fwht.h"
|
||||
|
||||
#include <vector>
|
||||
using namespace std;
|
||||
|
||||
/*
|
||||
* we count on that all integers are sufficiently large.
|
||||
* They should be, largest value occuring should be O(k*n) if initial vector is
|
||||
* consisted only from {0,1}^n, and we don't usually have codes of this size.
|
||||
*/
|
||||
|
||||
static void fwht (vector<int> x, vector<int>&r)
|
||||
{
|
||||
int bs, s;
|
||||
s = x.size();
|
||||
r.resize (s);
|
||||
bs = s >> 1;
|
||||
r.swap (x);
|
||||
while (bs) {
|
||||
x.swap (r);
|
||||
for (uint i = 0; i < s; ++i) {
|
||||
if ( (i / bs) & 1)
|
||||
r[i] = x[i - bs] - x[i];
|
||||
else
|
||||
r[i] = x[i] + x[i + bs];
|
||||
}
|
||||
bs >>= 1;
|
||||
}
|
||||
}
|
||||
|
||||
//we expect correct parameter size and preallocated out.
|
||||
void fwht_dyadic_multiply (const bvector& a, const bvector& b, bvector& out)
|
||||
{
|
||||
|
||||
//lift everyting to Z.
|
||||
vector<int> t, A, B;
|
||||
uint i;
|
||||
|
||||
t.resize (a.size() );
|
||||
A.resize (a.size() );
|
||||
B.resize (a.size() );
|
||||
|
||||
for (i = 0; i < a.size(); ++i) t[i] = a[i];
|
||||
fwht (t, A);
|
||||
|
||||
for (i = 0; i < b.size(); ++i) t[i] = b[i];
|
||||
fwht (t, B);
|
||||
|
||||
//multiply diagonals to A
|
||||
for (i = 0; i < A.size(); ++i) A[i] *= B[i];
|
||||
fwht (A, t);
|
||||
|
||||
uint bitpos = a.size(); //no problem as a.size() == 1<<m == 2^m
|
||||
for (i = 0; i < t.size(); ++i) out[i] = (t[i] & bitpos) ? 1 : 0;
|
||||
}
|
||||
|
152
lib/mce_qd.cpp
152
lib/mce_qd.cpp
|
@ -22,7 +22,7 @@ using namespace ccr;
|
|||
using namespace ccr::mce_qd;
|
||||
|
||||
#include "decoding.h"
|
||||
#include "fwht.h"
|
||||
#include "qd_utils.h"
|
||||
|
||||
#include <set>
|
||||
|
||||
|
@ -57,9 +57,9 @@ int mce_qd::generate (pubkey&pub, privkey&priv, prng&rng,
|
|||
|
||||
//convenience
|
||||
gf2m&fld = priv.fld;
|
||||
std::vector<uint>&Hsig = priv.Hsig;
|
||||
std::vector<uint>&essence = priv.essence;
|
||||
std::vector<uint>&support = priv.support;
|
||||
|
||||
std::vector<uint> support, Hsig;
|
||||
polynomial g;
|
||||
|
||||
//prepare for data
|
||||
|
@ -148,9 +148,8 @@ int mce_qd::generate (pubkey&pub, privkey&priv, prng&rng,
|
|||
|
||||
//now the blocks.
|
||||
uint block_size = 1 << T,
|
||||
h_block_count = (fld.n / 2) / block_size;
|
||||
uint& block_count = priv.block_count;
|
||||
block_count = h_block_count - block_discard;
|
||||
h_block_count = (fld.n / 2) / block_size,
|
||||
block_count = h_block_count - block_discard;
|
||||
|
||||
//assemble blocks to bl
|
||||
std::vector<polynomial> bl, blp;
|
||||
|
@ -181,50 +180,13 @@ int mce_qd::generate (pubkey&pub, privkey&priv, prng&rng,
|
|||
uint attempts = 0;
|
||||
for (attempts = 0; attempts < block_count; ++attempts) {
|
||||
|
||||
/*
|
||||
* try computing the redundancy block of G
|
||||
*
|
||||
* First, co-trace H, then compute G in form
|
||||
*
|
||||
* G^T = [I | X]
|
||||
*
|
||||
* Since H*G^T = [L | R] * [I | X] = L + R*X = 0
|
||||
* we have the solution: X = R^1 * L
|
||||
* (thanks to Rafael Misoczki)
|
||||
*
|
||||
* Inversion is done the quasi-dyadic way:
|
||||
*
|
||||
* - because for QD matrix m=delta(h) the product
|
||||
* m*m = sum(h) * I, binary QD matrix m is either
|
||||
* inversion of itself (m*m=I) or isn't invertible
|
||||
* and m*m=0. sum(h), the "count of ones in QD
|
||||
* signature mod 2", easily determines the result.
|
||||
*
|
||||
* - Using blockwise invertions/multiplications,
|
||||
* gaussian elimination needed to invert the right
|
||||
* square of H can be performed in O(m^2*block_count)
|
||||
* matrix operations. Matrix operations are either
|
||||
* addition (O(t) on QDs), multiplication(O(t log t)
|
||||
* on QDs) or inversion (O(t), as shown above).
|
||||
* Whole proces is therefore quite fast.
|
||||
*
|
||||
* Gaussian elimination on the QD signature should
|
||||
* result in something like this: (for m=3, t=4)
|
||||
*
|
||||
* 1010 0101 1001 1000 0000 0000
|
||||
* 0101 1100 1110 0000 1000 0000
|
||||
* 0111 1110 0100 0000 0000 1000
|
||||
*/
|
||||
|
||||
priv.hperm.generate_random (block_count, rng);
|
||||
permutation hpermInv;
|
||||
priv.hperm.compute_inversion (hpermInv);
|
||||
|
||||
std::vector<std::vector<bvector> > hblocks;
|
||||
bvector tmp;
|
||||
bool failed;
|
||||
uint i, j, k, l;
|
||||
|
||||
bvector col;
|
||||
uint i, j;
|
||||
|
||||
//prepare blocks of h
|
||||
hblocks.resize (block_count);
|
||||
|
@ -233,102 +195,16 @@ int mce_qd::generate (pubkey&pub, privkey&priv, prng&rng,
|
|||
|
||||
//fill them from Hsig
|
||||
for (i = 0; i < block_count; ++i) {
|
||||
tmp.from_poly_cotrace (bl[hpermInv[i]], fld);
|
||||
col.from_poly_cotrace (bl[hpermInv[i]], fld);
|
||||
for (j = 0; j < fld.m; ++j)
|
||||
tmp.get_block (j * block_size,
|
||||
col.get_block (j * block_size,
|
||||
block_size,
|
||||
hblocks[i][j]);
|
||||
}
|
||||
|
||||
/* do a modified QD-blockwise gaussian elimination on hblocks */
|
||||
failed = false;
|
||||
tmp.resize (block_size);
|
||||
for (i = 0; i < fld.m; ++i) { //gauss step
|
||||
//first, find a nonsingular matrix in the column
|
||||
for (j = i; j < fld.m; ++j)
|
||||
if (hblocks[block_count - fld.m + i][j]
|
||||
.hamming_weight() % 2) break;
|
||||
if (j >= fld.m) { //none found, break;
|
||||
failed = true;
|
||||
break;
|
||||
}
|
||||
|
||||
//bring it to correct position (swap it to i-th row)
|
||||
if (j > i) for (k = 0; k < block_count; ++k)
|
||||
hblocks[k][i].swap
|
||||
(hblocks[k][j]);
|
||||
|
||||
//now normalize the row
|
||||
for (j = i; j < fld.m; ++j) {
|
||||
l = hblocks [block_count - fld.m + i]
|
||||
[j].hamming_weight();
|
||||
if (l == 0) continue; //zero is just okay :]
|
||||
if (! (l % 2) ) //singular, make it regular by adding the i-th row
|
||||
for (k = 0;
|
||||
k < block_count;
|
||||
++k)
|
||||
hblocks[k][j].add
|
||||
(hblocks[k][i]);
|
||||
|
||||
//now a matrix is regular, we can easily make it I.
|
||||
//first, multiply the row
|
||||
for (k = 0; k < block_count; ++k) {
|
||||
//don't overwrite the matrix we're counting with
|
||||
if (k == block_count - fld.m + i) continue;
|
||||
fwht_dyadic_multiply
|
||||
(hblocks[block_count - fld.m + i][j],
|
||||
hblocks[k][j], tmp);
|
||||
hblocks[k][j] = tmp;
|
||||
}
|
||||
//change the block on the diagonal
|
||||
fwht_dyadic_multiply
|
||||
(hblocks[block_count - fld.m + i][j],
|
||||
hblocks[block_count - fld.m + i][j], tmp);
|
||||
hblocks[block_count - fld.m + i][j] = tmp;
|
||||
|
||||
//and zero the column below diagonal
|
||||
if (j > i) for (k = 0; k < block_count; ++k)
|
||||
hblocks[k][j].add
|
||||
(hblocks[k][i]);
|
||||
}
|
||||
}
|
||||
|
||||
if (failed) continue;
|
||||
|
||||
for (i = 0; i < fld.m; ++i) { //jordan step
|
||||
//normalize diagonal
|
||||
for (k = 0; k < block_count - i; ++k) {
|
||||
//we can safely rewrite the diagonal here (nothing's behind it)
|
||||
fwht_dyadic_multiply
|
||||
(hblocks[block_count - i - 1][fld.m - i - 1],
|
||||
hblocks[k][fld.m - i - 1], tmp);
|
||||
hblocks[k][fld.m - i - 1] = tmp;
|
||||
}
|
||||
|
||||
//now make zeroes above
|
||||
for (j = i + 1; j < fld.m; ++j) {
|
||||
l = hblocks[block_count - i - 1]
|
||||
[fld.m - j - 1].hamming_weight();
|
||||
if (l == 0) continue; //already zero
|
||||
if (! (l % 2) ) { //nonsingular, fix it by adding diagonal
|
||||
for (k = 0; k < block_count; ++k)
|
||||
hblocks[k][fld.m - j - 1].add
|
||||
(hblocks[k][fld.m - i - 1]);
|
||||
}
|
||||
for (k = 0; k < block_count - i; ++k) {
|
||||
//overwrite is also safe here
|
||||
fwht_dyadic_multiply
|
||||
(hblocks[block_count - i - 1]
|
||||
[fld.m - j - 1],
|
||||
hblocks[k][fld.m - j - 1], tmp);
|
||||
hblocks[k][fld.m - j - 1] = tmp;
|
||||
}
|
||||
//I+I=0
|
||||
for (k = 0; k < block_count; ++k)
|
||||
hblocks[k][fld.m - j - 1].add
|
||||
(hblocks[k][fld.m - i - 1]);
|
||||
}
|
||||
}
|
||||
/* do a modified QD-blockwise gaussian elimination on hblocks.
|
||||
* If it fails, retry. */
|
||||
if (!qd_to_right_echelon_form (hblocks) ) continue;
|
||||
|
||||
pub.qd_sigs.resize (block_count - fld.m);
|
||||
for (uint i = 0; i < block_count - fld.m; ++i) {
|
||||
|
@ -354,6 +230,9 @@ int mce_qd::generate (pubkey&pub, privkey&priv, prng&rng,
|
|||
int privkey::prepare()
|
||||
{
|
||||
uint s, i, j;
|
||||
std::vector<uint> Hsig, support;
|
||||
uint omega;
|
||||
|
||||
//compute H signature from essence
|
||||
Hsig.resize (fld.n / 2);
|
||||
Hsig[0] = fld.inv (essence[fld.m - 1]);
|
||||
|
@ -428,6 +307,7 @@ int privkey::prepare()
|
|||
// (so that it can be applied directly)
|
||||
uint block_size = 1 << T;
|
||||
uint pos, blk_perm;
|
||||
uint block_count = hperm.size();
|
||||
std::vector<uint> sbl1, sbl2, permuted_support;
|
||||
|
||||
sbl1.resize (block_size);
|
||||
|
|
196
lib/qd_utils.cpp
Normal file
196
lib/qd_utils.cpp
Normal file
|
@ -0,0 +1,196 @@
|
|||
|
||||
/*
|
||||
* This file is part of Codecrypt.
|
||||
*
|
||||
* Codecrypt is free software: you can redistribute it and/or modify it
|
||||
* under the terms of the GNU Lesser General Public License as published by
|
||||
* the Free Software Foundation, either version 3 of the License, or (at
|
||||
* your option) any later version.
|
||||
*
|
||||
* Codecrypt is distributed in the hope that it will be useful, but WITHOUT
|
||||
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
|
||||
* License for more details.
|
||||
*
|
||||
* You should have received a copy of the GNU Lesser General Public License
|
||||
* along with Codecrypt. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include "qd_utils.h"
|
||||
|
||||
#include <vector>
|
||||
using namespace std;
|
||||
|
||||
/*
|
||||
* we count on that all integers are sufficiently large.
|
||||
* They should be, largest value occuring should be O(k*n) if initial vector is
|
||||
* consisted only from {0,1}^n, and we don't usually have codes of this size.
|
||||
*/
|
||||
|
||||
static void fwht (vector<int> x, vector<int>&r)
|
||||
{
|
||||
int bs, s;
|
||||
s = x.size();
|
||||
r.resize (s);
|
||||
bs = s >> 1;
|
||||
r.swap (x);
|
||||
while (bs) {
|
||||
x.swap (r);
|
||||
for (uint i = 0; i < s; ++i) {
|
||||
if ( (i / bs) & 1)
|
||||
r[i] = x[i - bs] - x[i];
|
||||
else
|
||||
r[i] = x[i] + x[i + bs];
|
||||
}
|
||||
bs >>= 1;
|
||||
}
|
||||
}
|
||||
|
||||
//we expect correct parameter size and preallocated out.
|
||||
void fwht_dyadic_multiply (const bvector& a, const bvector& b, bvector& out)
|
||||
{
|
||||
|
||||
//lift everyting to Z.
|
||||
vector<int> t, A, B;
|
||||
uint i;
|
||||
|
||||
t.resize (a.size() );
|
||||
A.resize (a.size() );
|
||||
B.resize (a.size() );
|
||||
|
||||
for (i = 0; i < a.size(); ++i) t[i] = a[i];
|
||||
fwht (t, A);
|
||||
|
||||
for (i = 0; i < b.size(); ++i) t[i] = b[i];
|
||||
fwht (t, B);
|
||||
|
||||
//multiply diagonals to A
|
||||
for (i = 0; i < A.size(); ++i) A[i] *= B[i];
|
||||
fwht (A, t);
|
||||
|
||||
uint bitpos = a.size(); //no problem as a.size() == 1<<m == 2^m
|
||||
for (i = 0; i < t.size(); ++i) out[i] = (t[i] & bitpos) ? 1 : 0;
|
||||
}
|
||||
|
||||
bool qd_to_right_echelon_form (std::vector<std::vector<bvector> >&mat)
|
||||
{
|
||||
uint w = mat.size();
|
||||
if (!w) return false;
|
||||
uint h = mat[0].size();
|
||||
if (!h) return false;
|
||||
uint bs = mat[0][0].size();
|
||||
|
||||
uint i, j, k, l;
|
||||
|
||||
/*
|
||||
* Inversion is done the quasi-dyadic way:
|
||||
*
|
||||
* - because for QD matrix m=delta(h) the product
|
||||
* m*m = sum(h) * I, binary QD matrix m is either
|
||||
* inversion of itself (m*m=I) or isn't invertible
|
||||
* and m*m=0. sum(h), the "count of ones in QD
|
||||
* signature mod 2", easily determines the result.
|
||||
*
|
||||
* - Using blockwise invertions/multiplications,
|
||||
* gaussian elimination needed to invert the right
|
||||
* square of H can be performed in O(m^2*block_count)
|
||||
* matrix operations. Matrix operations are either
|
||||
* addition (O(t) on QDs), multiplication(O(t log t)
|
||||
* on QDs) or inversion (O(t), as shown above).
|
||||
* Whole proces is therefore quite fast.
|
||||
*
|
||||
* Gaussian elimination on the QD signature should
|
||||
* result in something like this: (for m=3, t=4)
|
||||
*
|
||||
* 1010 0101 1001 1000 0000 0000
|
||||
* 0101 1100 1110 0000 1000 0000
|
||||
* 0111 1110 0100 0000 0000 1000
|
||||
*/
|
||||
|
||||
bvector tmp;
|
||||
tmp.resize (bs);
|
||||
for (i = 0; i < h; ++i) { //gauss step
|
||||
//first, find a nonsingular matrix in the column
|
||||
for (j = i; j < h; ++j)
|
||||
if (mat[w - h + i][j]
|
||||
.hamming_weight() % 2) break;
|
||||
if (j >= h) //none found, die!
|
||||
return false;
|
||||
|
||||
//bring it to correct position (swap it to i-th row)
|
||||
if (j > i) for (k = 0; k < w; ++k)
|
||||
mat[k][i].swap
|
||||
(mat[k][j]);
|
||||
|
||||
//now normalize the row
|
||||
for (j = i; j < h; ++j) {
|
||||
l = mat [w - h + i]
|
||||
[j].hamming_weight();
|
||||
if (l == 0) continue; //zero is just okay :]
|
||||
if (! (l % 2) ) //singular, make it regular by adding the i-th row
|
||||
for (k = 0;
|
||||
k < w;
|
||||
++k)
|
||||
mat[k][j].add
|
||||
(mat[k][i]);
|
||||
|
||||
//now a matrix is regular, we can easily make it I.
|
||||
//first, multiply the row
|
||||
for (k = 0; k < w; ++k) {
|
||||
//don't overwrite the matrix we're counting with
|
||||
if (k == w - h + i) continue;
|
||||
fwht_dyadic_multiply
|
||||
(mat[w - h + i][j],
|
||||
mat[k][j], tmp);
|
||||
mat[k][j] = tmp;
|
||||
}
|
||||
//change the block on the diagonal
|
||||
fwht_dyadic_multiply
|
||||
(mat[w - h + i][j],
|
||||
mat[w - h + i][j], tmp);
|
||||
mat[w - h + i][j] = tmp;
|
||||
|
||||
//and zero the column below diagonal
|
||||
if (j > i) for (k = 0; k < w; ++k)
|
||||
mat[k][j].add
|
||||
(mat[k][i]);
|
||||
}
|
||||
}
|
||||
|
||||
for (i = 0; i < h; ++i) { //jordan step
|
||||
//normalize diagonal
|
||||
for (k = 0; k < w - i; ++k) {
|
||||
//we can safely rewrite the diagonal here (nothing's behind it)
|
||||
fwht_dyadic_multiply
|
||||
(mat[w - i - 1][h - i - 1],
|
||||
mat[k][h - i - 1], tmp);
|
||||
mat[k][h - i - 1] = tmp;
|
||||
}
|
||||
|
||||
//now make zeroes above
|
||||
for (j = i + 1; j < h; ++j) {
|
||||
l = mat[w - i - 1]
|
||||
[h - j - 1].hamming_weight();
|
||||
if (l == 0) continue; //already zero
|
||||
if (! (l % 2) ) { //nonsingular, fix it by adding diagonal
|
||||
for (k = 0; k < w; ++k)
|
||||
mat[k][h - j - 1].add
|
||||
(mat[k][h - i - 1]);
|
||||
}
|
||||
for (k = 0; k < w - i; ++k) {
|
||||
//overwrite is also safe here
|
||||
fwht_dyadic_multiply
|
||||
(mat[w - i - 1]
|
||||
[h - j - 1],
|
||||
mat[k][h - j - 1], tmp);
|
||||
mat[k][h - j - 1] = tmp;
|
||||
}
|
||||
//I+I=0
|
||||
for (k = 0; k < w; ++k)
|
||||
mat[k][h - j - 1].add
|
||||
(mat[k][h - i - 1]);
|
||||
}
|
||||
}
|
||||
|
||||
return true;
|
||||
}
|
|
@ -16,15 +16,18 @@
|
|||
* along with Codecrypt. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#ifndef _fwht_h_
|
||||
#define _fwht_h_
|
||||
#ifndef _qdutils_h_
|
||||
#define _qdutils_h_
|
||||
|
||||
#include "codecrypt.h"
|
||||
|
||||
using namespace ccr;
|
||||
|
||||
//parameters MUST be of 2^m size.
|
||||
//FWHT matrix mult in O(n log n). parameters MUST be of 2^m size.
|
||||
void fwht_dyadic_multiply (const bvector&, const bvector&, bvector&);
|
||||
|
||||
//create a generator using fwht
|
||||
bool qd_to_right_echelon_form (std::vector<std::vector<bvector> >&matrix);
|
||||
|
||||
#endif
|
||||
|
|
@ -44,7 +44,7 @@ int main()
|
|||
|
||||
ccr::mce_qd::privkey priv;
|
||||
ccr::mce_qd::pubkey pub;
|
||||
ccr::mce_qd::generate (pub, priv, r, 11, 5, 2);
|
||||
ccr::mce_qd::generate (pub, priv, r, 14, 8, 2);
|
||||
|
||||
cout << "cipher size: " << priv.cipher_size() << ' ' << pub.cipher_size() << endl;
|
||||
cout << "plain size: " << priv.plain_size() << ' ' << pub.plain_size() << endl;
|
||||
|
|
Loading…
Reference in a new issue