Add support for NIST256P1 elliptic curve

This enables SSH ECDSA public key authentication.
pull/25/head
Roman Zeyde 9 years ago
parent c58d4e03c5
commit 7c58fc11a4

5
.gitignore vendored

@ -8,3 +8,8 @@ build/
dist/
MANIFEST
TrezorCrypto.c
SConstruct
.sconsign.dblite
*.os
*.so
*.pyc

@ -30,7 +30,7 @@ ifdef SMALL
CFLAGS += -DUSE_PRECOMPUTED_IV=0 -DUSE_PRECOMPUTED_CP=0
endif
OBJS = bignum.o ecdsa.o secp256k1.o rand.o hmac.o bip32.o bip39.o pbkdf2.o base58.o
OBJS = bignum.o ecdsa.o secp256k1.o nist256p1.o rand.o hmac.o bip32.o bip39.o pbkdf2.o base58.o
OBJS += ripemd160.o
OBJS += sha2.o
OBJS += aescrypt.o aeskey.o aestab.o aes_modes.o

@ -140,18 +140,18 @@ void bn_rshift(bignum256 *a)
a->val[8] >>= 1;
}
// multiply x by 3/2 modulo prime.
// multiply x by 1/2 modulo prime.
// assumes x < 2*prime,
// guarantees x < 4*prime on exit.
void bn_mult_3_2(bignum256 * x, const bignum256 *prime)
void bn_mult_half(bignum256 * x, const bignum256 *prime)
{
int j;
uint32_t xodd = -(x->val[0] & 1);
// compute x = 3*x/2 mod prime
// if x is odd compute (3*x+prime)/2
uint32_t tmp1 = (3*x->val[0] + (prime->val[0] & xodd)) >> 1;
// compute x = x/2 mod prime
// if x is odd compute (x+prime)/2
uint32_t tmp1 = (x->val[0] + (prime->val[0] & xodd)) >> 1;
for (j = 0; j < 8; j++) {
uint32_t tmp2 = (3*x->val[j+1] + (prime->val[j+1] & xodd));
uint32_t tmp2 = (x->val[j+1] + (prime->val[j+1] & xodd));
tmp1 += (tmp2 & 1) << 29;
x->val[j] = tmp1 & 0x3fffffff;
tmp1 >>= 30;
@ -160,6 +160,20 @@ void bn_mult_3_2(bignum256 * x, const bignum256 *prime)
x->val[8] = tmp1;
}
// multiply x by k modulo prime.
// assumes x < prime,
// guarantees x < prime on exit.
void bn_mult_k(bignum256 *x, uint8_t k, const bignum256 *prime)
{
int j;
for (j = 0; j < 9; j++) {
x->val[j] = k * x->val[j];
}
bn_normalize(x);
bn_fast_mod(x, prime);
bn_mod(x, prime);
}
// assumes x < 2*prime, result < prime
void bn_mod(bignum256 *x, const bignum256 *prime)
{
@ -186,16 +200,10 @@ void bn_mod(bignum256 *x, const bignum256 *prime)
}
}
// Compute x := k * x (mod prime)
// both inputs must be smaller than 2 * prime.
// result is reduced to 0 <= x < 2 * prime
// This only works for primes between 2^256-2^196 and 2^256.
// this particular implementation accepts inputs up to 2^263 or 128*prime.
void bn_multiply(const bignum256 *k, bignum256 *x, const bignum256 *prime)
void bn_multiply_long(const bignum256 *k, const bignum256 *x, uint32_t res[18])
{
int i, j;
uint64_t temp = 0;
uint32_t res[18], coef;
// compute lower half of long multiplication
for (i = 0; i < 9; i++)
@ -216,43 +224,69 @@ void bn_multiply(const bignum256 *k, bignum256 *x, const bignum256 *prime)
temp >>= 30;
}
res[17] = temp;
}
void bn_multiply_reduce_step(uint32_t res[18], const bignum256 *prime, uint32_t i) {
// let k = i-8.
// invariants:
// res[0..(i+1)] = k * x (mod prime)
// 0 <= res < 2^(30k + 256) * (2^30 + 1)
// estimate (res / prime)
// coef = res / 2^(30k + 256) rounded down
// 0 <= coef <= 2^30
// subtract (coef * 2^(30k) * prime) from res
// note that we unrolled the first iteration
uint32_t j;
uint32_t coef = (res[i] >> 16) + (res[i + 1] << 14);
uint64_t temp = 0x1000000000000000ull + res[i - 8] - prime->val[0] * (uint64_t)coef;
res[i - 8] = temp & 0x3FFFFFFF;
for (j = 1; j < 9; j++) {
temp >>= 30;
temp += 0xFFFFFFFC0000000ull + res[i - 8 + j] - prime->val[j] * (uint64_t)coef;
res[i - 8 + j] = temp & 0x3FFFFFFF;
}
temp >>= 30;
temp += 0xFFFFFFFC0000000ull + res[i - 8 + j];
res[i - 8 + j] = temp & 0x3FFFFFFF;
// we rely on the fact that prime > 2^256 - 2^196
// res = oldres - coef*2^(30k) * prime;
// and
// coef * 2^(30k + 256) <= oldres < (coef+1) * 2^(30k + 256)
// Hence, 0 <= res < 2^30k (2^256 + coef * (2^256 - prime))
// Since coef * (2^256 - prime) < 2^226, we get
// 0 <= res < 2^(30k + 226) (2^30 + 1)
// Thus the invariant holds again.
}
void bn_multiply_reduce(bignum256 *x, uint32_t res[18], const bignum256 *prime)
{
int i;
// res = k * x is a normalized number (every limb < 2^30)
// 0 <= res < 2^526.
// compute modulo p division is only estimated so this may give result greater than prime but not bigger than 2 * prime
for (i = 16; i >= 8; i--) {
// let k = i-8.
// invariants:
// res[0..(i+1)] = k * x (mod prime)
// 0 <= res < 2^(30k + 256) * (2^30 + 1)
// estimate (res / prime)
coef = (res[i] >> 16) + (res[i + 1] << 14);
// coef = res / 2^(30k + 256) rounded down
// 0 <= coef <= 2^30
// subtract (coef * 2^(30k) * prime) from res
// note that we unrolled the first iteration
temp = 0x1000000000000000ull + res[i - 8] - prime->val[0] * (uint64_t)coef;
res[i - 8] = temp & 0x3FFFFFFF;
for (j = 1; j < 9; j++) {
temp >>= 30;
temp += 0xFFFFFFFC0000000ull + res[i - 8 + j] - prime->val[j] * (uint64_t)coef;
res[i - 8 + j] = temp & 0x3FFFFFFF;
}
// we don't clear res[i+1] but we never read it again.
// we rely on the fact that prime > 2^256 - 2^196
// res = oldres - coef*2^(30k) * prime;
// and
// coef * 2^(30k + 256) <= oldres < (coef+1) * 2^(30k + 256)
// Hence, 0 <= res < 2^30k (2^256 + coef * (2^256 - prime))
// Since coef * (2^256 - prime) < 2^226, we get
// 0 <= res < 2^(30k + 226) (2^30 + 1)
// Thus the invariant holds again.
bn_multiply_reduce_step(res, prime, i);
bn_multiply_reduce_step(res, prime, i); // apply twice, as a hack for NIST256P1 prime.
assert(res[i + 1] == 0);
}
// store the result
for (i = 0; i < 9; i++) {
x->val[i] = res[i];
}
}
// Compute x := k * x (mod prime)
// both inputs must be smaller than 2 * prime.
// result is reduced to 0 <= x < 2 * prime
// This only works for primes between 2^256-2^196 and 2^256.
// this particular implementation accepts inputs up to 2^263 or 128*prime.
void bn_multiply(const bignum256 *k, bignum256 *x, const bignum256 *prime)
{
uint32_t res[18] = {0};
bn_multiply_long(k, x, res);
bn_multiply_reduce(x, res, prime);
MEMSET_BZERO(res, sizeof(res));
}
@ -660,9 +694,9 @@ void bn_addmodi(bignum256 *a, uint32_t b, const bignum256 *prime) {
void bn_subtractmod(const bignum256 *a, const bignum256 *b, bignum256 *res, const bignum256 *prime)
{
int i;
uint32_t temp = 0;
uint32_t temp = 1;
for (i = 0; i < 9; i++) {
temp += a->val[i] + 2u * prime->val[i] - b->val[i];
temp += 0x3FFFFFFF + a->val[i] + 2u * prime->val[i] - b->val[i];
res->val[i] = temp & 0x3FFFFFFF;
temp >>= 30;
}

@ -57,7 +57,9 @@ void bn_lshift(bignum256 *a);
void bn_rshift(bignum256 *a);
void bn_mult_3_2(bignum256 *x, const bignum256 *prime);
void bn_mult_half(bignum256 *x, const bignum256 *prime);
void bn_mult_k(bignum256 *x, uint8_t k, const bignum256 *prime);
void bn_mod(bignum256 *x, const bignum256 *prime);

@ -32,6 +32,9 @@
#include "ripemd160.h"
#include "base58.h"
#include "macros.h"
#include "secp256k1.h"
static const ecdsa_curve *default_curve = &secp256k1;
int hdnode_from_xpub(uint32_t depth, uint32_t fingerprint, uint32_t child_num, const uint8_t *chain_code, const uint8_t *public_key, HDNode *out)
{
@ -56,7 +59,7 @@ int hdnode_from_xprv(uint32_t depth, uint32_t fingerprint, uint32_t child_num, c
if (bn_is_zero(&a)) { // == 0
failed = true;
} else {
if (!bn_is_less(&a, &order256k1)) { // >= order
if (!bn_is_less(&a, &default_curve->order)) { // >= order
failed = true;
}
MEMSET_BZERO(&a, sizeof(a));
@ -91,7 +94,7 @@ int hdnode_from_seed(const uint8_t *seed, int seed_len, HDNode *out)
if (bn_is_zero(&a)) { // == 0
failed = true;
} else {
if (!bn_is_less(&a, &order256k1)) { // >= order
if (!bn_is_less(&a, &default_curve->order)) { // >= order
failed = true;
}
MEMSET_BZERO(&a, sizeof(a));
@ -135,11 +138,11 @@ int hdnode_private_ckd(HDNode *inout, uint32_t i)
bool failed = false;
if (!bn_is_less(&b, &order256k1)) { // >= order
if (!bn_is_less(&b, &default_curve->order)) { // >= order
failed = true;
}
if (!failed) {
bn_addmod(&a, &b, &order256k1);
bn_addmod(&a, &b, &default_curve->order);
if (bn_is_zero(&a)) {
failed = true;
}
@ -182,7 +185,7 @@ int hdnode_public_ckd(HDNode *inout, uint32_t i)
memset(inout->private_key, 0, 32);
bool failed = false;
if (!ecdsa_read_pubkey(inout->public_key, &a)) {
if (!ecdsa_read_pubkey(default_curve, inout->public_key, &a)) {
failed = true;
}
@ -190,15 +193,15 @@ int hdnode_public_ckd(HDNode *inout, uint32_t i)
hmac_sha512(inout->chain_code, 32, data, sizeof(data), I);
memcpy(inout->chain_code, I + 32, 32);
bn_read_be(I, &c);
if (!bn_is_less(&c, &order256k1)) { // >= order
if (!bn_is_less(&c, &default_curve->order)) { // >= order
failed = true;
}
}
if (!failed) {
scalar_multiply(&c, &b); // b = c * G
point_add(&a, &b); // b = a + b
if (!ecdsa_validate_pubkey(&b)) {
scalar_multiply(default_curve, &c, &b); // b = c * G
point_add(default_curve, &a, &b); // b = a + b
if (!ecdsa_validate_pubkey(default_curve, &b)) {
failed = true;
}
}
@ -291,7 +294,7 @@ int hdnode_private_ckd_cached(HDNode *inout, const uint32_t *i, size_t i_count)
void hdnode_fill_public_key(HDNode *node)
{
ecdsa_get_public_key33(node->private_key, node->public_key);
ecdsa_get_public_key33(default_curve, node->private_key, node->public_key);
}
void hdnode_serialize(const HDNode *node, uint32_t version, char use_public, char *str, int strsize)

@ -26,6 +26,7 @@
#include <stdint.h>
#include <stdlib.h>
#include "ecdsa.h"
#include "options.h"
typedef struct {

@ -35,6 +35,9 @@
#include "base58.h"
#include "macros.h"
#include "secp256k1.h"
#include "nist256p1.h"
// Set cp2 = cp1
void point_copy(const curve_point *cp1, curve_point *cp2)
{
@ -42,10 +45,8 @@ void point_copy(const curve_point *cp1, curve_point *cp2)
}
// cp2 = cp1 + cp2
void point_add(const curve_point *cp1, curve_point *cp2)
void point_add(const ecdsa_curve *curve, const curve_point *cp1, curve_point *cp2)
{
int i;
uint32_t temp;
bignum256 lambda, inv, xr, yr;
if (point_is_infinity(cp1)) {
@ -56,7 +57,7 @@ void point_add(const curve_point *cp1, curve_point *cp2)
return;
}
if (point_is_equal(cp1, cp2)) {
point_double(cp2);
point_double(curve, cp2);
return;
}
if (point_is_negative_of(cp1, cp2)) {
@ -64,39 +65,34 @@ void point_add(const curve_point *cp1, curve_point *cp2)
return;
}
bn_subtractmod(&(cp2->x), &(cp1->x), &inv, &prime256k1);
bn_inverse(&inv, &prime256k1);
bn_subtractmod(&(cp2->y), &(cp1->y), &lambda, &prime256k1);
bn_multiply(&inv, &lambda, &prime256k1);
bn_subtractmod(&(cp2->x), &(cp1->x), &inv, &curve->prime);
bn_inverse(&inv, &curve->prime);
bn_subtractmod(&(cp2->y), &(cp1->y), &lambda, &curve->prime);
bn_multiply(&inv, &lambda, &curve->prime);
// xr = lambda^2 - x1 - x2
xr = lambda;
bn_multiply(&xr, &xr, &prime256k1);
temp = 1;
for (i = 0; i < 9; i++) {
temp += 0x3FFFFFFF + xr.val[i] + 2u * prime256k1.val[i] - cp1->x.val[i] - cp2->x.val[i];
xr.val[i] = temp & 0x3FFFFFFF;
temp >>= 30;
}
bn_fast_mod(&xr, &prime256k1);
bn_mod(&xr, &prime256k1);
bn_multiply(&xr, &xr, &curve->prime);
yr = cp1->x;
bn_addmod(&yr, &(cp2->x), &curve->prime);
bn_subtractmod(&xr, &yr, &xr, &curve->prime);
bn_fast_mod(&xr, &curve->prime);
bn_mod(&xr, &curve->prime);
// yr = lambda (x1 - xr) - y1
bn_subtractmod(&(cp1->x), &xr, &yr, &prime256k1);
bn_multiply(&lambda, &yr, &prime256k1);
bn_subtractmod(&yr, &(cp1->y), &yr, &prime256k1);
bn_fast_mod(&yr, &prime256k1);
bn_mod(&yr, &prime256k1);
bn_subtractmod(&(cp1->x), &xr, &yr, &curve->prime);
bn_multiply(&lambda, &yr, &curve->prime);
bn_subtractmod(&yr, &(cp1->y), &yr, &curve->prime);
bn_fast_mod(&yr, &curve->prime);
bn_mod(&yr, &curve->prime);
cp2->x = xr;
cp2->y = yr;
}
// cp = cp + cp
void point_double(curve_point *cp)
void point_double(const ecdsa_curve *curve, curve_point *cp)
{
int i;
uint32_t temp;
bignum256 lambda, xr, yr;
if (point_is_infinity(cp)) {
@ -107,31 +103,32 @@ void point_double(curve_point *cp)
return;
}
// lambda = 3/2 x^2 / y
// lambda = (3 x^2 + a) / (2 y)
lambda = cp->y;
bn_inverse(&lambda, &prime256k1);
bn_multiply(&cp->x, &lambda, &prime256k1);
bn_multiply(&cp->x, &lambda, &prime256k1);
bn_mult_3_2(&lambda, &prime256k1);
bn_mult_k(&lambda, 2, &curve->prime);
bn_inverse(&lambda, &curve->prime);
xr = cp->x;
bn_multiply(&xr, &xr, &curve->prime);
bn_mult_k(&xr, 3, &curve->prime);
bn_addmod(&xr, &curve->a, &curve->prime);
bn_multiply(&xr, &lambda, &curve->prime);
// xr = lambda^2 - 2*x
xr = lambda;
bn_multiply(&xr, &xr, &prime256k1);
temp = 1;
for (i = 0; i < 9; i++) {
temp += 0x3FFFFFFF + xr.val[i] + 2u * (prime256k1.val[i] - cp->x.val[i]);
xr.val[i] = temp & 0x3FFFFFFF;
temp >>= 30;
}
bn_fast_mod(&xr, &prime256k1);
bn_mod(&xr, &prime256k1);
bn_multiply(&xr, &xr, &curve->prime);
yr = cp->x;
bn_lshift(&yr);
bn_subtractmod(&xr, &yr, &xr, &curve->prime);
bn_fast_mod(&xr, &curve->prime);
bn_mod(&xr, &curve->prime);
// yr = lambda (x - xr) - y
bn_subtractmod(&(cp->x), &xr, &yr, &prime256k1);
bn_multiply(&lambda, &yr, &prime256k1);
bn_subtractmod(&yr, &(cp->y), &yr, &prime256k1);
bn_fast_mod(&yr, &prime256k1);
bn_mod(&yr, &prime256k1);
bn_subtractmod(&(cp->x), &xr, &yr, &curve->prime);
bn_multiply(&lambda, &yr, &curve->prime);
bn_subtractmod(&yr, &(cp->y), &yr, &curve->prime);
bn_fast_mod(&yr, &curve->prime);
bn_mod(&yr, &curve->prime);
cp->x = xr;
cp->y = yr;
@ -176,7 +173,7 @@ int point_is_negative_of(const curve_point *p, const curve_point *q)
// Negate a (modulo prime) if cond is 0xffffffff, keep it if cond is 0.
// The timing of this function does not depend on cond.
static void conditional_negate(uint32_t cond, bignum256 *a, const bignum256 *prime)
void conditional_negate(uint32_t cond, bignum256 *a, const bignum256 *prime)
{
int j;
uint32_t tmp = 1;
@ -193,7 +190,7 @@ typedef struct jacobian_curve_point {
bignum256 x, y, z;
} jacobian_curve_point;
static void curve_to_jacobian(const curve_point *p, jacobian_curve_point *jp) {
void curve_to_jacobian(const curve_point *p, jacobian_curve_point *jp, const bignum256 *prime) {
int i;
// randomize z coordinate
for (i = 0; i < 8; i++) {
@ -202,41 +199,39 @@ static void curve_to_jacobian(const curve_point *p, jacobian_curve_point *jp) {
jp->z.val[8] = (random32() & 0x7fff) + 1;
jp->x = jp->z;
bn_multiply(&jp->z, &jp->x, &prime256k1);
bn_multiply(&jp->z, &jp->x, prime);
// x = z^2
jp->y = jp->x;
bn_multiply(&jp->z, &jp->y, &prime256k1);
bn_multiply(&jp->z, &jp->y, prime);
// y = z^3
bn_multiply(&p->x, &jp->x, &prime256k1);
bn_multiply(&p->y, &jp->y, &prime256k1);
bn_mod(&jp->x, &prime256k1);
bn_mod(&jp->y, &prime256k1);
bn_multiply(&p->x, &jp->x, prime);
bn_multiply(&p->y, &jp->y, prime);
bn_mod(&jp->x, prime);
bn_mod(&jp->y, prime);
}
static void jacobian_to_curve(const jacobian_curve_point *jp, curve_point *p) {
void jacobian_to_curve(const jacobian_curve_point *jp, curve_point *p, const bignum256 *prime) {
p->y = jp->z;
bn_mod(&p->y, &prime256k1);
bn_inverse(&p->y, &prime256k1);
bn_mod(&p->y, prime);
bn_inverse(&p->y, prime);
// p->y = z^-1
p->x = p->y;
bn_multiply(&p->x, &p->x, &prime256k1);
bn_multiply(&p->x, &p->x, prime);
// p->x = z^-2
bn_multiply(&p->x, &p->y, &prime256k1);
bn_multiply(&p->x, &p->y, prime);
// p->y = z^-3
bn_multiply(&jp->x, &p->x, &prime256k1);
bn_multiply(&jp->x, &p->x, prime);
// p->x = jp->x * z^-2
bn_multiply(&jp->y, &p->y, &prime256k1);
bn_multiply(&jp->y, &p->y, prime);
// p->y = jp->y * z^-3
bn_mod(&p->x, &prime256k1);
bn_mod(&p->y, &prime256k1);
bn_mod(&p->x, prime);
bn_mod(&p->y, prime);
}
static void point_jacobian_add(const curve_point *p1, jacobian_curve_point *p2) {
void point_jacobian_add(const curve_point *p1, jacobian_curve_point *p2, const bignum256 *prime) {
bignum256 r, h;
bignum256 rsq, hcb, hcby2, hsqx2;
int j;
uint64_t tmp1;
/* usual algorithm:
*
@ -270,75 +265,70 @@ static void point_jacobian_add(const curve_point *p1, jacobian_curve_point *p2)
// h = x1 * z2^2 - x2;
// r = y1 * z2^3 - y2;
h = p2->z;
bn_multiply(&h, &h, &prime256k1); // h = z2^2
bn_multiply(&h, &h, prime); // h = z2^2
r = p2->z;
bn_multiply(&h, &r, &prime256k1); // r = z2^3
bn_multiply(&h, &r, prime); // r = z2^3
bn_multiply(&p1->x, &h, &prime256k1);
bn_subtractmod(&h, &p2->x, &h, &prime256k1);
bn_multiply(&p1->x, &h, prime);
bn_subtractmod(&h, &p2->x, &h, prime);
// h = x1 * z2^2 - x2;
bn_multiply(&p1->y, &r, &prime256k1);
bn_subtractmod(&r, &p2->y, &r, &prime256k1);
bn_multiply(&p1->y, &r, prime);
bn_subtractmod(&r, &p2->y, &r, prime);
// r = y1 * z2^3 - y2;
// hsqx2 = h^2
hsqx2 = h;
bn_multiply(&hsqx2, &hsqx2, &prime256k1);
bn_multiply(&hsqx2, &hsqx2, prime);
// hcb = h^3
hcb = h;
bn_multiply(&hsqx2, &hcb, &prime256k1);
bn_multiply(&hsqx2, &hcb, prime);
// hsqx2 = h^2 * x2
bn_multiply(&p2->x, &hsqx2, &prime256k1);
bn_multiply(&p2->x, &hsqx2, prime);
// hcby2 = h^3 * y2
hcby2 = hcb;
bn_multiply(&p2->y, &hcby2, &prime256k1);
bn_multiply(&p2->y, &hcby2, prime);
// rsq = r^2
rsq = r;
bn_multiply(&rsq, &rsq, &prime256k1);
bn_multiply(&rsq, &rsq, prime);
// z3 = h*z2
bn_multiply(&h, &p2->z, &prime256k1);
bn_mod(&p2->z, &prime256k1);
bn_multiply(&h, &p2->z, prime);
bn_mod(&p2->z, prime);
// x3 = r^2 - h^3 - 2h^2x2
tmp1 = 0;
for (j = 0; j < 9; j++) {
tmp1 += (uint64_t) rsq.val[j] + 4*prime256k1.val[j] - hcb.val[j] - 2*hsqx2.val[j];
assert(tmp1 < 5 * 0x40000000ull);
p2->x.val[j] = tmp1 & 0x3fffffff;
tmp1 >>= 30;
}
bn_fast_mod(&p2->x, &prime256k1);
bn_mod(&p2->x, &prime256k1);
bn_addmod(&hcb, &hsqx2, prime);
bn_addmod(&hcb, &hsqx2, prime);
bn_subtractmod(&rsq, &hcb, &p2->x, prime);
bn_fast_mod(&p2->x, prime);
bn_mod(&p2->x, prime);
// y3 = r*(h^2x2 - x3) - y2*h^3
bn_subtractmod(&hsqx2, &p2->x, &p2->y, &prime256k1);
bn_multiply(&r, &p2->y, &prime256k1);
bn_subtractmod(&p2->y, &hcby2, &p2->y, &prime256k1);
bn_fast_mod(&p2->y, &prime256k1);
bn_mod(&p2->y, &prime256k1);
bn_subtractmod(&hsqx2, &p2->x, &p2->y, prime);
bn_multiply(&r, &p2->y, prime);
bn_subtractmod(&p2->y, &hcby2, &p2->y, prime);
bn_fast_mod(&p2->y, prime);
bn_mod(&p2->y, prime);
}
static void point_jacobian_double(jacobian_curve_point *p) {
bignum256 m, msq, ysq, xysq;
int j;
uint32_t tmp1;
void point_jacobian_double(jacobian_curve_point *p, const ecdsa_curve *curve) {
bignum256 az4, m, msq, ysq, xysq;
const bignum256 *prime = &curve->prime;
/* usual algorithm:
*
* lambda = (3(x/z^2)^2 / 2y/z^3) = 3x^2/2yz
* lambda = (3((x/z^2)^2 + a) / 2y/z^3) = (3x^2 + az^4)/2yz
* x3/z3^2 = lambda^2 - 2x/z^2
* y3/z3^3 = lambda * (x/z^2 - x3/z3^2) - y/z^3
*
* to get rid of fraction we set
* m = 3/2 x^2
* m = (3 x^2 + az^4) / 2
* Hence,
* lambda = m / yz
* lambda = m / yz = m / z3
*
* With z3 = yz (the denominator of lambda)
* we get x3 = lambda^2*z3^2 - 2*x/z^2*z3^2
@ -347,57 +337,62 @@ static void point_jacobian_double(jacobian_curve_point *p) {
* = m * (xy^2 - x3) - y^4
*/
/* m = 3/2*x*x
/* m = (3*x^2 + a z^4) / 2
* x3 = m^2 - 2*xy^2
* y3 = m*(xy^2 - x3) - 8y^4
* z3 = y*z
*/
m = p->x;
bn_multiply(&m, &m, &prime256k1);
bn_mult_3_2(&m, &prime256k1);
bn_multiply(&m, &m, prime);
bn_mult_k(&m, 3, prime);
az4 = p->z;
bn_multiply(&az4, &az4, prime);
bn_multiply(&az4, &az4, prime);
bn_multiply(&curve->a, &az4, prime);
bn_addmod(&m, &az4, prime);
bn_mult_half(&m, prime);
// msq = m^2
msq = m;
bn_multiply(&msq, &msq, &prime256k1);
bn_multiply(&msq, &msq, prime);
// ysq = y^2
ysq = p->y;
bn_multiply(&ysq, &ysq, &prime256k1);
bn_multiply(&ysq, &ysq, prime);
// xysq = xy^2
xysq = p->x;
bn_multiply(&ysq, &xysq, &prime256k1);
bn_multiply(&ysq, &xysq, prime);
// z3 = yz
bn_multiply(&p->y, &p->z, &prime256k1);
bn_mod(&p->z, &prime256k1);
bn_multiply(&p->y, &p->z, prime);
bn_mod(&p->z, prime);
// x3 = m^2 - 2*xy^2
tmp1 = 0;
for (j = 0; j < 9; j++) {
tmp1 += msq.val[j] + 3*prime256k1.val[j] - 2*xysq.val[j];
p->x.val[j] = tmp1 & 0x3fffffff;
tmp1 >>= 30;
}
bn_fast_mod(&p->x, &prime256k1);
bn_mod(&p->x, &prime256k1);
p->x = xysq;
bn_mod(&p->x, prime);
bn_lshift(&p->x);
bn_subtractmod(&msq, &p->x, &p->x, prime);
bn_fast_mod(&p->x, prime);
bn_mod(&p->x, prime);
// y3 = m*(xy^2 - x3) - y^4
bn_subtractmod(&xysq, &p->x, &p->y, &prime256k1);
bn_multiply(&m, &p->y, &prime256k1);
bn_multiply(&ysq, &ysq, &prime256k1);
bn_subtractmod(&p->y, &ysq, &p->y, &prime256k1);
bn_fast_mod(&p->y, &prime256k1);
bn_mod(&p->y, &prime256k1);
bn_subtractmod(&xysq, &p->x, &p->y, prime);
bn_multiply(&m, &p->y, prime);
bn_multiply(&ysq, &ysq, prime);
bn_subtractmod(&p->y, &ysq, &p->y, prime);
bn_fast_mod(&p->y, prime);
bn_mod(&p->y, prime);
}
// res = k * p
void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
void point_multiply(const ecdsa_curve *curve, const bignum256 *k, const curve_point *p, curve_point *res)
{
// this algorithm is loosely based on
// Katsuyuki Okeya and Tsuyoshi Takagi, The Width-w NAF Method Provides
// Small Memory and Fast Elliptic Scalar Multiplications Secure against
// Side Channel Attacks.
assert (bn_is_less(k, &order256k1));
assert (bn_is_less(k, &curve->order));
int i, j;
int pos, shift;
@ -406,21 +401,22 @@ void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
uint32_t bits, sign, nsign;
jacobian_curve_point jres;
curve_point pmult[8];
const bignum256 *prime = &curve->prime;
// is_even = 0xffffffff if k is even, 0 otherwise.
// add 2^256.
// make number odd: subtract order256k1 if even
// make number odd: subtract curve->order if even
uint32_t tmp = 1;
uint32_t is_non_zero = 0;
for (j = 0; j < 8; j++) {
is_non_zero |= k->val[j];
tmp += 0x3fffffff + k->val[j] - (order256k1.val[j] & is_even);
tmp += 0x3fffffff + k->val[j] - (curve->order.val[j] & is_even);
a.val[j] = tmp & 0x3fffffff;
tmp >>= 30;
}
is_non_zero |= k->val[j];
a.val[j] = tmp + 0xffff + k->val[j] - (order256k1.val[j] & is_even);
a.val[j] = tmp + 0xffff + k->val[j] - (curve->order.val[j] & is_even);
assert((a.val[0] & 1) != 0);
// special case 0*p: just return zero. We don't care about constant time.
@ -429,7 +425,7 @@ void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
return;
}
// Now a = k + 2^256 (mod order256k1) and a is odd.
// Now a = k + 2^256 (mod curve->order) and a is odd.
//
// The idea is to bring the new a into the form.
// sum_{i=0..64} a[i] 16^i, where |a[i]| < 16 and a[i] is odd.
@ -437,7 +433,7 @@ void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
// add 1 to it and subtract 16 from a[i-1]. Afterwards,
// a[64] = 1, which is the 2^256 that we added before.
//
// Since k = a - 2^256 (mod order256k1), we can compute
// Since k = a - 2^256 (mod curve->order), we can compute
// k*p = sum_{i=0..63} a[i] 16^i * p
//
// We compute |a[i]| * p in advance for all possible
@ -445,12 +441,12 @@ void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
// We compute p, 3*p, ..., 15*p and store it in the table pmult.
// store p^2 temporarily in pmult[7]
pmult[7] = *p;
point_double(&pmult[7]);
point_double(curve, &pmult[7]);
// compute 3*p, etc by repeatedly adding p^2.
pmult[0] = *p;
for (i = 1; i < 8; i++) {
pmult[i] = pmult[7];
point_add(&pmult[i-1], &pmult[i]);
point_add(curve, &pmult[i-1], &pmult[i]);
}
// now compute res = sum_{i=0..63} a[i] * 16^i * p step by step,
@ -464,15 +460,15 @@ void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
sign = (bits >> 4) - 1;
bits ^= sign;
bits &= 15;
curve_to_jacobian(&pmult[bits>>1], &jres);
curve_to_jacobian(&pmult[bits>>1], &jres, prime);
for (i = 62; i >= 0; i--) {
// sign = sign(a[i+1]) (0xffffffff for negative, 0 for positive)
// invariant jres = (-1)^sign sum_{j=i+1..63} (a[j] * 16^{j-i-1} * p)
point_jacobian_double(&jres);
point_jacobian_double(&jres);
point_jacobian_double(&jres);
point_jacobian_double(&jres);
point_jacobian_double(&jres, curve);
point_jacobian_double(&jres, curve);
point_jacobian_double(&jres, curve);
point_jacobian_double(&jres, curve);
// get lowest 5 bits of a >> (i*4).
pos = i*4/30; shift = i*4 % 30;
@ -483,44 +479,45 @@ void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res)
// negate last result to make signs of this round and the
// last round equal.
conditional_negate(sign ^ nsign, &jres.z, &prime256k1);
conditional_negate(sign ^ nsign, &jres.z, prime);
// add odd factor
point_jacobian_add(&pmult[bits >> 1], &jres);
point_jacobian_add(&pmult[bits >> 1], &jres, prime);
sign = nsign;
}
conditional_negate(sign, &jres.z, &prime256k1);
jacobian_to_curve(&jres, res);
conditional_negate(sign, &jres.z, prime);
jacobian_to_curve(&jres, res, prime);
}
#if USE_PRECOMPUTED_CP
// res = k * G
// k must be a normalized number with 0 <= k < order256k1
void scalar_multiply(const bignum256 *k, curve_point *res)
// k must be a normalized number with 0 <= k < curve->order
void scalar_multiply(const ecdsa_curve *curve, const bignum256 *k, curve_point *res)
{
assert (bn_is_less(k, &order256k1));
assert (bn_is_less(k, &curve->order));
int i, j;
bignum256 a;
uint32_t is_even = (k->val[0] & 1) - 1;
uint32_t lowbits;
jacobian_curve_point jres;
const bignum256 *prime = &curve->prime;
// is_even = 0xffffffff if k is even, 0 otherwise.
// add 2^256.
// make number odd: subtract order256k1 if even
// make number odd: subtract curve->order if even
uint32_t tmp = 1;
uint32_t is_non_zero = 0;
for (j = 0; j < 8; j++) {
is_non_zero |= k->val[j];
tmp += 0x3fffffff + k->val[j] - (order256k1.val[j] & is_even);
tmp += 0x3fffffff + k->val[j] - (curve->order.val[j] & is_even);
a.val[j] = tmp & 0x3fffffff;
tmp >>= 30;
}
is_non_zero |= k->val[j];
a.val[j] = tmp + 0xffff + k->val[j] - (order256k1.val[j] & is_even);
a.val[j] = tmp + 0xffff + k->val[j] - (curve->order.val[j] & is_even);
assert((a.val[0] & 1) != 0);
// special case 0*G: just return zero. We don't care about constant time.
@ -529,7 +526,7 @@ void scalar_multiply(const bignum256 *k, curve_point *res)
return;
}
// Now a = k + 2^256 (mod order256k1) and a is odd.
// Now a = k + 2^256 (mod curve->order) and a is odd.
//
// The idea is to bring the new a into the form.
// sum_{i=0..64} a[i] 16^i, where |a[i]| < 16 and a[i] is odd.
@ -537,12 +534,12 @@ void scalar_multiply(const bignum256 *k, curve_point *res)
// add 1 to it and subtract 16 from a[i-1]. Afterwards,
// a[64] = 1, which is the 2^256 that we added before.
//
// Since k = a - 2^256 (mod order256k1), we can compute
// Since k = a - 2^256 (mod curve->order), we can compute
// k*G = sum_{i=0..63} a[i] 16^i * G
//
// We have a big table secp256k1_cp that stores all possible
// We have a big table curve->cp that stores all possible
// values of |a[i]| 16^i * G.
// secp256k1_cp[i][j] = (2*j+1) * 16^i * G
// curve->cp[i][j] = (2*j+1) * 16^i * G
// now compute res = sum_{i=0..63} a[i] * 16^i * G step by step.
// initial res = |a[0]| * G. Note that a[0] = a & 0xf if (a&0x10) != 0
@ -552,7 +549,7 @@ void scalar_multiply(const bignum256 *k, curve_point *res)
lowbits = a.val[0] & ((1 << 5) - 1);
lowbits ^= (lowbits >> 4) - 1;
lowbits &= 15;
curve_to_jacobian(&secp256k1_cp[0][lowbits >> 1], &jres);
curve_to_jacobian(&curve->cp[0][lowbits >> 1], &jres, prime);
for (i = 1; i < 64; i ++) {
// invariant res = sign(a[i-1]) sum_{j=0..i-1} (a[j] * 16^j * G)
@ -569,26 +566,26 @@ void scalar_multiply(const bignum256 *k, curve_point *res)
lowbits &= 15;
// negate last result to make signs of this round and the
// last round equal.
conditional_negate((lowbits & 1) - 1, &jres.y, &prime256k1);
conditional_negate((lowbits & 1) - 1, &jres.y, prime);
// add odd factor
point_jacobian_add(&secp256k1_cp[i][lowbits >> 1], &jres);
point_jacobian_add(&curve->cp[i][lowbits >> 1], &jres, prime);
}
conditional_negate(((a.val[0] >> 4) & 1) - 1, &jres.y, &prime256k1);
jacobian_to_curve(&jres, res);
conditional_negate(((a.val[0] >> 4) & 1) - 1, &jres.y, prime);
jacobian_to_curve(&jres, res, prime);
}
#else
void scalar_multiply(const bignum256 *k, curve_point *res)
void scalar_multiply(const ecdsa_curve *curve, const bignum256 *k, curve_point *res)
{
point_multiply(k, &G256k1, res);
point_multiply(curve, k, &curve->G, res);
}
#endif
// generate random K for signing
int generate_k_random(bignum256 *k) {
int generate_k_random(const ecdsa_curve *curve, bignum256 *k) {
int i, j;
for (j = 0; j < 10000; j++) {
for (i = 0; i < 8; i++) {
@ -596,7 +593,7 @@ int generate_k_random(bignum256 *k) {
}
k->val[8] = random32() & 0xFFFF;
// if k is too big or too small, we don't like it
if ( !bn_is_zero(k) && bn_is_less(k, &order256k1) ) {
if ( !bn_is_zero(k) && bn_is_less(k, &curve->order) ) {
return 0; // good number - no error
}
}
@ -606,7 +603,7 @@ int generate_k_random(bignum256 *k) {
// generate K in a deterministic way, according to RFC6979
// http://tools.ietf.org/html/rfc6979
int generate_k_rfc6979(bignum256 *secret, const uint8_t *priv_key, const uint8_t *hash)
int generate_k_rfc6979(const ecdsa_curve *curve, bignum256 *secret, const uint8_t *priv_key, const uint8_t *hash)
{
int i;
uint8_t v[32], k[32], bx[2*32], buf[32 + 1 + sizeof(bx)];
@ -614,7 +611,7 @@ int generate_k_rfc6979(bignum256 *secret, const uint8_t *priv_key, const uint8_t
memcpy(bx, priv_key, 32);
bn_read_be(hash, &z1);
bn_mod(&z1, &order256k1);
bn_mod(&z1, &curve->order);
bn_write_be(&z1, bx + 32);
memset(v, 1, sizeof(v));
@ -635,7 +632,7 @@ int generate_k_rfc6979(bignum256 *secret, const uint8_t *priv_key, const uint8_t
for (i = 0; i < 10000; i++) {
hmac_sha256(k, sizeof(k), v, sizeof(v), v);
bn_read_be(v, secret);
if ( !bn_is_zero(secret) && bn_is_less(secret, &order256k1) ) {
if ( !bn_is_zero(secret) && bn_is_less(secret, &curve->order) ) {
return 0; // good number -> no error
}
memcpy(buf, v, sizeof(v));
@ -649,11 +646,11 @@ int generate_k_rfc6979(bignum256 *secret, const uint8_t *priv_key, const uint8_t
// msg is a data to be signed
// msg_len is the message length
int ecdsa_sign(const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby)
int ecdsa_sign(const ecdsa_curve *curve, const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby)
{
uint8_t hash[32];
sha256_Raw(msg, msg_len, hash);
int res = ecdsa_sign_digest(priv_key, hash, sig, pby);
int res = ecdsa_sign_digest(curve, priv_key, hash, sig, pby);
MEMSET_BZERO(hash, sizeof(hash));
return res;
@ -661,12 +658,12 @@ int ecdsa_sign(const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, ui
// msg is a data to be signed
// msg_len is the message length
int ecdsa_sign_double(const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby)
int ecdsa_sign_double(const ecdsa_curve *curve, const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby)
{
uint8_t hash[32];
sha256_Raw(msg, msg_len, hash);
sha256_Raw(hash, 32, hash);
int res = ecdsa_sign_digest(priv_key, hash, sig, pby);
int res = ecdsa_sign_digest(curve, priv_key, hash, sig, pby);
MEMSET_BZERO(hash, sizeof(hash));
return res;
}
@ -675,7 +672,7 @@ int ecdsa_sign_double(const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_
// priv_key is a 32 byte big endian stored number
// sig is 64 bytes long array for the signature
// digest is 32 bytes of digest
int ecdsa_sign_digest(const uint8_t *priv_key, const uint8_t *digest, uint8_t *sig, uint8_t *pby)
int ecdsa_sign_digest(const ecdsa_curve *curve, const uint8_t *priv_key, const uint8_t *digest, uint8_t *sig, uint8_t *pby)
{
uint32_t i;
curve_point R;
@ -686,24 +683,24 @@ int ecdsa_sign_digest(const uint8_t *priv_key, const uint8_t *digest, uint8_t *s
#if USE_RFC6979
// generate K deterministically
if (generate_k_rfc6979(&k, priv_key, digest) != 0) {
if (generate_k_rfc6979(curve, &k, priv_key, digest) != 0) {
result = 1;
}
#else
// generate random number k
if (generate_k_random(&k) != 0) {
if (generate_k_random(curve, &k) != 0) {
result = 1;
}
#endif
if (result == 0) {
// compute k*G
scalar_multiply(&k, &R);
scalar_multiply(curve, &k, &R);
if (pby) {
*pby = R.y.val[0] & 1;
}
// r = (rx mod n)
bn_mod(&R.x, &order256k1);
bn_mod(&R.x, &curve->order);
// if r is zero, we fail
if (bn_is_zero(&R.x))
{
@ -712,17 +709,17 @@ int ecdsa_sign_digest(const uint8_t *priv_key, const uint8_t *digest, uint8_t *s
}
if (result == 0) {
bn_inverse(&k, &order256k1);
bn_inverse(&k, &curve->order);
bn_read_be(priv_key, da);
bn_multiply(&R.x, da, &order256k1);
bn_multiply(&R.x, da, &curve->order);
for (i = 0; i < 8; i++) {
da->val[i] += z.val[i];
da->val[i + 1] += (da->val[i] >> 30);
da->val[i] &= 0x3FFFFFFF;
}
da->val[8] += z.val[8];
bn_multiply(da, &k, &order256k1);
bn_mod(&k, &order256k1);
bn_multiply(da, &k, &curve->order);
bn_mod(&k, &curve->order);
// if k is zero, we fail
if (bn_is_zero(&k)) {
result = 3;
@ -730,13 +727,6 @@ int ecdsa_sign_digest(const uint8_t *priv_key, const uint8_t *digest, uint8_t *s
}
if (result == 0) {
// if S > order/2 => S = -S
if (bn_is_less(&order256k1_half, &k)) {
bn_subtract(&order256k1, &k, &k);
if (pby) {
*pby = !*pby;
}
}
// we are done, R.x and k is the result signature
bn_write_be(&R.x, sig);
bn_write_be(&k, sig + 32);
@ -748,28 +738,28 @@ int ecdsa_sign_digest(const uint8_t *priv_key, const uint8_t *digest, uint8_t *s
return 0;
}
void ecdsa_get_public_key33(const uint8_t *priv_key, uint8_t *pub_key)
void ecdsa_get_public_key33(const ecdsa_curve *curve, const uint8_t *priv_key, uint8_t *pub_key)
{
curve_point R;
bignum256 k;
bn_read_be(priv_key, &k);
// compute k*G
scalar_multiply(&k, &R);
scalar_multiply(curve, &k, &R);
pub_key[0] = 0x02 | (R.y.val[0] & 0x01);
bn_write_be(&R.x, pub_key + 1);
MEMSET_BZERO(&R, sizeof(R));
MEMSET_BZERO(&k, sizeof(k));
}
void ecdsa_get_public_key65(const uint8_t *priv_key, uint8_t *pub_key)
void ecdsa_get_public_key65(const ecdsa_curve *curve, const uint8_t *priv_key, uint8_t *pub_key)
{
curve_point R;
bignum256 k;
bn_read_be(priv_key, &k);
// compute k*G
scalar_multiply(&k, &R);
scalar_multiply(curve, &k, &R);
pub_key[0] = 0x04;
bn_write_be(&R.x, pub_key + 1);
bn_write_be(&R.y, pub_key + 33);
@ -825,30 +815,30 @@ int ecdsa_address_decode(const char *addr, uint8_t *out)
return base58_decode_check(addr, out, 21) == 21;
}
void uncompress_coords(uint8_t odd, const bignum256 *x, bignum256 *y)
void uncompress_coords(const ecdsa_curve *curve, uint8_t odd, const bignum256 *x, bignum256 *y)
{
// y^2 = x^3 + 0*x + 7
memcpy(y, x, sizeof(bignum256)); // y is x
bn_multiply(x, y, &prime256k1); // y is x^2
bn_multiply(x, y, &prime256k1); // y is x^3
bn_addmodi(y, 7, &prime256k1); // y is x^3 + 7
bn_sqrt(y, &prime256k1); // y = sqrt(y)
bn_multiply(x, y, &curve->prime); // y is x^2
bn_multiply(x, y, &curve->prime); // y is x^3
bn_addmodi(y, 7, &curve->prime); // y is x^3 + 7
bn_sqrt(y, &curve->prime); // y = sqrt(y)
if ((odd & 0x01) != (y->val[0] & 1)) {
bn_subtract(&prime256k1, y, y); // y = -y
bn_subtract(&curve->prime, y, y); // y = -y
}
}
int ecdsa_read_pubkey(const uint8_t *pub_key, curve_point *pub)
int ecdsa_read_pubkey(const ecdsa_curve *curve, const uint8_t *pub_key, curve_point *pub)
{
if (pub_key[0] == 0x04) {
bn_read_be(pub_key + 1, &(pub->x));
bn_read_be(pub_key + 33, &(pub->y));
return ecdsa_validate_pubkey(pub);
return ecdsa_validate_pubkey(curve, pub);
}
if (pub_key[0] == 0x02 || pub_key[0] == 0x03) { // compute missing y coords
bn_read_be(pub_key + 1, &(pub->x));
uncompress_coords(pub_key[0], &(pub->x), &(pub->y));
return ecdsa_validate_pubkey(pub);
uncompress_coords(curve, pub_key[0], &(pub->x), &(pub->y));
return ecdsa_validate_pubkey(curve, pub);
}
// error
return 0;
@ -859,7 +849,7 @@ int ecdsa_read_pubkey(const uint8_t *pub_key, curve_point *pub)
// - pub->x and pub->y are in range [0,p-1].
// - pub is on the curve.
int ecdsa_validate_pubkey(const curve_point *pub)
int ecdsa_validate_pubkey(const ecdsa_curve *curve, const curve_point *pub)
{
bignum256 y_2, x_3_b;
@ -867,7 +857,7 @@ int ecdsa_validate_pubkey(const curve_point *pub)
return 0;
}
if (!bn_is_less(&(pub->x), &prime256k1) || !bn_is_less(&(pub->y), &prime256k1)) {
if (!bn_is_less(&(pub->x), &curve->prime) || !bn_is_less(&(pub->y), &curve->prime)) {
return 0;
}
@ -875,13 +865,13 @@ int ecdsa_validate_pubkey(const curve_point *pub)
memcpy(&x_3_b, &(pub->x), sizeof(bignum256));
// y^2
bn_multiply(&(pub->y), &y_2, &prime256k1);
bn_mod(&y_2, &prime256k1);
bn_multiply(&(pub->y), &y_2, &curve->prime);
bn_mod(&y_2, &curve->prime);
// x^3 + b
bn_multiply(&(pub->x), &x_3_b, &prime256k1);
bn_multiply(&(pub->x), &x_3_b, &prime256k1);
bn_addmodi(&x_3_b, 7, &prime256k1);
bn_multiply(&(pub->x), &x_3_b, &curve->prime);
bn_multiply(&(pub->x), &x_3_b, &curve->prime);
bn_addmodi(&x_3_b, 7, &curve->prime);
if (!bn_is_equal(&x_3_b, &y_2)) {
return 0;
@ -896,32 +886,32 @@ int ecdsa_validate_pubkey(const curve_point *pub)
// msg is a data that was signed
// msg_len is the message length
int ecdsa_verify(const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len)
int ecdsa_verify(const ecdsa_curve *curve, const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len)
{
uint8_t hash[32];
sha256_Raw(msg, msg_len, hash);
int res = ecdsa_verify_digest(pub_key, sig, hash);
int res = ecdsa_verify_digest(curve, pub_key, sig, hash);
MEMSET_BZERO(hash, sizeof(hash));
return res;
}
int ecdsa_verify_double(const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len)
int ecdsa_verify_double(const ecdsa_curve *curve, const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len)
{
uint8_t hash[32];
sha256_Raw(msg, msg_len, hash);
sha256_Raw(hash, 32, hash);
int res = ecdsa_verify_digest(pub_key, sig, hash);
int res = ecdsa_verify_digest(curve, pub_key, sig, hash);
MEMSET_BZERO(hash, sizeof(hash));
return res;
}
// returns 0 if verification succeeded
int ecdsa_verify_digest(const uint8_t *pub_key, const uint8_t *sig, const uint8_t *digest)
int ecdsa_verify_digest(const ecdsa_curve *curve, const uint8_t *pub_key, const uint8_t *sig, const uint8_t *digest)
{
curve_point pub, res;
bignum256 r, s, z;
if (!ecdsa_read_pubkey(pub_key, &pub)) {
if (!ecdsa_read_pubkey(curve, pub_key, &pub)) {
return 1;
}
@ -931,14 +921,14 @@ int ecdsa_verify_digest(const uint8_t *pub_key, const uint8_t *sig, const uint8_
bn_read_be(digest, &z);
if (bn_is_zero(&r) || bn_is_zero(&s) ||
(!bn_is_less(&r, &order256k1)) ||
(!bn_is_less(&s, &order256k1))) return 2;
(!bn_is_less(&r, &curve->order)) ||
(!bn_is_less(&s, &curve->order))) return 2;
bn_inverse(&s, &order256k1); // s^-1
bn_multiply(&s, &z, &order256k1); // z*s^-1
bn_mod(&z, &order256k1);
bn_multiply(&r, &s, &order256k1); // r*s^-1
bn_mod(&s, &order256k1);
bn_inverse(&s, &curve->order); // s^-1
bn_multiply(&s, &z, &curve->order); // z*s^-1
bn_mod(&z, &curve->order);
bn_multiply(&r, &s, &curve->order); // r*s^-1
bn_mod(&s, &curve->order);
int result = 0;
if (bn_is_zero(&z)) {
@ -946,14 +936,14 @@ int ecdsa_verify_digest(const uint8_t *pub_key, const uint8_t *sig, const uint8_
// I don't expect this to happen any time soon
result = 3;
} else {
scalar_multiply(&z, &res);
scalar_multiply(curve, &z, &res);
}
if (result == 0) {
// both pub and res can be infinity, can have y = 0 OR can be equal -> false negative
point_multiply(&s, &pub, &pub);
point_add(&pub, &res);
bn_mod(&(res.x), &order256k1);
point_multiply(curve, &s, &pub, &pub);
point_add(curve, &pub, &res);
bn_mod(&(res.x), &curve->order);
// signature does not match
if (!bn_is_equal(&res.x, &r)) {
result = 5;
@ -1006,3 +996,17 @@ int ecdsa_sig_to_der(const uint8_t *sig, uint8_t *der)
*len = *len1 + *len2 + 4;
return *len + 2;
}
const ecdsa_curve *get_curve_by_name(const char *curve_name) {
if (curve_name == 0) {
return 0;
}
if (strcmp(curve_name, "secp256k1") == 0) {
return &secp256k1;
}
if (strcmp(curve_name, "nist256p1") == 0) {
return &nist256p1;
}
return 0;
}

@ -26,38 +26,57 @@
#include <stdint.h>
#include "options.h"
#include "secp256k1.h"
#include "bignum.h"
// curve point x and y
typedef struct {
bignum256 x, y;
} curve_point;
typedef struct {
bignum256 prime; // prime order of the finite field
curve_point G; // initial curve point
bignum256 order; // order of G
bignum256 a; // coefficient 'a' of the elliptic curve
#if USE_PRECOMPUTED_CP
const curve_point cp[64][8];
#endif
} ecdsa_curve;
void point_copy(const curve_point *cp1, curve_point *cp2);
void point_add(const curve_point *cp1, curve_point *cp2);
void point_double(curve_point *cp);
void point_multiply(const bignum256 *k, const curve_point *p, curve_point *res);
void point_add(const ecdsa_curve *curve, const curve_point *cp1, curve_point *cp2);
void point_double(const ecdsa_curve *curve, curve_point *cp);
void point_multiply(const ecdsa_curve *curve, const bignum256 *k, const curve_point *p, curve_point *res);
void point_set_infinity(curve_point *p);
int point_is_infinity(const curve_point *p);
int point_is_equal(const curve_point *p, const curve_point *q);
int point_is_negative_of(const curve_point *p, const curve_point *q);
void scalar_multiply(const bignum256 *k, curve_point *res);
void uncompress_coords(uint8_t odd, const bignum256 *x, bignum256 *y);
int ecdsa_sign(const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby);
int ecdsa_sign_double(const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby);
int ecdsa_sign_digest(const uint8_t *priv_key, const uint8_t *digest, uint8_t *sig, uint8_t *pby);
void ecdsa_get_public_key33(const uint8_t *priv_key, uint8_t *pub_key);
void ecdsa_get_public_key65(const uint8_t *priv_key, uint8_t *pub_key);
void scalar_multiply(const ecdsa_curve *curve, const bignum256 *k, curve_point *res);
void uncompress_coords(const ecdsa_curve *curve, uint8_t odd, const bignum256 *x, bignum256 *y);
int ecdsa_sign(const ecdsa_curve *curve, const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby);
int ecdsa_sign_double(const ecdsa_curve *curve, const uint8_t *priv_key, const uint8_t *msg, uint32_t msg_len, uint8_t *sig, uint8_t *pby);
int ecdsa_sign_digest(const ecdsa_curve *curve, const uint8_t *priv_key, const uint8_t *digest, uint8_t *sig, uint8_t *pby);
void ecdsa_get_public_key33(const ecdsa_curve *curve, const uint8_t *priv_key, uint8_t *pub_key);
void ecdsa_get_public_key65(const ecdsa_curve *curve, const uint8_t *priv_key, uint8_t *pub_key);
void ecdsa_get_pubkeyhash(const uint8_t *pub_key, uint8_t *pubkeyhash);
void ecdsa_get_address_raw(const uint8_t *pub_key, uint8_t version, uint8_t *addr_raw);
void ecdsa_get_address(const uint8_t *pub_key, uint8_t version, char *addr, int addrsize);
void ecdsa_get_wif(const uint8_t *priv_key, uint8_t version, char *wif, int wifsize);
int ecdsa_address_decode(const char *addr, uint8_t *out);
int ecdsa_read_pubkey(const uint8_t *pub_key, curve_point *pub);
int ecdsa_validate_pubkey(const curve_point *pub);
int ecdsa_verify(const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len);
int ecdsa_verify_double(const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len);
int ecdsa_verify_digest(const uint8_t *pub_key, const uint8_t *sig, const uint8_t *digest);
int ecdsa_read_pubkey(const ecdsa_curve *curve, const uint8_t *pub_key, curve_point *pub);
int ecdsa_validate_pubkey(const ecdsa_curve *curve, const curve_point *pub);
int ecdsa_verify(const ecdsa_curve *curve, const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len);
int ecdsa_verify_double(const ecdsa_curve *curve, const uint8_t *pub_key, const uint8_t *sig, const uint8_t *msg, uint32_t msg_len);
int ecdsa_verify_digest(const ecdsa_curve *curve, const uint8_t *pub_key, const uint8_t *sig, const uint8_t *digest);
int ecdsa_sig_to_der(const uint8_t *sig, uint8_t *der);
const ecdsa_curve *get_curve_by_name(const char *curve_name);
// Private
int generate_k_rfc6979(bignum256 *secret, const uint8_t *priv_key, const uint8_t *hash);
int generate_k_random(bignum256 *k);
int generate_k_rfc6979(const ecdsa_curve *curve, bignum256 *secret, const uint8_t *priv_key, const uint8_t *hash);
int generate_k_random(const ecdsa_curve *curve, bignum256 *k);
#endif

@ -0,0 +1,50 @@
/**
* Copyright (c) 2013-2014 Tomas Dzetkulic
* Copyright (c) 2013-2014 Pavol Rusnak
*
* Permission is hereby granted, free of charge, to any person obtaining
* a copy of this software and associated documentation files (the "Software"),
* to deal in the Software without restriction, including without limitation
* the rights to use, copy, modify, merge, publish, distribute, sublicense,
* and/or sell copies of the Software, and to permit persons to whom the
* Software is furnished to do so, subject to the following conditions:
*
* The above copyright notice and this permission notice shall be included
* in all copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
* OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
* FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
* THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES
* OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
* ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
* OTHER DEALINGS IN THE SOFTWARE.
*/
#include "nist256p1.h"
const ecdsa_curve nist256p1 = {
/* .prime */ {
/*.val =*/ {0x3fffffff, 0x3fffffff, 0x3fffffff, 0x3f, 0x0, 0x0, 0x1000, 0x3fffc000, 0xffff}
},
/* G */ {
/*.x =*/{/*.val =*/{0x1898c296, 0x1284e517, 0x1eb33a0f, 0xdf604b, 0x2440f277, 0x339b958e, 0x4247f8b, 0x347cb84b, 0x6b17}},
/*.y =*/{/*.val =*/{0x37bf51f5, 0x2ed901a0, 0x3315ecec, 0x338cd5da, 0xf9e162b, 0x1fad29f0, 0x27f9b8ee, 0x10b8bf86, 0x4fe3}}
},
/* order */ {
/*.val =*/{0x3c632551, 0xee72b0b, 0x3179e84f, 0x39beab69, 0x3fffffbc, 0x3fffffff, 0xfff, 0x3fffc000, 0xffff}
},
/* a */ {
/*.val =*/{0x3ffffffc, 0x3fffffff, 0x3fffffff, 0x3f, 0x0, 0x0, 0x1000, 0x3fffc000, 0xffff}
}
#if USE_PRECOMPUTED_CP
,
/* cp */ {
#include "nist256p1.table"
}
#endif
};

@ -0,0 +1,33 @@
/**
* Copyright (c) 2013-2014 Tomas Dzetkulic
* Copyright (c) 2013-2014 Pavol Rusnak
*
* Permission is hereby granted, free of charge, to any person obtaining
* a copy of this software and associated documentation files (the "Software"),
* to deal in the Software without restriction, including without limitation
* the rights to use, copy, modify, merge, publish, distribute, sublicense,
* and/or sell copies of the Software, and to permit persons to whom the
* Software is furnished to do so, subject to the following conditions:
*
* The above copyright notice and this permission notice shall be included
* in all copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
* OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
* FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
* THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES
* OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
* ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
* OTHER DEALINGS IN THE SOFTWARE.
*/
#ifndef __NIST256P1_H__
#define __NIST256P1_H__
#include <stdint.h>
#include "ecdsa.h"
extern const ecdsa_curve nist256p1;
#endif

File diff suppressed because it is too large Load Diff

File diff suppressed because it is too large Load Diff

@ -26,30 +26,8 @@
#include <stdint.h>
#include "bignum.h"
#include "ecdsa.h"
// curve point x and y
typedef struct {
bignum256 x, y;
} curve_point;
// secp256k1 prime
extern const bignum256 prime256k1;
// secp256k1 initial curve point
extern const curve_point G256k1;
// secp256k1 order of G
extern const bignum256 order256k1;
// secp256k1 order of G / 2
extern const bignum256 order256k1_half;
// 3/2 in G_p
extern const bignum256 three_over_two256k1;
#if USE_PRECOMPUTED_CP
extern const curve_point secp256k1_cp[64][8];
#endif
extern const ecdsa_curve secp256k1;
#endif

File diff suppressed because it is too large Load Diff

@ -29,6 +29,13 @@
#include "ecdsa.h"
#include "rand.h"
#include "secp256k1.h"
#define CURVE (&secp256k1)
#define prime256k1 (secp256k1.prime)
#define G256k1 (secp256k1.G)
#define order256k1 (secp256k1.order)
#define secp256k1_cp (secp256k1.cp)
int main(int argc, char *argv[])
{
@ -81,21 +88,21 @@ int main(int argc, char *argv[])
}
// use our ECDSA signer to sign the message with the key
if (ecdsa_sign(priv_key, msg, msg_len, sig, 0) != 0) {
if (ecdsa_sign(CURVE, priv_key, msg, msg_len, sig, 0) != 0) {
printf("trezor-crypto signing failed\n");
break;
}
// generate public key from private key
ecdsa_get_public_key33(priv_key, pub_key33);
ecdsa_get_public_key65(priv_key, pub_key65);
ecdsa_get_public_key33(&secp256k1, priv_key, pub_key33);
ecdsa_get_public_key65(&secp256k1, priv_key, pub_key65);
// use our ECDSA verifier to verify the message signature
if (ecdsa_verify(pub_key65, sig, msg, msg_len) != 0) {
if (ecdsa_verify(CURVE, pub_key65, sig, msg, msg_len) != 0) {
printf("trezor-crypto verification failed (pub_key_len = 65)\n");
break;
}
if (ecdsa_verify(pub_key33, sig, msg, msg_len) != 0) {
if (ecdsa_verify(CURVE, pub_key33, sig, msg, msg_len) != 0) {
printf("trezor-crypto verification failed (pub_key_len = 33)\n");
break;
}

@ -0,0 +1,304 @@
import ctypes as c
import random
import ecdsa
import hashlib
import subprocess
import binascii
import pytest
import os
def bytes2num(s):
res = 0
for i, b in enumerate(reversed(bytearray(s))):
res += b << (i * 8)
return res
curves = {
'nist256p1': ecdsa.curves.NIST256p,
'secp256k1': ecdsa.curves.SECP256k1
}
random_iters = int(os.environ.get('ITERS', 1))
scons_file = '''
srcs = 'ecdsa bignum secp256k1 nist256p1 sha2 rand hmac ripemd160 base58'
srcs = [(s + '.c') for s in srcs.split()]
flags = ('-Os -g -W -Wall -Wextra -Wimplicit-function-declaration '
'-Wredundant-decls -Wstrict-prototypes -Wundef -Wshadow '
'-Wpointer-arith -Wformat -Wreturn-type -Wsign-compare -Wmultichar '
'-Wformat-nonliteral -Winit-self -Wuninitialized -Wformat-security '
'-Werror -Wno-sequence-point ')
SharedLibrary('ecdsa', srcs, CCFLAGS=flags)
'''
open('SConstruct', 'w').write(scons_file)
subprocess.check_call('scons -s', shell=True)
lib = c.cdll.LoadLibrary('./libecdsa.so')
lib.get_curve_by_name.restype = c.c_void_p
BIGNUM = c.c_uint32 * 9
class Random(random.Random):
def randbytes(self, n):
buf = (c.c_uint8 * n)()
for i in range(n):
buf[i] = self.randrange(0, 256)
return buf
def randpoint(self, curve):
k = self.randrange(0, curve.order)
return k * curve.generator
def int2bn(x, bn_type=BIGNUM):
b = bn_type()
b._int = x
for i in range(len(b)):
b[i] = x % (1 << 30)
x = x >> 30
return b
def bn2int(b):
x = 0
for i in range(len(b)):
x += (b[i] << (30 * i))
return x
@pytest.fixture(params=range(random_iters))
def r(request):
seed = request.param
return Random(seed + int(os.environ.get('SEED', 0)))
@pytest.fixture(params=list(sorted(curves)))
def curve(request):
name = request.param
curve_ptr = lib.get_curve_by_name(name)
assert curve_ptr, 'curve {} not found'.format(name)
curve_obj = curves[name]
curve_obj.ptr = c.c_void_p(curve_ptr)
curve_obj.p = curve_obj.curve.p() # shorthand
return curve_obj
def test_inverse(curve, r):
x = r.randrange(1, curve.p)
y = int2bn(x)
lib.bn_inverse(y, int2bn(curve.p))
y = bn2int(y)
y_ = ecdsa.numbertheory.inverse_mod(x, curve.p)
assert y == y_
def test_inverse(curve, r):
x = r.randrange(0, 2*curve.p)
y = int2bn(x)
lib.bn_mult_half(y, int2bn(curve.p))
y = bn2int(y)
if y > curve.p:
y -= curve.p
half = ecdsa.numbertheory.inverse_mod(2, curve.p)
assert y == (x * half) % curve.p
def test_subtractmod(curve, r):
x = r.randrange(0, 2 ** 256)
y = r.randrange(0, 2 ** 256)
z = int2bn(0)
lib.bn_subtractmod(int2bn(x), int2bn(y), z, int2bn(curve.p))
z = bn2int(z)
z_ = x + 2*curve.p - y
assert z == z_
def test_subtract2(r):
x = r.randrange(0, 2 ** 256)
y = r.randrange(0, 2 ** 256)
x, y = max(x, y), min(x, y)
z = int2bn(0)
lib.bn_subtract(int2bn(x), int2bn(y), z)
z = bn2int(z)
z_ = x - y
assert z == z_
def test_addmod(curve, r):
x = r.randrange(0, 2 ** 256)
y = r.randrange(0, 2 ** 256)
z_ = (x + y) % curve.p
z = int2bn(x)
lib.bn_addmod(z, int2bn(y), int2bn(curve.p))
z = bn2int(z)
assert z == z_
def test_multiply(curve, r):
k = r.randrange(0, 2 * curve.p)
x = r.randrange(0, 2 * curve.p)
z = (k * x) % curve.p
k = int2bn(k)
z_ = int2bn(x)
p_ = int2bn(curve.p)
lib.bn_multiply(k, z_, p_)
z_ = bn2int(z_)
assert z_ < 2*curve.p
if z_ >= curve.p:
z_ = z_ - curve.p
assert z_ == z
def test_multiply1(curve, r):
k = r.randrange(0, 2 * curve.p)
x = r.randrange(0, 2 * curve.p)
kx = k * x
res = int2bn(0, bn_type=(c.c_uint32 * 18))
lib.bn_multiply_long(int2bn(k), int2bn(x), res)
res = bn2int(res)
assert res == kx
def test_multiply2(curve, r):
x = int2bn(0)
s = r.randrange(0, 2 ** 526)
res = int2bn(s, bn_type=(c.c_uint32 * 18))
prime = int2bn(curve.p)
lib.bn_multiply_reduce(x, res, prime)
x = bn2int(x)
x_ = s % curve.p
assert x == x_
def test_fast_mod(curve, r):
x = r.randrange(0, 128*curve.p)
y = int2bn(x)
lib.bn_fast_mod(y, int2bn(curve.p))
y = bn2int(y)
assert y < 2*curve.p
if y >= curve.p:
y -= curve.p
assert x % curve.p == y
def test_mod(curve, r):
x = r.randrange(0, 2*curve.p)
y = int2bn(x)
lib.bn_mod(y, int2bn(curve.p))
assert bn2int(y) == x % curve.p
POINT = BIGNUM * 2
to_POINT = lambda p: POINT(int2bn(p.x()), int2bn(p.y()))
from_POINT = lambda p: (bn2int(p[0]), bn2int(p[1]))
JACOBIAN = BIGNUM * 3
to_JACOBIAN = lambda jp: JACOBIAN(int2bn(jp[0]), int2bn(jp[1]), int2bn(jp[2]))
from_JACOBIAN = lambda p: (bn2int(p[0]), bn2int(p[1]), bn2int(p[2]))
def test_point_multiply(curve, r):
p = r.randpoint(curve)
k = r.randrange(0, 2 ** 256)
kp = k * p
res = POINT(int2bn(0), int2bn(0))
lib.point_multiply(curve.ptr, int2bn(k), to_POINT(p), res)
res = from_POINT(res)
assert res == (kp.x(), kp.y())
def test_point_add(curve, r):
p1 = r.randpoint(curve)
p2 = r.randpoint(curve)
#print '-' * 80
q = p1 + p2
q1 = to_POINT(p1)
q2 = to_POINT(p2)
lib.point_add(curve.ptr, q1, q2)
q_ = from_POINT(q2)
assert q_ == (q.x(), q.y())
def test_point_double(curve, r):
p = r.randpoint(curve)
q = p.double()
q_ = to_POINT(p)
lib.point_double(curve.ptr, q_)
q_ = from_POINT(q_)
assert q_ == (q.x(), q.y())
def test_point_to_jacobian(curve, r):
p = r.randpoint(curve)
jp = JACOBIAN()
lib.curve_to_jacobian(to_POINT(p), jp, int2bn(curve.p))
jx, jy, jz = from_JACOBIAN(jp)
assert jx == (p.x() * jz ** 2) % curve.p
assert jy == (p.y() * jz ** 3) % curve.p
q = POINT()
lib.jacobian_to_curve(jp, q, int2bn(curve.p))
q = from_POINT(q)
assert q == (p.x(), p.y())
def test_cond_negate(curve, r):
x = r.randrange(0, curve.p)
a = int2bn(x)
lib.conditional_negate(0, a, int2bn(curve.p))
assert bn2int(a) == x
lib.conditional_negate(-1, a, int2bn(curve.p))
assert bn2int(a) == curve.p - x
def test_jacobian_add(curve, r):
p1 = r.randpoint(curve)
p2 = r.randpoint(curve)
prime = int2bn(curve.p)
q = POINT()
jp2 = JACOBIAN()
lib.curve_to_jacobian(to_POINT(p2), jp2, prime)
lib.point_jacobian_add(to_POINT(p1), jp2, prime)
lib.jacobian_to_curve(jp2, q, prime)
q = from_POINT(q)
p_ = p1 + p2
assert (p_.x(), p_.y()) == q
def test_jacobian_double(curve, r):
p = r.randpoint(curve)
p2 = p.double()
prime = int2bn(curve.p)
q = POINT()
jp = JACOBIAN()
lib.curve_to_jacobian(to_POINT(p), jp, prime)
lib.point_jacobian_double(jp, curve.ptr)
lib.jacobian_to_curve(jp, q, prime)
q = from_POINT(q)
assert (p2.x(), p2.y()) == q
def sigdecode(sig, _):
return map(bytes2num, [sig[:32], sig[32:]])
def test_sign(curve, r):
priv = r.randbytes(32)
digest = r.randbytes(32)
sig = r.randbytes(64)
lib.ecdsa_sign_digest(curve.ptr, priv, digest, sig, c.c_void_p(0))
exp = bytes2num(priv)
sk = ecdsa.SigningKey.from_secret_exponent(exp, curve,
hashfunc=hashlib.sha256)
vk = sk.get_verifying_key()
sig_ref = sk.sign_digest_deterministic(digest, hashfunc=hashlib.sha256)
assert binascii.hexlify(sig) == binascii.hexlify(sig_ref)
assert vk.verify_digest(sig, digest, sigdecode)

@ -37,6 +37,13 @@
#include "rand.h"
#include "sha2.h"
#include "options.h"
#include "secp256k1.h"
#define CURVE (&secp256k1)
#define prime256k1 (secp256k1.prime)
#define G256k1 (secp256k1.G)
#define order256k1 (secp256k1.order)
#define secp256k1_cp (secp256k1.cp)
uint8_t *fromhex(const char *str)
{
@ -502,7 +509,7 @@ END_TEST
#define test_deterministic(KEY, MSG, K) do { \
sha256_Raw((uint8_t *)MSG, strlen(MSG), buf); \
res = generate_k_rfc6979(&k, fromhex(KEY), buf); \
res = generate_k_rfc6979(CURVE, &k, fromhex(KEY), buf); \
ck_assert_int_eq(res, 0); \
bn_write_be(&k, buf); \
ck_assert_mem_eq(buf, fromhex(K), 32); \
@ -537,13 +544,13 @@ START_TEST(test_sign_speed)
memcpy(priv_key, fromhex("c55ece858b0ddd5263f96810fe14437cd3b5e1fbd7c6a2ec1e031f05e86d8bd5"), 32);
for (i = 0 ; i < 250; i++) {
res = ecdsa_sign(priv_key, msg, sizeof(msg), sig, 0);
res = ecdsa_sign(CURVE, priv_key, msg, sizeof(msg), sig, 0);
ck_assert_int_eq(res, 0);
}
memcpy(priv_key, fromhex("509a0382ff5da48e402967a671bdcde70046d07f0df52cff12e8e3883b426a0a"), 32);
for (i = 0 ; i < 250; i++) {
res = ecdsa_sign(priv_key, msg, sizeof(msg), sig, 0);
res = ecdsa_sign(CURVE, priv_key, msg, sizeof(msg), sig, 0);
ck_assert_int_eq(res, 0);
}
@ -568,9 +575,9 @@ START_TEST(test_verify_speed)
memcpy(pub_key65, fromhex("044054fd18aeb277aeedea01d3f3986ff4e5be18092a04339dcf4e524e2c0a09746c7083ed2097011b1223a17a644e81f59aa3de22dac119fd980b36a8ff29a244"), 65);
for (i = 0 ; i < 25; i++) {
res = ecdsa_verify(pub_key65, sig, msg, sizeof(msg));
res = ecdsa_verify(CURVE, pub_key65, sig, msg, sizeof(msg));
ck_assert_int_eq(res, 0);
res = ecdsa_verify(pub_key33, sig, msg, sizeof(msg));
res = ecdsa_verify(CURVE, pub_key33, sig, msg, sizeof(msg));
ck_assert_int_eq(res, 0);
}
@ -579,9 +586,9 @@ START_TEST(test_verify_speed)
memcpy(pub_key65, fromhex("04ff45a5561a76be930358457d113f25fac790794ec70317eff3b97d7080d457196235193a15778062ddaa44aef7e6901b781763e52147f2504e268b2d572bf197"), 65);
for (i = 0 ; i < 25; i++) {
res = ecdsa_verify(pub_key65, sig, msg, sizeof(msg));
res = ecdsa_verify(CURVE, pub_key65, sig, msg, sizeof(msg));
ck_assert_int_eq(res, 0);
res = ecdsa_verify(pub_key33, sig, msg, sizeof(msg));
res = ecdsa_verify(CURVE, pub_key33, sig, msg, sizeof(msg));
ck_assert_int_eq(res, 0);
}
@ -1034,43 +1041,43 @@ START_TEST(test_pubkey_validity)
int res;
memcpy(pub_key, fromhex("0226659c1cf7321c178c07437150639ff0c5b7679c7ea195253ed9abda2e081a37"), 33);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("025b1654a0e78d28810094f6c5a96b8efb8a65668b578f170ac2b1f83bc63ba856"), 33);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("03433f246a12e6486a51ff08802228c61cf895175a9b49ed4766ea9a9294a3c7fe"), 33);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("03aeb03abeee0f0f8b4f7a5d65ce31f9570cef9f72c2dd8a19b4085a30ab033d48"), 33);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("0496e8f2093f018aff6c2e2da5201ee528e2c8accbf9cac51563d33a7bb74a016054201c025e2a5d96b1629b95194e806c63eb96facaedc733b1a4b70ab3b33e3a"), 65);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("0498010f8a687439ff497d3074beb4519754e72c4b6220fb669224749591dde416f3961f8ece18f8689bb32235e436874d2174048b86118a00afbd5a4f33a24f0f"), 65);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("04f80490839af36d13701ec3f9eebdac901b51c362119d74553a3c537faff31b17e2a59ebddbdac9e87b816307a7ed5b826b8f40b92719086238e1bebf19b77a4d"), 65);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 1);
memcpy(pub_key, fromhex("04f80490839af36d13701ec3f9eebdac901b51c362119d74553a3c537faff31b17e2a59ebddbdac9e87b816307a7ed5b826b8f40b92719086238e1bebf00000000"), 65);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 0);
memcpy(pub_key, fromhex("04f80490839af36d13701ec3f9eebdac901b51c362119d74553a3c537faff31b17e2a59ebddbdac9e87b816307a7ed5b8211111111111111111111111111111111"), 65);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 0);
memcpy(pub_key, fromhex("00"), 1);
res = ecdsa_read_pubkey(pub_key, &pub);
res = ecdsa_read_pubkey(CURVE, pub_key, &pub);
ck_assert_int_eq(res, 0);
}
END_TEST
@ -1214,21 +1221,21 @@ START_TEST(test_secp256k1_cp) {
bn_normalize(&a);
// note that this is not a trivial test. We add 64 curve
// points in the table to get that particular curve point.
scalar_multiply(&a, &p);
scalar_multiply(CURVE, &a, &p);
ck_assert_mem_eq(&p, &secp256k1_cp[i][j], sizeof(curve_point));
bn_zero(&p.y); // test that point_multiply is not a noop
point_multiply(&a, &G256k1, &p);
bn_zero(&p.y); // test that point_multiply CURVE, is not a noop
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert_mem_eq(&p, &secp256k1_cp[i][j], sizeof(curve_point));
// even/odd has different behaviour;
// increment by one and test again
p1 = p;
point_add(&G256k1, &p1);
point_add(CURVE, &G256k1, &p1);
bn_addmodi(&a, 1, &order256k1);
scalar_multiply(&a, &p);
scalar_multiply(CURVE, &a, &p);
ck_assert_mem_eq(&p, &p1, sizeof(curve_point));
bn_zero(&p.y); // test that point_multiply is not a noop
point_multiply(&a, &G256k1, &p);
bn_zero(&p.y); // test that point_multiply CURVE, is not a noop
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert_mem_eq(&p, &p1, sizeof(curve_point));
}
}
@ -1240,43 +1247,43 @@ START_TEST(test_mult_border_cases) {
curve_point p;
curve_point expected;
bn_zero(&a); // a == 0
scalar_multiply(&a, &p);
scalar_multiply(CURVE, &a, &p);
ck_assert(point_is_infinity(&p));
point_multiply(&a, &p, &p);
point_multiply(CURVE, &a, &p, &p);
ck_assert(point_is_infinity(&p));
point_multiply(&a, &G256k1, &p);
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert(point_is_infinity(&p));
bn_addmodi(&a, 1, &order256k1); // a == 1
scalar_multiply(&a, &p);
scalar_multiply(CURVE, &a, &p);
ck_assert_mem_eq(&p, &G256k1, sizeof(curve_point));
point_multiply(&a, &G256k1, &p);
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert_mem_eq(&p, &G256k1, sizeof(curve_point));
bn_subtract(&order256k1, &a, &a); // a == -1
expected = G256k1;
bn_subtract(&prime256k1, &expected.y, &expected.y);
scalar_multiply(&a, &p);
scalar_multiply(CURVE, &a, &p);
ck_assert_mem_eq(&p, &expected, sizeof(curve_point));
point_multiply(&a, &G256k1, &p);
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert_mem_eq(&p, &expected, sizeof(curve_point));
bn_subtract(&order256k1, &a, &a);
bn_addmodi(&a, 1, &order256k1); // a == 2
expected = G256k1;
point_add(&expected, &expected);
scalar_multiply(&a, &p);
point_add(CURVE, &expected, &expected);
scalar_multiply(CURVE, &a, &p);
ck_assert_mem_eq(&p, &expected, sizeof(curve_point));
point_multiply(&a, &G256k1, &p);
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert_mem_eq(&p, &expected, sizeof(curve_point));
bn_subtract(&order256k1, &a, &a); // a == -2
expected = G256k1;
point_add(&expected, &expected);
point_add(CURVE, &expected, &expected);
bn_subtract(&prime256k1, &expected.y, &expected.y);
scalar_multiply(&a, &p);
scalar_multiply(CURVE, &a, &p);
ck_assert_mem_eq(&p, &expected, sizeof(curve_point));
point_multiply(&a, &G256k1, &p);
point_multiply(CURVE, &a, &G256k1, &p);
ck_assert_mem_eq(&p, &expected, sizeof(curve_point));
}
END_TEST
@ -1289,11 +1296,11 @@ START_TEST(test_scalar_mult) {
curve_point p1, p2, p3;
for (i = 0; i < 1000; i++) {
/* test distributivity: (a + b)G = aG + bG */
scalar_multiply(&a, &p1);
scalar_multiply(&b, &p2);
scalar_multiply(CURVE, &a, &p1);
scalar_multiply(CURVE, &b, &p2);
bn_addmod(&a, &b, &order256k1);
scalar_multiply(&a, &p3);
point_add(&p1, &p2);
scalar_multiply(CURVE, &a, &p3);
point_add(CURVE, &p1, &p2);
ck_assert_mem_eq(&p2, &p3, sizeof(curve_point));
// new "random" numbers
a = p3.x;
@ -1311,11 +1318,11 @@ START_TEST(test_point_mult) {
curve_point p1, p2, p3;
for (i = 0; i < 200; i++) {
/* test distributivity: (a + b)P = aP + bP */
point_multiply(&a, &p, &p1);
point_multiply(&b, &p, &p2);
point_multiply(CURVE, &a, &p, &p1);
point_multiply(CURVE, &b, &p, &p2);
bn_addmod(&a, &b, &order256k1);
point_multiply(&a, &p, &p3);
point_add(&p1, &p2);
point_multiply(CURVE, &a, &p, &p3);
point_add(CURVE, &p1, &p2);
ck_assert_mem_eq(&p2, &p3, sizeof(curve_point));
// new "random" numbers and a "random" point
a = p1.x;
@ -1335,16 +1342,16 @@ START_TEST(test_scalar_point_mult) {
/* test commutativity and associativity:
* a(bG) = (ab)G = b(aG)
*/
scalar_multiply(&a, &p1);
point_multiply(&b, &p1, &p1);
scalar_multiply(CURVE, &a, &p1);
point_multiply(CURVE, &b, &p1, &p1);
scalar_multiply(&b, &p2);
point_multiply(&a, &p2, &p2);
scalar_multiply(CURVE, &b, &p2);
point_multiply(CURVE, &a, &p2, &p2);
ck_assert_mem_eq(&p1, &p2, sizeof(curve_point));
bn_multiply(&a, &b, &order256k1);
scalar_multiply(&b, &p2);
scalar_multiply(CURVE, &b, &p2);
ck_assert_mem_eq(&p1, &p2, sizeof(curve_point));

2
tools/.gitignore vendored

@ -1,3 +1,3 @@
xpubaddrgen
mksecptable
mktable
bip39bruteforce

@ -23,9 +23,9 @@ CFLAGS += $(OPTFLAGS) \
-Werror \
-I..
all: xpubaddrgen mksecptable bip39bruteforce
all: xpubaddrgen mktable bip39bruteforce
OBJS = ../bip32.o ../bip39.o ../ecdsa.o ../sha2.o ../bignum.o ../base58.o ../secp256k1.o ../ripemd160.o ../hmac.o ../rand.o ../pbkdf2.o
OBJS = ../bip32.o ../bip39.o ../ecdsa.o ../sha2.o ../bignum.o ../base58.o ../secp256k1.o ../nist256p1.o ../ripemd160.o ../hmac.o ../rand.o ../pbkdf2.o
%.o: %.c %.h options.h
$(CC) $(CFLAGS) -o $@ -c $<
@ -33,11 +33,11 @@ OBJS = ../bip32.o ../bip39.o ../ecdsa.o ../sha2.o ../bignum.o ../base58.o ../sec
xpubaddrgen: xpubaddrgen.o $(OBJS)
$(CC) xpubaddrgen.o $(OBJS) -o xpubaddrgen
mksecptable: mksecptable.o $(OBJS)
$(CC) mksecptable.o $(OBJS) -o mksecptable
mktable: mktable.o $(OBJS)
$(CC) mktable.o $(OBJS) -o mktable
bip39bruteforce: bip39bruteforce.o $(OBJS)
$(CC) bip39bruteforce.o $(OBJS) -o bip39bruteforce
clean:
rm -f *.o xpubaddrgen mksecptable bip39bruteforce
rm -f *.o xpubaddrgen mktable bip39bruteforce

@ -45,9 +45,10 @@ It will print ```<jobid> error``` when there was an error processing job jobid.
It will print ```error``` when it encountered a malformed line.
mksecptable
mktable
-----------
mksecptable computes the points of the form `(2*j+1)*16^i*G` and prints them in the format to be included in `secp256k1.c`. These points are used by the fast ECC multiplication.
mktable computes the points of the form `(2*j+1)*16^i*G` and prints them in the format to be included in `secp256k1.c` and `nist256p1.c`.
These points are used by the fast ECC multiplication.
It is only meant to be run if the `scalar_mult` algorithm changes.

@ -10,10 +10,21 @@
* The entry secp256k1_cp[i][j] contains the number (2*j+1)*16^i*G,
* where G is the generator of secp256k1.
*/
int main(int __attribute__((unused)) argc, char __attribute__((unused)) **argv) {
int main(int argc, char **argv) {
int i,j,k;
curve_point ng = G256k1;
curve_point pow2ig = G256k1;
if (argc != 2) {
printf("Usage: %s CURVE_NAME\n", argv[0]);
return 1;
}
const char *name = argv[1];
const ecdsa_curve *curve = get_curve_by_name(name);
if (curve == 0) {
printf("Unknown curve '%s'\n", name);
return 1;
}
curve_point ng = curve->G;
curve_point pow2ig = curve->G;
for (i = 0; i < 64; i++) {
// invariants:
// pow2ig = 16^i * G
@ -29,7 +40,7 @@ int main(int __attribute__((unused)) argc, char __attribute__((unused)) **argv)
bn_zero(&a);
a.val[(4*i) / 30] = ((uint32_t) 2*j+1) << ((4*i) % 30);
bn_normalize(&a);
point_multiply(&a, &G256k1, &checkresult);
point_multiply(curve, &a, &curve->G, &checkresult);
assert(point_is_equal(&checkresult, &ng));
#endif
printf("\t\t/* %2d*16^%d*G: */\n\t\t{{{", 2*j + 1, i);
@ -46,9 +57,9 @@ int main(int __attribute__((unused)) argc, char __attribute__((unused)) **argv)
printf("}}}\n\t},\n");
} else {
printf("}}},\n");
point_add(&pow2ig, &ng);
point_add(curve, &pow2ig, &ng);
}
point_add(&pow2ig, &ng);
point_add(curve, &pow2ig, &ng);
}
pow2ig = ng;
}
Loading…
Cancel
Save