pub struct MatZq { /* private fields */ }Expand description
MatZq is a matrix with entries of type Zq.
Attributes:
§Examples
§Matrix usage
use qfall_math::{
integer::Z,
integer_mod_q::MatZq,
traits::{MatrixGetEntry, MatrixSetEntry},
};
use std::str::FromStr;
// instantiate new matrix
let id_mat = MatZq::from_str("[[1, 0],[0, 1]] mod 2").unwrap();
// clone object, set and get entry
let mut clone = id_mat.clone();
clone.set_entry(0, 0, 2);
let entry:Z = clone.get_entry(1, 1).unwrap();
assert_eq!(entry, Z::ONE);
// to_string incl. (de-)serialization
assert_eq!("[[1, 0],[0, 1]] mod 2", &id_mat.to_string());§Vector usage
use qfall_math::{
integer::Z,
integer_mod_q::MatZq,
};
use std::str::FromStr;
let row_vec = MatZq::from_str("[[1, 1, 1]] mod 2").unwrap();
let col_vec = MatZq::from_str("[[1],[-1],[0]] mod 2").unwrap();
// check if matrix instance is vector
assert!(row_vec.is_row_vector());
assert!(col_vec.is_column_vector());Implementations§
Source§impl MatZq
impl MatZq
Sourcepub fn add_safe(&self, other: &Self) -> Result<MatZq, MathError>
pub fn add_safe(&self, other: &Self) -> Result<MatZq, MathError>
Implements addition for two MatZq matrices.
Parameters:
other: specifies the value to add toself
Returns the sum of both matrices as a MatZq or an
error if the matrix dimensions mismatch.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a: MatZq = MatZq::from_str("[[1, 2, 3],[3, 4, 5]] mod 7").unwrap();
let b: MatZq = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c: MatZq = a.add_safe(&b).unwrap();§Errors and Failures
- Returns a
MathErrorof typeMathError::MismatchingMatrixDimensionif the matrix dimensions mismatch. - Returns a
MathErrorof typeMathError::MismatchingModulusif the moduli mismatch.
Source§impl MatZq
impl MatZq
Sourcepub fn mul_safe(&self, other: &Self) -> Result<Self, MathError>
pub fn mul_safe(&self, other: &Self) -> Result<Self, MathError>
Implements multiplication for two MatZq values.
Parameters:
other: specifies the value to multiply withself
Returns the product of self and other as a MatZq or an error if the dimensions of self
and other do not match for multiplication or the moduli mismatch.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a = MatZq::from_str("[[2, 1],[1, 2]] mod 7").unwrap();
let b = MatZq::from_str("[[1, 0],[0, 1]] mod 7").unwrap();
let c: MatZq = a.mul_safe(&b).unwrap();§Errors and Failures
- Returns a
MathErrorof typeMathError::MismatchingMatrixDimensionif the dimensions ofselfandotherdo not match for multiplication. - Returns a
MathErrorof typeMathError::MismatchingModulusif the moduli mismatch.
Source§impl MatZq
impl MatZq
Sourcepub fn mul_scalar_safe(&self, scalar: &Zq) -> Result<Self, MathError>
pub fn mul_scalar_safe(&self, scalar: &Zq) -> Result<Self, MathError>
Implements multiplication for a MatZq matrix with a Zq.
Parameters:
scalar: specifies the scalar by which the matrix is multiplied
Returns the product of self and scalar as a MatZq or
an error if the moduli mismatch.
§Examples
use qfall_math::integer_mod_q::{MatZq, Zq};
use std::str::FromStr;
let mat_1 = MatZq::from_str("[[42, 17],[8, 6]] mod 61").unwrap();
let integer = Zq::from((2, 61));
let mat_2 = &mat_1.mul_scalar_safe(&integer).unwrap();§Errors and Failures
- Returns a
MathErrorof typeMismatchingModulusif the moduli mismatch.
Source§impl MatZq
impl MatZq
Sourcepub fn sub_safe(&self, other: &Self) -> Result<MatZq, MathError>
pub fn sub_safe(&self, other: &Self) -> Result<MatZq, MathError>
Implements subtraction for two MatZq matrices.
Parameters:
other: specifies the value to subtract fromself
Returns the result of the subtraction as a MatZq or an
error if the matrix dimensions or moduli mismatch.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a: MatZq = MatZq::from_str("[[1, 2, 3],[3, 4, 5]] mod 7").unwrap();
let b: MatZq = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c: MatZq = a.sub_safe(&b).unwrap();§Errors and Failures
- Returns a
MathErrorof typeMathError::MismatchingMatrixDimensionif the matrix dimensions mismatch. - Returns a
MathErrorof typeMathError::MismatchingModulusif the moduli mismatch.
Source§impl MatZq
impl MatZq
Sourcepub fn new(
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
modulus: impl Into<Modulus>,
) -> Self
pub fn new( num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, modulus: impl Into<Modulus>, ) -> Self
Creates a new matrix with num_rows rows, num_cols columns,
zeros as entries and modulus as the modulus.
Parameters:
num_rows: number of rows the new matrix should havenum_cols: number of columns the new matrix should havemodulus: the common modulus of the matrix entries
Returns a new MatZq instance of the provided dimensions.
§Examples
use qfall_math::integer_mod_q::MatZq;
let matrix = MatZq::new(5, 10, 7);§Panics …
- if the number of rows or columns is negative,
0, or does not fit into ani64. - if
modulusis smaller than2.
Sourcepub fn identity(
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
modulus: impl Into<Modulus>,
) -> Self
pub fn identity( num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, modulus: impl Into<Modulus>, ) -> Self
Generate a num_rows times num_columns matrix with 1 on the
diagonal and 0 anywhere else with a given modulus.
Parameters:
rum_rows: the number of rows of the identity matrixnum_columns: the number of columns of the identity matrixmodulus: the modulus of the matrix
Returns a matrix with 1 across the diagonal and 0 anywhere else.
§Examples
use qfall_math::integer_mod_q::MatZq;
let matrix = MatZq::identity(2, 3, 3);
let identity = MatZq::identity(10, 10, 3);§Panics …
- if the provided number of rows and columns or the modulus are not suited to create a matrix.
For further information see
MatZq::new.
Source§impl MatZq
impl MatZq
Sourcepub fn from_utf8(
message: &str,
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
modulus: impl Into<Modulus>,
) -> Result<Self, MathError>
pub fn from_utf8( message: &str, num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, modulus: impl Into<Modulus>, ) -> Result<Self, MathError>
Create a MatZq from a String, i.e. its UTF8-Encoding.
This function can only construct positive or zero integers, but not negative ones.
If the number of bytes and number of entries does not line up, we pad the message
with '0's.
The inverse of this function is MatZq::to_utf8.
WARNING: This implementation requires the modulus to be larger than
any single entry in the matrix. This function will denote the same number of bytes
to every entry and sequentially move through your message to encode them.
If a decimal presentation of these bytes is ever larger than the specified modulus,
the function will return an error.
Parameters:
message: specifies the message that is transformed via its UTF8-Encoding to a newMatZqinstance.num_rows: number of rows the new matrix should havenum_cols: number of columns the new matrix should havemodulus: specifies the modulus of the matrix, it is required to be larger than any entry of the matrix
Returns a MatZq with corresponding entries to the message’s UTF8-Encoding or
a ConversionError if the modulus isn’t larger than
every single entry of the matrix after distributing the (potentially padded) UTF8-Bytes
equally over the matrix.
§Examples
use qfall_math::integer_mod_q::MatZq;
let message = "hello!";
let matrix = MatZq::from_utf8(&message, 3, 2, 257).unwrap();§Errors and Failures
- Returns a
MathErrorof typeConversionErrorif the modulus isn’t larger than the largest entry of the matrix after equally distributing the (potentially padded) UTF8-Conversion over the matrix.
§Panics …
- if the provided number of rows and columns are not suited to create a matrix.
For further information see
MatZq::new.
Source§impl MatZq
impl MatZq
Sourcepub fn get_representative_least_nonnegative_residue(&self) -> MatZ
pub fn get_representative_least_nonnegative_residue(&self) -> MatZ
Creates a MatZ where each entry is a representative of the
equivalence class of each entry from a MatZq.
The values in the output matrix are in the range of [0, modulus).
Use MatZq::get_representative_least_absolute_residue if they should be
in the range [-modulus/2, modulus/2].
Returns the matrix as a MatZ.
§Examples
use qfall_math::integer::MatZ;
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat_zq = MatZq::from_str("[[1, 2],[3, -1]] mod 5").unwrap();
let mat_z = mat_zq.get_representative_least_nonnegative_residue();
assert_eq!(mat_z.to_string(), "[[1, 2],[3, 4]]");Sourcepub fn get_representative_least_absolute_residue(&self) -> MatZ
pub fn get_representative_least_absolute_residue(&self) -> MatZ
Creates a MatZ where each entry is a representative of the
equivalence class of each entry from a MatZq with the
representatives close to 0.
The values in the output matrix are in the range of [-modulus/2, modulus/2].
For even moduli, the positive representative is chosen for the element modulus / 2.
Use MatZq::get_representative_least_nonnegative_residue if they should be
in the range [0, modulus).
Returns an MatZ representation of the given matrix with
representatives chosen close to 0.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat_zq_1 = MatZq::from_str("[[1,2],[3,4]] mod 5").unwrap();
let mat_zq_2 = MatZq::from_str("[[1,2],[3,4]] mod 4").unwrap();
let mat_z_1 = mat_zq_1.get_representative_least_absolute_residue();
let mat_z_2 = mat_zq_2.get_representative_least_absolute_residue();
assert_eq!(mat_z_1.to_string(), "[[1, 2],[-2, -1]]");
assert_eq!(mat_z_2.to_string(), "[[1, 2],[-1, 0]]");Source§impl MatZq
impl MatZq
Sourcepub fn inverse(&self) -> Option<MatZq>
pub fn inverse(&self) -> Option<MatZq>
Returns the inverse of the matrix if it exists (is square and
has a determinant co-prime to the modulus) and None otherwise.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mut matrix = MatZq::from_str("[[1, 2],[3, 4]] mod 7").unwrap();
let matrix_invert = matrix.inverse().unwrap();
let id = &matrix_invert * &matrix;
assert_eq!("[[5, 1],[5, 3]] mod 7", matrix_invert.to_string());
assert!(id.is_identity());Sourcepub fn inverse_prime(&self) -> Option<MatZq>
pub fn inverse_prime(&self) -> Option<MatZq>
Returns the inverse of the matrix if it exists (is square and
has a determinant co-prime to the modulus) and None otherwise.
Note that the modulus is assumed to be prime, otherwise the function panics.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mut matrix = MatZq::from_str("[[1, 2],[3, 4]] mod 7").unwrap();
let matrix_invert = matrix.inverse_prime().unwrap();
let id = &matrix_invert * &matrix;
assert_eq!("[[5, 1],[5, 3]] mod 7", matrix_invert.to_string());
assert!(id.is_identity());§Panics …
- if the modulus is not prime.
Sourcepub fn gaussian_elimination_prime(self) -> MatZq
pub fn gaussian_elimination_prime(self) -> MatZq
Returns the row echelon form of the matrix using gaussian elimination.
Note that the modulus is assumed to be prime, otherwise the function panics.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mut matrix = MatZq::from_str("[[1, 2],[3, 4]] mod 7").unwrap();
let matrix_gauss = matrix.gaussian_elimination_prime();
assert_eq!("[[1, 0],[0, 1]] mod 7", matrix_gauss.to_string());§Panics …
- if the modulus is not prime.
Source§impl MatZq
impl MatZq
Sourcepub fn norm_l_2_infty_sqrd(&self) -> Z
pub fn norm_l_2_infty_sqrd(&self) -> Z
Outputs the squared l_{2, ∞}-norm, i.e. it computes the squared Euclidean norm of each column of the matrix and returns the largest one.
§Examples
use qfall_math::{integer_mod_q::MatZq, integer::Z};
use std::str::FromStr;
let mat = MatZq::from_str("[[2, 3],[2, 0]] mod 7").unwrap();
let eucl_norm = mat.norm_l_2_infty_sqrd();
// 3^2 + 0^2 = 9
assert_eq!(Z::from(9), eucl_norm);Sourcepub fn norm_l_2_infty(&self) -> Q
pub fn norm_l_2_infty(&self) -> Q
Outputs the l_{2, ∞}-norm, i.e. it computes the Euclidean norm of each column of the matrix and returns the largest one.
§Examples
use qfall_math::{integer_mod_q::MatZq, rational::Q};
use std::str::FromStr;
let mat = MatZq::from_str("[[2, 3],[2, 0],[3, 4],[3, 4]] mod 5").unwrap();
let eucl_norm = mat.norm_l_2_infty();
// sqrt(4 * 2^2) = 4
assert_eq!(Q::from(4), eucl_norm);Sourcepub fn norm_l_infty_infty(&self) -> Z
pub fn norm_l_infty_infty(&self) -> Z
Outputs the l_{∞, ∞}-norm, i.e. it computes the ∞-norm of each column of the matrix and returns the largest one.
§Examples
use qfall_math::{integer_mod_q::MatZq, integer::Z};
use std::str::FromStr;
let mat = MatZq::from_str("[[2, 4],[2, 0]] mod 7").unwrap();
let eucl_norm = mat.norm_l_infty_infty();
// max{2, 3} = 3
assert_eq!(Z::from(3), eucl_norm);Source§impl MatZq
impl MatZq
Sourcepub fn is_identity(&self) -> bool
pub fn is_identity(&self) -> bool
Checks if a MatZq is the identity matrix.
Returns true if every diagonal entry of the upper square matrix is 1
and all other entries are 0.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let value = MatZq::from_str("[[1, 0],[0, 1]] mod 17").unwrap();
assert!(value.is_identity());use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let value = MatZq::from_str("[[1, 0],[0, 1],[0, 0]] mod 17").unwrap();
assert!(value.is_identity());Sourcepub fn is_symmetric(&self) -> bool
pub fn is_symmetric(&self) -> bool
Source§impl MatZq
impl MatZq
Sourcepub fn sample_binomial(
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
modulus: impl Into<Modulus>,
n: impl Into<Z>,
p: impl Into<Q>,
) -> Result<Self, MathError>
pub fn sample_binomial( num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, modulus: impl Into<Modulus>, n: impl Into<Z>, p: impl Into<Q>, ) -> Result<Self, MathError>
Outputs a MatZq instance with entries chosen according to the binomial
distribution parameterized by n and p.
Parameters:
num_rows: specifies the number of rows the new matrix should havenum_cols: specifies the number of columns the new matrix should havemodulus: specifies theModulusof the newMatZqinstancen: specifies the number of trialsp: specifies the probability of success
Returns a new MatZq instance with entries chosen
according to the binomial distribution or a MathError
if n < 0, p ∉ (0,1), n does not fit into an i64,
or the dimensions of the matrix were chosen too small.
§Examples
use qfall_math::integer_mod_q::MatZq;
let sample = MatZq::sample_binomial(2, 2, 7, 2, 0.5).unwrap();§Errors and Failures
- Returns a
MathErrorof typeInvalidIntegerInputifn < 0. - Returns a
MathErrorof typeInvalidIntervalifp ∉ (0,1). - Returns a
MathErrorof typeConversionErrorifndoes not fit into ani64.
§Panics …
- if the provided number of rows and columns are not suited to create a matrix.
For further information see
MatZq::new. - if
modulusis smaller than2.
Sourcepub fn sample_binomial_with_offset(
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
offset: impl Into<Z>,
modulus: impl Into<Modulus>,
n: impl Into<Z>,
p: impl Into<Q>,
) -> Result<Self, MathError>
pub fn sample_binomial_with_offset( num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, offset: impl Into<Z>, modulus: impl Into<Modulus>, n: impl Into<Z>, p: impl Into<Q>, ) -> Result<Self, MathError>
Outputs a MatZq instance with entries chosen according to the binomial
distribution parameterized by n and p with given offset.
Parameters:
num_rows: specifies the number of rows the new matrix should havenum_cols: specifies the number of columns the new matrix should haveoffset: specifies an offset applied to each sample collected from the binomial distributionmodulus: specifies theModulusof the newMatZqinstancen: specifies the number of trialsp: specifies the probability of success
Returns a new MatZq instance with entries chosen
according to the binomial distribution or a MathError
if n < 0, p ∉ (0,1), n does not fit into an i64,
or the dimensions of the matrix were chosen too small.
§Examples
use qfall_math::integer_mod_q::MatZq;
let sample = MatZq::sample_binomial_with_offset(2, 2, -1, 7, 2, 0.5).unwrap();§Errors and Failures
- Returns a
MathErrorof typeInvalidIntegerInputifn < 0. - Returns a
MathErrorof typeInvalidIntervalifp ∉ (0,1). - Returns a
MathErrorof typeConversionErrorifndoes not fit into ani64.
§Panics …
- if the provided number of rows and columns are not suited to create a matrix.
For further information see
MatZq::new. - if
modulusis smaller than2.
Source§impl MatZq
impl MatZq
Sourcepub fn sample_discrete_gauss(
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
modulus: impl Into<Modulus>,
center: impl Into<Q>,
s: impl Into<Q>,
) -> Result<MatZq, MathError>
pub fn sample_discrete_gauss( num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, modulus: impl Into<Modulus>, center: impl Into<Q>, s: impl Into<Q>, ) -> Result<MatZq, MathError>
Initializes a new matrix with dimensions num_rows x num_columns and with each entry
sampled independently according to the discrete Gaussian distribution,
using Z::sample_discrete_gauss.
Parameters:
num_rows: specifies the number of rows the new matrix should havenum_cols: specifies the number of columns the new matrix should havecenter: specifies the positions of the center with peak probabilitys: specifies the Gaussian parameter, which is proportional to the standard deviationsigma * sqrt(2 * pi) = s
Returns a matrix with each entry sampled independently from the
specified discrete Gaussian distribution or an error if s < 0.
§Examples
use qfall_math::integer_mod_q::MatZq;
let sample = MatZq::sample_discrete_gauss(3, 1, 83, 0, 1.25f32).unwrap();§Errors and Failures
- Returns a
MathErrorof typeInvalidIntegerInputifs < 0.
§Panics …
- if the provided number of rows and columns or the modulus are not suited to create a matrix.
For further information see
MatZq::new. - if the provided
modulus < 2.
Sourcepub fn sample_d(
basis: &MatZq,
center: &MatQ,
s: impl Into<Q>,
) -> Result<Self, MathError>
pub fn sample_d( basis: &MatZq, center: &MatQ, s: impl Into<Q>, ) -> Result<Self, MathError>
SampleD samples a discrete Gaussian from the lattice with a provided basis.
We do not check whether basis is actually a basis. Hence, the callee is
responsible for making sure that basis provides a suitable basis.
Parameters:
basis: specifies a basis for the lattice from which is sampledcenter: specifies the positions of the center with peak probabilitys: specifies the Gaussian parameter, which is proportional to the standard deviationsigma * sqrt(2 * pi) = s
Returns a lattice vector sampled according to the discrete Gaussian distribution
or an error if s < 0, the number of rows of the basis and center differ,
or if center is not a column vector.
§Examples
use qfall_math::{integer_mod_q::MatZq, rational::MatQ};
let basis = MatZq::identity(5, 5, 17);
let center = MatQ::new(5, 1);
let sample = MatZq::sample_d(&basis, ¢er, 1.25f32).unwrap();§Errors and Failures
- Returns a
MathErrorof typeInvalidIntegerInputifs < 0. - Returns a
MathErrorof typeMismatchingMatrixDimensionif the number of rows of thebasisandcenterdiffer. - Returns a
MathErrorof typeStringConversionErrorifcenteris not a column vector.
This function implements SampleD according to:
- [1] Gentry, Craig and Peikert, Chris and Vaikuntanathan, Vinod (2008). Trapdoors for hard lattices and new cryptographic constructions. In: Proceedings of the fortieth annual ACM symposium on Theory of computing. https://dl.acm.org/doi/pdf/10.1145/1374376.1374407
Sourcepub fn sample_d_precomputed_gso(
basis: &MatZq,
basis_gso: &MatQ,
center: &MatQ,
s: impl Into<Q>,
) -> Result<Self, MathError>
pub fn sample_d_precomputed_gso( basis: &MatZq, basis_gso: &MatQ, center: &MatQ, s: impl Into<Q>, ) -> Result<Self, MathError>
SampleD samples a discrete Gaussian from the lattice with a provided basis.
We do not check whether basis is actually a basis or whether basis_gso is
actually the gso of basis. Hence, the callee is responsible for making sure
that basis provides a suitable basis and basis_gso is a corresponding GSO.
Parameters:
basis: specifies a basis for the lattice from which is sampledbasis_gso: specifies the precomputed gso forbasiscenter: specifies the positions of the center with peak probabilitys: specifies the Gaussian parameter, which is proportional to the standard deviationsigma * sqrt(2 * pi) = s
Returns a lattice vector sampled according to the discrete Gaussian distribution
or an error if s < 0, the number of rows of the basis and center differ,
or if center is not a column vector.
§Examples
use qfall_math::{integer::MatZ, integer_mod_q::MatZq, rational::MatQ};
let basis = MatZq::identity(5, 5, 17);
let center = MatQ::new(5, 1);
let basis_gso = MatQ::from(&basis.get_representative_least_nonnegative_residue()).gso();
let sample = MatZq::sample_d_precomputed_gso(&basis, &basis_gso, ¢er, 1.25f32).unwrap();§Errors and Failures
- Returns a
MathErrorof typeInvalidIntegerInputifs < 0. - Returns a
MathErrorof typeMismatchingMatrixDimensionif the number of rows of thebasisandcenterdiffer. - Returns a
MathErrorof typeStringConversionErrorifcenteris not a column vector.
§Panics …
- if the number of rows/columns of
basis_gsoandbasismismatch.
This function implements SampleD according to:
- [1] Gentry, Craig and Peikert, Chris and Vaikuntanathan, Vinod (2008). Trapdoors for hard lattices and new cryptographic constructions. In: Proceedings of the fortieth annual ACM symposium on Theory of computing. https://dl.acm.org/doi/pdf/10.1145/1374376.1374407
Source§impl MatZq
impl MatZq
Sourcepub fn sample_uniform(
num_rows: impl TryInto<i64> + Display,
num_cols: impl TryInto<i64> + Display,
modulus: impl Into<Z>,
) -> Self
pub fn sample_uniform( num_rows: impl TryInto<i64> + Display, num_cols: impl TryInto<i64> + Display, modulus: impl Into<Z>, ) -> Self
Outputs a MatZq instance with entries chosen uniform at random
in [0, modulus).
The internally used uniform at random chosen bytes are generated
by ThreadRng, which uses ChaCha12 and
is considered cryptographically secure.
Parameters:
num_rows: specifies the number of rows the new matrix should havenum_cols: specifies the number of columns the new matrix should havemodulus: specifies the modulus of the matrix and defines the interval over which is sampled
Returns a new MatZq instance with entries chosen
uniformly at random in [0, modulus).
§Examples
use qfall_math::integer_mod_q::MatZq;
let matrix = MatZq::sample_uniform(3, 3, 17);§Panics …
- if the provided number of rows and columns or the modulus are not suited to create a matrix.
For further information see
MatZq::new.
Source§impl MatZq
impl MatZq
Sourcepub fn reverse_columns(&mut self)
pub fn reverse_columns(&mut self)
Swaps the i-th column with the n-i-th column for all i <= n/2
of the specified matrix with n columns.
§Examples
use qfall_math::integer_mod_q::MatZq;
let mut matrix = MatZq::new(4, 3, 5);
matrix.reverse_columns();Sourcepub fn reverse_rows(&mut self)
pub fn reverse_rows(&mut self)
Swaps the i-th row with the n-i-th row for all i <= n/2
of the specified matrix with n rows.
§Examples
use qfall_math::integer_mod_q::MatZq;
let mut matrix = MatZq::new(4, 3, 5);
matrix.reverse_rows();Sourcepub fn change_modulus(&mut self, modulus: impl Into<Modulus>)
pub fn change_modulus(&mut self, modulus: impl Into<Modulus>)
Changes the modulus of the given matrix to the new modulus. It takes the representation of each coefficient in [0, q) as the new matrix entries and reduces them by the new modulus automatically.
Parameters:
modulus: the new modulus of the matrix
§Examples
use qfall_math::integer_mod_q::{MatZq, Modulus};
use std::str::FromStr;
let mut mat = MatZq::from_str("[[1, 2]] mod 3").unwrap();
mat.change_modulus(2);§Panics …
- if
modulusis smaller than2.
Source§impl MatZq
impl MatZq
Sourcepub fn solve_gaussian_elimination(&self, y: &MatZq) -> Option<MatZq>
pub fn solve_gaussian_elimination(&self, y: &MatZq) -> Option<MatZq>
Computes a solution for a system of linear equations under a certain modulus.
It solves Ax = y for x with A being a MatZq value.
If no solution is found, None is returned.
The function uses Gaussian elimination together with Factor refinement
to split the modulus and the Chinese remainder theorem and Hensel lifting
to combine solutions under the split modulus.
For Hensel lifting we use the method from [[1]].
Note that this function does not compute a solution whenever there is one.
If the matrix has not full rank under a modulus that divides the given one,
None may be returned even if the system may be solvable.
If the number of columns exceeds the number of rows by a factor of 2,
this is very unlikely to happen.
Parameters:
y: the syndrome for which a solution has to be computed.
Returns a solution for the linear system or None, if none could be computed.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat = MatZq::from_str("[[2, 2, 3],[2, 5, 7]] mod 8").unwrap();
let y = MatZq::from_str("[[3],[5]] mod 8").unwrap();
let x = mat.solve_gaussian_elimination(&y).unwrap();
assert_eq!(y, mat*x);§Panics …
- if the the number of rows of the matrix and the syndrome are different.
- if the syndrome is not a column vector.
- if the moduli mismatch.
§Reference
- [1] John D. Dixon. “Exact Solution of Linear Equations Using P-Adic Expansions” https://link.springer.com/article/10.1007/BF01459082
Source§impl MatZq
impl MatZq
Sourcepub fn sort_by_column<T: Ord>(
&self,
cond_func: fn(&Self) -> Result<T, MathError>,
) -> Result<Self, MathError>
pub fn sort_by_column<T: Ord>( &self, cond_func: fn(&Self) -> Result<T, MathError>, ) -> Result<Self, MathError>
Sorts the columns of the matrix based on some condition defined by cond_func in an ascending order.
This condition is usually a norm with the described input-output behaviour.
Parameters:
cond_func: computes values implementingOrdover the columns of the specified matrix. These values are then used to re-order / sort the rows of the matrix.
Returns an empty Ok if the action could be performed successfully.
A MathError is returned if the execution of cond_func returned an error.
§Examples
§Use a build-in function as condition
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat = MatZq::from_str("[[3, 2, 1]] mod 7").unwrap();
let cmp = MatZq::from_str("[[1, 2, 3]] mod 7").unwrap();
let sorted = mat.sort_by_column(MatZq::norm_eucl_sqrd).unwrap();
assert_eq!(cmp, sorted);§Use a custom function as condition
This function needs to take a column vector as input and output a type implementing PartialOrd
use qfall_math::{integer_mod_q::MatZq, integer::Z, error::MathError, traits::{MatrixDimensions, MatrixGetEntry}};
use std::str::FromStr;
let mat = MatZq::from_str("[[3, 2, 1]] mod 7").unwrap();
let cmp = MatZq::from_str("[[1, 2, 3]] mod 7").unwrap();
fn custom_cond_func(matrix: &MatZq) -> Result<Z, MathError> {
let mut sum = Z::ZERO;
for entry in MatrixGetEntry::<Z>::get_entries_rowwise(matrix){
sum += entry;
}
Ok(sum)
}
let sorted = mat.sort_by_column(custom_cond_func).unwrap();
assert_eq!(cmp, sorted);§Errors and Failures
- Returns a
MathErrorof the same type ascond_funcif the execution ofcond_funcfails.
Sourcepub fn sort_by_row<T: Ord>(
&self,
cond_func: fn(&Self) -> Result<T, MathError>,
) -> Result<Self, MathError>
pub fn sort_by_row<T: Ord>( &self, cond_func: fn(&Self) -> Result<T, MathError>, ) -> Result<Self, MathError>
Sorts the rows of the matrix based on some condition defined by cond_func in an ascending order.
This condition is usually a norm with the described input-output behaviour.
Parameters:
cond_func: computes values implementingOrdover the columns of the specified matrix. These values are then used to re-order / sort the columns of the matrix.
Returns an empty Ok if the action could be performed successfully.
A MathError is returned if the execution of cond_func returned an error.
§Examples
§Use a build-in function as condition
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat = MatZq::from_str("[[3],[2],[1]] mod 7").unwrap();
let cmp = MatZq::from_str("[[1],[2],[3]] mod 7").unwrap();
let sorted = mat.sort_by_row(MatZq::norm_infty).unwrap();
assert_eq!(cmp, sorted);§Use a custom function as condition
This function needs to take a row vector as input and output a type implementing PartialOrd
use qfall_math::{integer_mod_q::MatZq, integer::Z, error::MathError, traits::{MatrixDimensions, MatrixGetEntry}};
use std::str::FromStr;
let mat = MatZq::from_str("[[3],[2],[1]] mod 7").unwrap();
let cmp = MatZq::from_str("[[1],[2],[3]] mod 7").unwrap();
fn custom_cond_func(matrix: &MatZq) -> Result<Z, MathError> {
let mut sum = Z::ZERO;
for entry in MatrixGetEntry::<Z>::get_entries_rowwise(matrix){
sum += entry;
}
Ok(sum)
}
let sorted = mat.sort_by_row(custom_cond_func).unwrap();
assert_eq!(cmp, sorted);§Errors and Failures
- Returns a
MathErrorof the same type ascond_funcif the execution ofcond_funcfails.
Source§impl MatZq
impl MatZq
Sourcepub fn tensor_product_safe(&self, other: &Self) -> Result<Self, MathError>
pub fn tensor_product_safe(&self, other: &Self) -> Result<Self, MathError>
Computes the tensor product of self with other.
Parameters:
other: the value with which the tensor product is computed.
Returns the tensor product of self with other or an error if the
moduli of the provided matrices mismatch.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat_1 = MatZq::from_str("[[1, 1],[2, 2]] mod 7").unwrap();
let mat_2 = MatZq::from_str("[[1, 2],[3, 4]] mod 7").unwrap();
let mat_ab = mat_1.tensor_product_safe(&mat_2).unwrap();
let mat_ba = mat_2.tensor_product_safe(&mat_1).unwrap();
let res_ab = "[[1, 2, 1, 2],[3, 4, 3, 4],[2, 4, 2, 4],[6, 1, 6, 1]] mod 7";
let res_ba = "[[1, 1, 2, 2],[2, 2, 4, 4],[3, 3, 4, 4],[6, 6, 1, 1]] mod 7";
assert_eq!(mat_ab, MatZq::from_str(res_ab).unwrap());
assert_eq!(mat_ba, MatZq::from_str(res_ba).unwrap());§Errors and Failures
- Returns a
MathErrorof typeMismatchingModulusif the moduli of the provided matrices mismatch.
Source§impl MatZq
impl MatZq
Sourcepub fn to_utf8(&self) -> Result<String, FromUtf8Error>
pub fn to_utf8(&self) -> Result<String, FromUtf8Error>
Enables conversion to a UTF8-Encoded String for MatZq values.
Every entry is padded with 00s s.t. all entries contain the same number of bytes.
Afterwards, they are appended row-by-row and converted.
The inverse to this function is MatZq::from_utf8 for valid UTF8-Encodings.
Warning: Not every byte-sequence forms a valid UTF8-Encoding.
In these cases, an error is returned. Please check the format of your message again.
The matrix entries are evaluated row by row, i.e. in the order of the output of mat_zq.to_string().
Returns the corresponding UTF8-encoded String or a
FromUtf8Error if the byte sequence contains an invalid UTF8-character.
§Examples
use qfall_math::integer::MatZ;
use std::str::FromStr;
let matrix = MatZ::from_str("[[104, 101, 108],[108, 111, 33]]").unwrap();
let message = matrix.to_utf8().unwrap();
assert_eq!("hello!", message);§Errors and Failures
- Returns a
FromUtf8Errorif the integer’s byte sequence contains invalid UTF8-characters.
Source§impl MatZq
impl MatZq
Sourcepub fn pretty_string(
&self,
nr_printed_rows: u64,
nr_printed_columns: u64,
) -> String
pub fn pretty_string( &self, nr_printed_rows: u64, nr_printed_columns: u64, ) -> String
Outputs the matrix as a String, where the upper leftmost nr_printed_rows x nr_printed_columns
submatrix is output entirely as well as the corresponding entries in the last column and row of the matrix.
Parameters:
nr_printed_rows: defines the number of rows of the upper leftmost matrix that are printed entirelynr_printed_columns: defines the number of columns of the upper leftmost matrix that are printed entirely
Returns a String representing the abbreviated matrix.
§Example
use qfall_math::integer::MatZ;
let matrix = MatZ::identity(10, 10);
println!("Matrix: {}", matrix.pretty_string(2, 2));
// outputs the following:
// Matrix: [
// [1, 0, , ..., 0],
// [0, 1, , ..., 0],
// [...],
// [0, 0, , ..., 1]
// ]Source§impl MatZq
impl MatZq
Sourcepub fn trace(&self) -> Result<Zq, MathError>
pub fn trace(&self) -> Result<Zq, MathError>
Returns the trace of a matrix and an error, if the matrix is not square.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let matrix = MatZq::from_str("[[1, 2],[3, 4]] mod 5").unwrap();
let trace = matrix.trace().unwrap();§Errors and Failures
- Returns a
MathErrorof typeNoSquareMatrixif the matrix is not a square matrix
Source§impl MatZq
impl MatZq
Sourcepub fn transpose(&self) -> Self
pub fn transpose(&self) -> Self
Returns the transposed form of the given matrix, i.e. rows get transformed to columns and vice versa.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let mat = MatZq::from_str("[[2, 1],[2, 1],[2, 1]] mod 4").unwrap();
let cmp = MatZq::from_str("[[2, 2, 2],[1, 1, 1]] mod 4").unwrap();
assert_eq!(mat.transpose(), cmp);Source§impl MatZq
impl MatZq
Sourcepub unsafe fn get_fmpz_mod_mat_struct(&mut self) -> &mut fmpz_mod_mat_struct
pub unsafe fn get_fmpz_mod_mat_struct(&mut self) -> &mut fmpz_mod_mat_struct
Returns a mutable reference to the field matrix of type fmpz_mod_mat_struct.
WARNING: The returned struct is part of flint_sys.
Any changes to this object are unsafe and may introduce memory leaks.
This function is a passthrough to enable users of this library to use flint_sys
and with that FLINT functions that might not be covered in our library yet.
If this is the case, please consider contributing to this open-source project
by opening a Pull Request at qfall_math
to provide this feature in the future.
§Safety
Any flint_sys struct and function is part of a FFI to the C-library FLINT.
As FLINT is a C-library, it does not provide all memory safety features
that Rust and our Wrapper provide.
Thus, using functions of flint_sys can introduce memory leaks.
Source§impl MatZq
impl MatZq
Sourcepub unsafe fn get_fmpz_mod_ctx(&mut self) -> &mut fmpz_mod_ctx
pub unsafe fn get_fmpz_mod_ctx(&mut self) -> &mut fmpz_mod_ctx
Returns a mutable reference to the underlying fmpz_mod_ctx by calling get_fmpz_mod_ctx on modulus.
WARNING: The returned struct is part of flint_sys.
Any changes to this object are unsafe and may introduce memory leaks.
In case you are calling this function to a modulus struct,
please be aware that most moduli are shared across multiple instances and all
modifications of this struct will affect any other instance with a reference to this object.
This function is a passthrough to enable users of this library to use flint_sys
and with that FLINT functions that might not be covered in our library yet.
If this is the case, please consider contributing to this open-source project
by opening a Pull Request at qfall_math
to provide this feature in the future.
§Safety
Any flint_sys struct and function is part of a FFI to the C-library FLINT.
As FLINT is a C-library, it does not provide all memory safety features
that Rust and our Wrapper provide.
Thus, using functions of flint_sys can introduce memory leaks.
Source§impl MatZq
impl MatZq
Sourcepub unsafe fn set_fmpz_mod_mat_struct(
&mut self,
flint_struct: fmpz_mod_mat_struct,
)
pub unsafe fn set_fmpz_mod_mat_struct( &mut self, flint_struct: fmpz_mod_mat_struct, )
Sets the field matrix of type fmpz_mod_mat_struct to flint_struct.
Parameters:
flint_struct: value to set the attribute to
This function is a passthrough to enable users of this library to use flint_sys
and with that FLINT functions that might not be covered in our library yet.
If this is the case, please consider contributing to this open-source project
by opening a Pull Request at qfall_math
to provide this feature in the future.
§Safety
Ensure that the old struct does not share any memory with any other structs that might be used in the future. The memory of the old struct is freed using this function.
Any flint_sys struct and function is part of a FFI to the C-library FLINT.
As FLINT is a C-library, it does not provide all memory safety features
that Rust and our Wrapper provide.
Thus, using functions of flint_sys can introduce memory leaks.
Source§impl MatZq
impl MatZq
Sourcepub unsafe fn set_fmpz_mod_ctx(&mut self, flint_struct: fmpz_mod_ctx)
pub unsafe fn set_fmpz_mod_ctx(&mut self, flint_struct: fmpz_mod_ctx)
Sets the field fmpz_mod_ctx to flint_struct by calling set_fmpz_mod_ctx on modulus.
Parameters:
flint_struct: value to set the attribute to
This function is a passthrough to enable users of this library to use flint_sys
and with that FLINT functions that might not be covered in our library yet.
If this is the case, please consider contributing to this open-source project
by opening a Pull Request at qfall_math
to provide this feature in the future.
§Safety
Ensure that the old struct does not share any memory with any other structs that might be used in the future. The memory of the old struct is freed using this function.
Any flint_sys struct and function is part of a FFI to the C-library FLINT.
As FLINT is a C-library, it does not provide all memory safety features
that Rust and our Wrapper provide.
Thus, using functions of flint_sys can introduce memory leaks.
Source§impl MatZq
impl MatZq
Sourcepub fn dot_product(&self, other: &Self) -> Result<Zq, MathError>
pub fn dot_product(&self, other: &Self) -> Result<Zq, MathError>
Returns the dot product of two vectors of type MatZq.
Note that the dimensions of the two vectors are irrelevant for the dot product.
Parameters:
other: specifies the other vector the dot product is calculated over
Returns the resulting dot product as a Zq or an error
if the given MatZq instances aren’t vectors, have different
numbers of entries, or mismatching moduli.
§Examples
use qfall_math::integer_mod_q::{MatZq, Zq};
use std::str::FromStr;
let vec_1 = MatZq::from_str("[[1],[2],[3]] mod 5").unwrap();
let vec_2 = MatZq::from_str("[[1, 3, 2]] mod 5").unwrap();
let dot_prod = vec_1.dot_product(&vec_2).unwrap();
// 1*1 + 2*3 + 3*2 = 3 mod 5
assert_eq!(Zq::from((3, 5)), dot_prod);§Errors and Failures
- Returns a
MathErrorof typeVectorFunctionCalledOnNonVectorif the givenMatZqinstance is not a (row or column) vector. - Returns a
MathErrorof typeMismatchingMatrixDimensionif the given vectors have different lengths. - Returns a
MathErrorof typeMismatchingModulusif the provided matrices have different moduli.
Source§impl MatZq
impl MatZq
Sourcepub fn is_row_vector(&self) -> bool
pub fn is_row_vector(&self) -> bool
Returns true if the provided MatZq has only one row,
i.e. is a row vector. Otherwise, returns false.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let col_vec = MatZq::from_str("[[1],[2],[3]] mod 4").unwrap();
let row_vec = MatZq::from_str("[[1, 2, 3]] mod 4").unwrap();
assert!(row_vec.is_row_vector());
assert!(!col_vec.is_row_vector());Sourcepub fn is_column_vector(&self) -> bool
pub fn is_column_vector(&self) -> bool
Returns true if the provided MatZq has only one column,
i.e. is a column vector. Otherwise, returns false.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let col_vec = MatZq::from_str("[[1],[2],[3]] mod 4").unwrap();
let row_vec = MatZq::from_str("[[1, 2, 3]] mod 4").unwrap();
assert!(col_vec.is_column_vector());
assert!(!row_vec.is_column_vector());Sourcepub fn is_vector(&self) -> bool
pub fn is_vector(&self) -> bool
Returns true if the provided MatZq has only one column or one row,
i.e. is a vector. Otherwise, returns false.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let col_vec = MatZq::from_str("[[1],[2],[3]] mod 4").unwrap();
let row_vec = MatZq::from_str("[[1, 2, 3]] mod 4").unwrap();
assert!(col_vec.is_vector());
assert!(row_vec.is_vector());Sourcepub fn has_single_entry(&self) -> bool
pub fn has_single_entry(&self) -> bool
Source§impl MatZq
impl MatZq
Sourcepub fn norm_eucl_sqrd(&self) -> Result<Z, MathError>
pub fn norm_eucl_sqrd(&self) -> Result<Z, MathError>
Returns the squared Euclidean norm or squared 2-norm of the given (row or column) vector or an error if the given matrix is not a vector.
Each length of an entry is defined as the shortest distance to the next zero instance in the ring.
§Examples
use qfall_math::{
integer::Z,
integer_mod_q::MatZq,
};
use std::str::FromStr;
let vec = MatZq::from_str("[[1],[2],[3]] mod 4").unwrap();
let sqrd_2_norm = vec.norm_eucl_sqrd().unwrap();
// 1*1 + 2*2 + 1*1 = 6
assert_eq!(Z::from(6), sqrd_2_norm);§Errors and Failures
- Returns a
MathErrorof typeMathError::VectorFunctionCalledOnNonVectorif the givenMatZqinstance is not a (row or column) vector.
Sourcepub fn norm_eucl(&self) -> Result<Q, MathError>
pub fn norm_eucl(&self) -> Result<Q, MathError>
Returns the Euclidean norm or 2-norm of the given (row or column) vector or an error if the given matrix is not a vector.
Each length of an entry is defined as the shortest distance to the next zero instance in the ring.
§Examples
use qfall_math::{
rational::Q,
integer_mod_q::MatZq,
};
use std::str::FromStr;
let vec = MatZq::from_str("[[2],[2],[2],[2]] mod 4").unwrap();
let eucl_norm = vec.norm_eucl().unwrap();
// sqrt(4 * 2^2) = 4
assert_eq!(Q::from(4), eucl_norm);§Errors and Failures
- Returns a
MathErrorof typeMathError::VectorFunctionCalledOnNonVectorif the givenMatZqinstance is not a (row or column) vector.
Sourcepub fn norm_infty(&self) -> Result<Z, MathError>
pub fn norm_infty(&self) -> Result<Z, MathError>
Returns the infinity norm or ∞-norm of the given (row or column) vector or an error if the given matrix is not a vector.
Each length of an entry is defined as the shortest distance to the next zero instance in the ring.
§Examples
use qfall_math::{
integer::Z,
integer_mod_q::MatZq,
};
use std::str::FromStr;
let vec = MatZq::from_str("[[1],[2],[3]] mod 3").unwrap();
let infty_norm = vec.norm_infty().unwrap();
// max{1, 1, 0} = 1
assert_eq!(Z::ONE, infty_norm);§Errors and Failures
- Returns a
MathErrorof typeMathError::VectorFunctionCalledOnNonVectorif the givenMatZqinstance is not a (row or column) vector.
Trait Implementations§
Source§impl Add<&MatZ> for &MatZq
impl Add<&MatZ> for &MatZq
Source§fn add(self, other: &MatZ) -> Self::Output
fn add(self, other: &MatZ) -> Self::Output
Implements the Add trait for a MatZ and a MatZq matrix.
Add is implemented for any combination of MatZ and MatZq and vice versa.
Parameters:
other: specifies the value to add toself
Returns the sum of both numbers as a MatZq.
§Examples
use qfall_math::{integer::MatZ, integer_mod_q::MatZq};
use std::str::FromStr;
let a = MatZ::from_str("[[1, 2, 3],[3, 4, 5]]").unwrap();
let b = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c = &a + &b;
let d = a.clone() + b.clone();
let e = &b + &a;
let f = b + a;§Panics …
- if the dimensions of both matrices mismatch.
Source§impl Add for &MatZq
impl Add for &MatZq
Source§fn add(self, other: Self) -> Self::Output
fn add(self, other: Self) -> Self::Output
Implements the Add trait for two MatZq values.
Add is implemented for any combination of MatZq and borrowed MatZq.
Parameters:
other: specifies the value to add toself
Returns the sum of both numbers as a MatZq.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a: MatZq = MatZq::from_str("[[1, 2, 3],[3, 4, 5]] mod 7").unwrap();
let b: MatZq = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c: MatZq = &a + &b;
let d: MatZq = a + b;
let e: MatZq = &c + d;
let f: MatZq = c + &e;§Panics …
- if the dimensions of both matrices mismatch.
- if the moduli mismatch.
Source§impl AddAssign<&MatZ> for MatZq
impl AddAssign<&MatZ> for MatZq
Source§fn add_assign(&mut self, other: &MatZ)
fn add_assign(&mut self, other: &MatZ)
Documentation at MatZq::add_assign.
Source§impl AddAssign<&MatZq> for MatZq
impl AddAssign<&MatZq> for MatZq
Source§fn add_assign(&mut self, other: &Self)
fn add_assign(&mut self, other: &Self)
Computes the addition of self and other reusing
the memory of self.
AddAssign can be used on MatZq in combination with
MatZq and MatZ.
Parameters:
other: specifies the value to add toself
§Examples
use qfall_math::{integer_mod_q::MatZq, integer::MatZ};
let mut a = MatZq::identity(2, 2, 7);
let b = MatZq::new(2, 2, 7);
let c = MatZ::new(2, 2);
a += &b;
a += b;
a += &c;
a += c;§Panics …
- if the matrix dimensions mismatch.
- if the moduli of the matrices mismatch.
Source§impl AddAssign<MatZ> for MatZq
impl AddAssign<MatZ> for MatZq
Source§fn add_assign(&mut self, other: MatZ)
fn add_assign(&mut self, other: MatZ)
Documentation at MatZq::add_assign.
Source§impl AddAssign for MatZq
impl AddAssign for MatZq
Source§fn add_assign(&mut self, other: MatZq)
fn add_assign(&mut self, other: MatZq)
Documentation at MatZq::add_assign.
Source§impl Clone for MatZq
impl Clone for MatZq
Source§impl CompareBase<&MatZ> for MatZq
impl CompareBase<&MatZ> for MatZq
Source§impl CompareBase<&MatZq> for MatNTTPolynomialRingZq
impl CompareBase<&MatZq> for MatNTTPolynomialRingZq
Source§fn compare_base(&self, other: &&MatZq) -> bool
fn compare_base(&self, other: &&MatZq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &&MatZq) -> Option<MathError>
fn call_compare_base_error(&self, other: &&MatZq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl CompareBase<&MatZq> for MatPolynomialRingZq
impl CompareBase<&MatZq> for MatPolynomialRingZq
Source§fn compare_base(&self, other: &&MatZq) -> bool
fn compare_base(&self, other: &&MatZq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &&MatZq) -> Option<MathError>
fn call_compare_base_error(&self, other: &&MatZq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl CompareBase<&MatZq> for MatZq
impl CompareBase<&MatZq> for MatZq
Source§fn compare_base(&self, other: &&MatZq) -> bool
fn compare_base(&self, other: &&MatZq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &&MatZq) -> Option<MathError>
fn call_compare_base_error(&self, other: &&MatZq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl CompareBase<&Zq> for MatZq
impl CompareBase<&Zq> for MatZq
Source§fn compare_base(&self, other: &&Zq) -> bool
fn compare_base(&self, other: &&Zq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &&Zq) -> Option<MathError>
fn call_compare_base_error(&self, other: &&Zq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl<Integer: Into<Z>> CompareBase<Integer> for MatZq
impl<Integer: Into<Z>> CompareBase<Integer> for MatZq
Source§impl CompareBase<MatZ> for MatZq
impl CompareBase<MatZ> for MatZq
Source§impl CompareBase<MatZq> for MatNTTPolynomialRingZq
impl CompareBase<MatZq> for MatNTTPolynomialRingZq
Source§fn compare_base(&self, other: &MatZq) -> bool
fn compare_base(&self, other: &MatZq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &MatZq) -> Option<MathError>
fn call_compare_base_error(&self, other: &MatZq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl CompareBase<MatZq> for MatPolynomialRingZq
impl CompareBase<MatZq> for MatPolynomialRingZq
Source§fn compare_base(&self, other: &MatZq) -> bool
fn compare_base(&self, other: &MatZq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &MatZq) -> Option<MathError>
fn call_compare_base_error(&self, other: &MatZq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl CompareBase<Zq> for MatZq
impl CompareBase<Zq> for MatZq
Source§fn compare_base(&self, other: &Zq) -> bool
fn compare_base(&self, other: &Zq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &Zq) -> Option<MathError>
fn call_compare_base_error(&self, other: &Zq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl CompareBase for MatZq
impl CompareBase for MatZq
Source§fn compare_base(&self, other: &MatZq) -> bool
fn compare_base(&self, other: &MatZq) -> bool
Compares the moduli of the two elements.
Parameters:
other: The other object whose base is compared toself
Returns true if the moduli match and false otherwise.
Source§fn call_compare_base_error(&self, other: &MatZq) -> Option<MathError>
fn call_compare_base_error(&self, other: &MatZq) -> Option<MathError>
Returns an error that gives a small explanation of how the moduli are incomparable.
Parameters:
other: The other object whose base is compared toself
Returns a MathError of type MismatchingModulus.
Source§impl Concatenate for &MatZq
impl Concatenate for &MatZq
Source§fn concat_vertical(self, other: Self) -> Result<Self::Output, MathError>
fn concat_vertical(self, other: Self) -> Result<Self::Output, MathError>
Concatenates self with other vertically, i.e. other is added below.
Parameters:
other: the other matrix to concatenate withself
Returns a vertical concatenation of the two matrices or a an error, if the matrices can not be concatenated vertically.
§Examples
use qfall_math::traits::*;
use qfall_math::integer_mod_q::MatZq;
let mat_1 = MatZq::new(13, 5, 19);
let mat_2 = MatZq::new(17, 5, 19);
let mat_vert = mat_1.concat_vertical(&mat_2).unwrap();§Errors and Failures
- Returns a
MathErrorof typeMismatchingMatrixDimensionif the matrices can not be concatenated due to mismatching dimensions. - Returns a
MathErrorof typeMismatchingModulusif the matrices can not be concatenated due to mismatching moduli.
Source§fn concat_horizontal(self, other: Self) -> Result<Self::Output, MathError>
fn concat_horizontal(self, other: Self) -> Result<Self::Output, MathError>
Concatenates self with other horizontally, i.e. other is added on the right.
Parameters:
other: the other matrix to concatenate withself
Returns a horizontal concatenation of the two matrices or a an error, if the matrices can not be concatenated horizontally.
§Examples
use qfall_math::traits::*;
use qfall_math::integer_mod_q::MatZq;
let mat_1 = MatZq::new(17, 5, 19);
let mat_2 = MatZq::new(17, 6, 19);
let mat_vert = mat_1.concat_horizontal(&mat_2).unwrap();§Errors and Failures
- Returns a
MathErrorof typeMismatchingMatrixDimensionif the matrices can not be concatenated due to mismatching dimensions. - Returns a
MathErrorof typeMismatchingModulusif the matrices can not be concatenated due to mismatching moduli.
type Output = MatZq
Source§impl<'de> Deserialize<'de> for MatZq
impl<'de> Deserialize<'de> for MatZq
Source§fn deserialize<D>(deserializer: D) -> Result<Self, D::Error>where
D: Deserializer<'de>,
fn deserialize<D>(deserializer: D) -> Result<Self, D::Error>where
D: Deserializer<'de>,
Implements the deserialize option. This allows to create a MatZq from a given Json-object.
Source§impl Display for MatZq
impl Display for MatZq
Source§fn fmt(&self, f: &mut Formatter<'_>) -> Result
fn fmt(&self, f: &mut Formatter<'_>) -> Result
Allows to convert a matrix of type MatZq into a String.
Returns the Matrix in form of a String. For matrix [[1, 2, 3],[4, 5, 6]] mod 4
the String looks like this [[1, 2, 3],[0, 1, 2]] mod 4.
§Examples
use qfall_math::integer_mod_q::MatZq;
use core::fmt;
use std::str::FromStr;
let matrix = MatZq::from_str("[[1, 2, 3],[4, 5, 6]] mod 4").unwrap();
println!("{matrix}");use qfall_math::integer_mod_q::MatZq;
use core::fmt;
use std::str::FromStr;
let matrix = MatZq::from_str("[[1, 2, 3],[4, 5, 6]] mod 4").unwrap();
let matrix_string = matrix.to_string();Source§impl Drop for MatZq
impl Drop for MatZq
Source§fn drop(&mut self)
fn drop(&mut self)
Drops the given MatZq value and frees the allocated memory.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let str_1 = "[[1, 2, 3],[4, 5, 6]] mod 4";
{
let a = MatZq::from_str(str_1).unwrap();
} // as a's scope ends here, it get's droppeduse qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let str_1 = "[[1, 2, 3],[4, 5, 6]] mod 4";
let a = MatZq::from_str(str_1).unwrap();
drop(a); // explicitly drops a's valueSource§impl From<&MatZq> for String
impl From<&MatZq> for String
Source§fn from(value: &MatZq) -> Self
fn from(value: &MatZq) -> Self
Converts a MatZq into its String representation.
Parameters:
value: specifies the matrix that will be represented as aString
Returns a String of the form "[[row_0],[row_1],...[row_n]] mod q".
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let matrix = MatZq::from_str("[[6, 1],[5, 2]] mod 4").unwrap();
let string: String = matrix.into();Source§impl<Mod: Into<Modulus>> From<(&MatZ, Mod)> for MatZq
impl<Mod: Into<Modulus>> From<(&MatZ, Mod)> for MatZq
Source§fn from((matrix, modulus): (&MatZ, Mod)) -> Self
fn from((matrix, modulus): (&MatZ, Mod)) -> Self
Creates a MatZq from a MatZ and a value that implements Into<Modulus>.
Parameters:
matrix: the matrix from which the entries are takenmodulus: the modulus of the matrix
Returns a MatZq.
§Examples
use qfall_math::integer::MatZ;
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let m = MatZ::from_str("[[1, 2],[3, -1]]").unwrap();
let a = MatZq::from((&m, 17));Source§impl<Mod: Into<Modulus>> From<(MatZ, Mod)> for MatZq
impl<Mod: Into<Modulus>> From<(MatZ, Mod)> for MatZq
Source§fn from((matrix, modulus): (MatZ, Mod)) -> Self
fn from((matrix, modulus): (MatZ, Mod)) -> Self
Creates a MatZq from a MatZ and a value that implements Into<Modulus>.
Parameters:
matrix: the matrix from which the entries are takenmodulus: the modulus of the matrix
Returns a new MatZq matrix with entries from the MatZ instance modulo modulus.
§Examples
use qfall_math::integer::MatZ;
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let m = MatZ::from_str("[[1, 2],[3, -1]]").unwrap();
let a = MatZq::from((m, 17));Source§impl FromCoefficientEmbedding<&MatZq> for ModulusPolynomialRingZq
impl FromCoefficientEmbedding<&MatZq> for ModulusPolynomialRingZq
Source§fn from_coefficient_embedding(embedding: &MatZq) -> Self
fn from_coefficient_embedding(embedding: &MatZq) -> Self
Computes a polynomial from a vector.
The first i-th entry of the column vector is taken
as the coefficient of the polynomial.
It inverts the operation of
ModulusPolynomialRingZq::into_coefficient_embedding.
Parameters:
embedding: the column vector that encodes the embedding
Returns a polynomial that corresponds to the embedding.
§Examples
use std::str::FromStr;
use qfall_math::{
integer_mod_q::{MatZq, ModulusPolynomialRingZq},
traits::FromCoefficientEmbedding,
};
let vector = MatZq::from_str("[[17],[3],[-5]] mod 19").unwrap();
let poly = ModulusPolynomialRingZq::from_coefficient_embedding(&vector);
let cmp_poly = ModulusPolynomialRingZq::from_str("3 17 3 -5 mod 19").unwrap();
assert_eq!(cmp_poly, poly);§Panics …
- if the provided embedding is not a column vector.
Source§impl FromCoefficientEmbedding<&MatZq> for PolyOverZq
impl FromCoefficientEmbedding<&MatZq> for PolyOverZq
Source§fn from_coefficient_embedding(embedding: &MatZq) -> Self
fn from_coefficient_embedding(embedding: &MatZq) -> Self
Computes a polynomial from a vector.
The first i-th entry of the column vector is taken
as the coefficient of the polynomial.
It inverts the operation of
PolyOverZq::into_coefficient_embedding.
Parameters:
embedding: the column vector that encodes the embedding
Returns a polynomial that corresponds to the embedding.
§Examples
use std::str::FromStr;
use qfall_math::{
integer_mod_q::{MatZq, PolyOverZq},
traits::FromCoefficientEmbedding,
};
let vector = MatZq::from_str("[[17],[3],[-5]] mod 19").unwrap();
let poly = PolyOverZq::from_coefficient_embedding(&vector);
let cmp_poly = PolyOverZq::from_str("3 17 3 -5 mod 19").unwrap();
assert_eq!(cmp_poly, poly);§Panics …
- if the provided embedding is not a column vector.
Source§impl FromStr for MatZq
impl FromStr for MatZq
Source§fn from_str(string: &str) -> Result<Self, MathError>
fn from_str(string: &str) -> Result<Self, MathError>
Creates a MatZq matrix with entries in Zq from a String.
Parameters:
string: the matrix of form:"[[1, 2, 3],[4, 5, 6]] mod 4"for a 2x3 matrix with entries 1, 2, 3 in the first row, 4, 5, 6 in the second row and 4 as modulus.
Note that the strings for entries and the modulus are trimmed, i.e. all whitespaces around all values are ignored.
Returns a MatZq or an error if the matrix is not formatted in a suitable way,
the number of rows or columns is too large (must fit into i64),
the number of entries in rows is unequal or if the modulus or an entry is not formatted correctly.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let matrix = MatZq::from_str("[[1, 2, 3],[4, 5, 6]] mod 4").unwrap();use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let str_1 = "[[1, 2, 3],[4, 5, 6]] mod 4";
let matrix = MatZq::from_str(str_1).unwrap();use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let string = String::from("[[1, 2, 3],[4, 5, 6]] mod 4");
let matrix = MatZq::from_str(&string).unwrap();§Errors and Failures
- Returns a
MathErrorof typeStringConversionError- if the matrix is not formatted in a suitable way,
- if the number of rows or columns is too large (must fit into i64),
- if the number of entries in rows is unequal,
- if the delimiter
modcould not be found, or - if the modulus or an entry is not formatted correctly.
For further information see
Z::from_str.
§Panics …
- if the provided number of rows and columns or the modulus are not suited to create a matrix.
For further information see
MatZq::new. - if the modulus is smaller than
2.
Source§impl IntoCoefficientEmbedding<MatZq> for &ModulusPolynomialRingZq
impl IntoCoefficientEmbedding<MatZq> for &ModulusPolynomialRingZq
Source§fn into_coefficient_embedding(self, size: impl Into<i64>) -> MatZq
fn into_coefficient_embedding(self, size: impl Into<i64>) -> MatZq
Computes the coefficient embedding of the polynomial
in a MatZq as a column vector, where the i-th entry
of the vector corresponds to the i-th coefficient.
It inverts the operation of ModulusPolynomialRingZq::from_coefficient_embedding.
Parameters:
size: determines the number of rows of the embedding. It has to be larger than the degree of the polynomial.
Returns a coefficient embedding as a column vector if size is large enough.
§Examples
use std::str::FromStr;
use qfall_math::{
integer_mod_q::{MatZq, ModulusPolynomialRingZq},
traits::IntoCoefficientEmbedding,
};
let poly = ModulusPolynomialRingZq::from_str("3 17 3 -5 mod 19").unwrap();
let vector = poly.into_coefficient_embedding(4);
let cmp_vector = MatZq::from_str("[[17],[3],[-5],[0]] mod 19").unwrap();
assert_eq!(cmp_vector, vector);§Panics …
- if
sizeis not larger than the degree of the polynomial, i.e. not all coefficients can be embedded.
Source§impl IntoCoefficientEmbedding<MatZq> for &PolyOverZq
impl IntoCoefficientEmbedding<MatZq> for &PolyOverZq
Source§fn into_coefficient_embedding(self, size: impl Into<i64>) -> MatZq
fn into_coefficient_embedding(self, size: impl Into<i64>) -> MatZq
Computes the coefficient embedding of the polynomial
in a MatZq as a column vector, where the i-th entry
of the vector corresponds to the i-th coefficient.
It inverts the operation of PolyOverZq::from_coefficient_embedding.
Parameters:
size: determines the number of rows of the embedding. It has to be larger than the degree of the polynomial.
Returns a coefficient embedding as a column vector if size is large enough.
§Examples
use std::str::FromStr;
use qfall_math::{
integer_mod_q::{MatZq, PolyOverZq},
traits::IntoCoefficientEmbedding,
};
let poly = PolyOverZq::from_str("3 17 3 -5 mod 19").unwrap();
let vector = poly.into_coefficient_embedding(4);
let cmp_vector = MatZq::from_str("[[17],[3],[-5],[0]] mod 19").unwrap();
assert_eq!(cmp_vector, vector);§Panics …
- if
sizeis not larger than the degree of the polynomial, i.e. not all coefficients can be embedded.
Source§impl MatrixDimensions for MatZq
impl MatrixDimensions for MatZq
Source§fn get_num_rows(&self) -> i64
fn get_num_rows(&self) -> i64
Source§impl MatrixGetEntry<Z> for MatZq
impl MatrixGetEntry<Z> for MatZq
Source§unsafe fn get_entry_unchecked(&self, row: i64, column: i64) -> Z
unsafe fn get_entry_unchecked(&self, row: i64, column: i64) -> Z
Outputs the Z value of a specific matrix entry
without checking whether it’s part of the matrix.
Parameters:
row: specifies the row in which the entry is locatedcolumn: specifies the column in which the entry is located
Returns the Z value of the matrix at the position of the given
row and column.
§Safety
To use this function safely, make sure that the selected entry is part of the matrix. If it is not, memory leaks, unexpected panics, etc. might occur.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::traits::MatrixGetEntry;
use qfall_math::integer::Z;
use std::str::FromStr;
let matrix = MatZq::from_str("[[1, 2, 3],[4, 5, 6],[7, 8, 9]] mod 10").unwrap();
let entry_1 :Z = unsafe { matrix.get_entry_unchecked(0, 2) };
let entry_2 :Z = unsafe { matrix.get_entry_unchecked(2, 1) };
let entry_3 :Z = unsafe { matrix.get_entry_unchecked(2, 1) };
assert_eq!(3, entry_1);
assert_eq!(8, entry_2);
assert_eq!(8, entry_3);Source§fn get_entry(
&self,
row: impl TryInto<i64> + Display,
column: impl TryInto<i64> + Display,
) -> Result<T, MathError>
fn get_entry( &self, row: impl TryInto<i64> + Display, column: impl TryInto<i64> + Display, ) -> Result<T, MathError>
Source§fn get_entries(&self) -> Vec<Vec<T>>
fn get_entries(&self) -> Vec<Vec<T>>
Vec<Vec<T>> containing all entries of the matrix s.t.
any entry in row i and column j can be accessed via entries[i][j]
if entries = matrix.get_entries. Read moreSource§fn get_entries_rowwise(&self) -> Vec<T>
fn get_entries_rowwise(&self) -> Vec<T>
Source§impl MatrixGetEntry<Zq> for MatZq
impl MatrixGetEntry<Zq> for MatZq
Source§unsafe fn get_entry_unchecked(&self, row: i64, column: i64) -> Zq
unsafe fn get_entry_unchecked(&self, row: i64, column: i64) -> Zq
Outputs the Zq value of a specific matrix entry
without checking whether it’s part of the matrix.
Parameters:
row: specifies the row in which the entry is locatedcolumn: specifies the column in which the entry is located
Returns the Zq value of the matrix at the position of the given
row and column.
§Safety
To use this function safely, make sure that the selected entry is part of the matrix. If it is not, memory leaks, unexpected panics, etc. might occur.
§Examples
use qfall_math::integer_mod_q::{MatZq, Zq};
use qfall_math::traits::MatrixGetEntry;
use std::str::FromStr;
let matrix = MatZq::from_str("[[1, 2, 3],[4, 5, 6],[7, 8, 9]] mod 10").unwrap();
assert_eq!(Zq::from((3, 10)), unsafe { matrix.get_entry_unchecked(0, 2) } );
assert_eq!(Zq::from((8, 10)), unsafe { matrix.get_entry_unchecked(2, 1) } );
assert_eq!(Zq::from((8, 10)), unsafe { matrix.get_entry_unchecked(2, 1) } );Source§fn get_entry(
&self,
row: impl TryInto<i64> + Display,
column: impl TryInto<i64> + Display,
) -> Result<T, MathError>
fn get_entry( &self, row: impl TryInto<i64> + Display, column: impl TryInto<i64> + Display, ) -> Result<T, MathError>
Source§fn get_entries(&self) -> Vec<Vec<T>>
fn get_entries(&self) -> Vec<Vec<T>>
Vec<Vec<T>> containing all entries of the matrix s.t.
any entry in row i and column j can be accessed via entries[i][j]
if entries = matrix.get_entries. Read moreSource§fn get_entries_rowwise(&self) -> Vec<T>
fn get_entries_rowwise(&self) -> Vec<T>
Source§impl MatrixGetSubmatrix for MatZq
impl MatrixGetSubmatrix for MatZq
Source§unsafe fn get_submatrix_unchecked(
&self,
row_1: i64,
row_2: i64,
col_1: i64,
col_2: i64,
) -> Self
unsafe fn get_submatrix_unchecked( &self, row_1: i64, row_2: i64, col_1: i64, col_2: i64, ) -> Self
Returns a deep copy of the submatrix defined by the given parameters and does not check the provided dimensions. There is also a safe version of this function that checks the input.
Parameters:
row_1: the starting row of the submatrix
row_2: the ending row of the submatrix
col_1: the starting column of the submatrix
col_2: the ending column of the submatrix
Returns the submatrix from (row_1, col_1) to (row_2, col_2)(exclusively).
§Examples
use qfall_math::{integer_mod_q::MatZq, traits::MatrixGetSubmatrix};
use std::str::FromStr;
let mat = MatZq::identity(3, 3, 17);
let sub_mat_1 = mat.get_submatrix(0, 2, 1, 1).unwrap();
let sub_mat_2 = mat.get_submatrix(0, -1, 1, -2).unwrap();
let sub_mat_3 = unsafe{mat.get_submatrix_unchecked(0, 3, 1, 2)};
let e_2 = MatZq::from_str("[[0],[1],[0]] mod 17").unwrap();
assert_eq!(e_2, sub_mat_1);
assert_eq!(e_2, sub_mat_2);
assert_eq!(e_2, sub_mat_3);§Safety
To use this function safely, make sure that the selected submatrix is part of the matrix. If it is not, memory leaks, unexpected panics, etc. might occur.
Source§fn get_row(
&self,
row: impl TryInto<i64> + Display + Clone,
) -> Result<Self, MathError>
fn get_row( &self, row: impl TryInto<i64> + Display + Clone, ) -> Result<Self, MathError>
Source§unsafe fn get_row_unchecked(&self, row: i64) -> Self
unsafe fn get_row_unchecked(&self, row: i64) -> Self
Source§fn get_column(
&self,
column: impl TryInto<i64> + Display + Clone,
) -> Result<Self, MathError>
fn get_column( &self, column: impl TryInto<i64> + Display + Clone, ) -> Result<Self, MathError>
Source§unsafe fn get_column_unchecked(&self, column: i64) -> Self
unsafe fn get_column_unchecked(&self, column: i64) -> Self
Source§fn get_submatrix(
&self,
row_1: impl TryInto<i64> + Display,
row_2: impl TryInto<i64> + Display,
col_1: impl TryInto<i64> + Display,
col_2: impl TryInto<i64> + Display,
) -> Result<Self, MathError>
fn get_submatrix( &self, row_1: impl TryInto<i64> + Display, row_2: impl TryInto<i64> + Display, col_1: impl TryInto<i64> + Display, col_2: impl TryInto<i64> + Display, ) -> Result<Self, MathError>
(row_1, col_1) to (row_2, col_2)(inclusively) are collected in
a new matrix.
Note that row_1 >= row_2 and col_1 >= col_2 must hold after converting negative indices.
Otherwise the function will panic. Read moreSource§impl MatrixSetEntry<&Zq> for MatZq
impl MatrixSetEntry<&Zq> for MatZq
Source§unsafe fn set_entry_unchecked(&mut self, row: i64, column: i64, value: &Zq)
unsafe fn set_entry_unchecked(&mut self, row: i64, column: i64, value: &Zq)
Sets the value of a specific matrix entry according to a given value of type Zq
without checking whether the coordinate is part of the matrix,
if the moduli match or the entry is reduced.
Parameters:
row: specifies the row in which the entry is locatedcolumn: specifies the column in which the entry is locatedvalue: specifies the value to which the entry is set
§Safety
To use this function safely, make sure that the selected entry is part of the matrix. If it is not, memory leaks, unexpected panics, etc. might occur.
§Examples
use qfall_math::integer_mod_q::{MatZq, Zq};
use qfall_math::traits::*;
let mut matrix = MatZq::new(3, 3, 10);
let value = Zq::from((5, 10));
unsafe {
matrix.set_entry_unchecked(0, 1, &value);
matrix.set_entry_unchecked(2, 2, Zq::from((19, 10)));
}
assert_eq!("[[0, 5, 0],[0, 0, 0],[0, 0, 9]] mod 10", matrix.to_string());Source§impl<Integer: Into<Z>> MatrixSetEntry<Integer> for MatZq
impl<Integer: Into<Z>> MatrixSetEntry<Integer> for MatZq
Source§fn set_entry(
&mut self,
row: impl TryInto<i64> + Display,
column: impl TryInto<i64> + Display,
value: Integer,
) -> Result<(), MathError>
fn set_entry( &mut self, row: impl TryInto<i64> + Display, column: impl TryInto<i64> + Display, value: Integer, ) -> Result<(), MathError>
Sets the value of a specific matrix entry according to a given value
that implements Into<Z>.
Parameters:
row: specifies the row in which the entry is locatedcolumn: specifies the column in which the entry is locatedvalue: specifies the value to which the entry is set
Negative indices can be used to index from the back, e.g., -1 for
the last element.
Returns an empty Ok if the action could be performed successfully.
Otherwise, a MathError is returned if the specified entry is not part of the matrix.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::traits::*;
let mut matrix = MatZq::new(3, 3, 10);
matrix.set_entry(0, 1, 5).unwrap();
matrix.set_entry(-1, 2, 19).unwrap();
assert_eq!("[[0, 5, 0],[0, 0, 0],[0, 0, 9]] mod 10", matrix.to_string());§Errors and Failures
- Returns a
MathErrorof typeMathError::OutOfBoundsifroworcolumnare greater than the matrix size.
Source§unsafe fn set_entry_unchecked(&mut self, row: i64, column: i64, value: Integer)
unsafe fn set_entry_unchecked(&mut self, row: i64, column: i64, value: Integer)
Sets the value of a specific matrix entry according to a given value
that implements Into<Z> without checking whether the coordinate is part of the matrix,
if the moduli match or the entry is reduced.
Parameters:
row: specifies the row in which the entry is locatedcolumn: specifies the column in which the entry is locatedvalue: specifies the value to which the entry is set
§Safety
To use this function safely, make sure that the selected entry is part of the matrix. If it is not, memory leaks, unexpected panics, etc. might occur.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::traits::*;
let mut matrix = MatZq::new(3, 3, 10);
unsafe {
matrix.set_entry_unchecked(0, 1, 5);
matrix.set_entry_unchecked(2, 2, 19);
}
assert_eq!("[[0, 5, 0],[0, 0, 0],[0, 0, 19]] mod 10", matrix.to_string());Source§impl MatrixSetEntry<Zq> for MatZq
impl MatrixSetEntry<Zq> for MatZq
Source§fn set_entry(
&mut self,
row: impl TryInto<i64> + Display,
column: impl TryInto<i64> + Display,
value: Zq,
) -> Result<(), MathError>
fn set_entry( &mut self, row: impl TryInto<i64> + Display, column: impl TryInto<i64> + Display, value: Zq, ) -> Result<(), MathError>
Documentation can be found at MatZq::set_entry for &Zq.
Source§unsafe fn set_entry_unchecked(&mut self, row: i64, column: i64, value: Zq)
unsafe fn set_entry_unchecked(&mut self, row: i64, column: i64, value: Zq)
Documentation can be found at MatZq::set_entry for &Zq.
Source§impl MatrixSetSubmatrix for MatZq
impl MatrixSetSubmatrix for MatZq
Source§unsafe fn set_submatrix_unchecked(
&mut self,
row_self_start: i64,
col_self_start: i64,
row_self_end: i64,
col_self_end: i64,
other: &Self,
row_other_start: i64,
col_other_start: i64,
row_other_end: i64,
col_other_end: i64,
)
unsafe fn set_submatrix_unchecked( &mut self, row_self_start: i64, col_self_start: i64, row_self_end: i64, col_self_end: i64, other: &Self, row_other_start: i64, col_other_start: i64, row_other_end: i64, col_other_end: i64, )
Sets the matrix entries in self to entries defined in other.
The entries in self starting from (row_self_start, col_self_start) up to
(row_self_end, col_self_end)are set to be
the entries from the submatrix from other defined by (row_other_start, col_other_start)
to (row_other_end, col_other_end) (exclusively).
Parameters:
row_self_start: the starting row of the matrix in which to set a submatrix
col_self_start: the starting column of the matrix in which to set a submatrix
other: the matrix from where to take the submatrix to set
row_other_start: the starting row of the specified submatrix
col_other_start: the starting column of the specified submatrix
row_other_end: the ending row of the specified submatrix
col_other_end:the ending column of the specified submatrix
§Examples
use qfall_math::integer_mod_q::{MatZq, Modulus};
use qfall_math::integer::MatZ;
use qfall_math::traits::MatrixSetSubmatrix;
use std::str::FromStr;
let mat = MatZ::identity(3, 3);
let modulus = Modulus::from(17);
let mut mat = MatZq::from((&mat, &modulus));
mat.set_submatrix(0, 1, &mat.clone(), 0, 0, 1, 1).unwrap();
// [[1,1,0],[0,0,1],[0,0,1]]
let mat_cmp = MatZ::from_str("[[1, 1, 0],[0, 0, 1],[0, 0, 1]]").unwrap();
assert_eq!(mat, MatZq::from((&mat_cmp, &modulus)));
unsafe{ mat.set_submatrix_unchecked(2, 0, 3, 2, &mat.clone(), 0, 0, 1, 2) };
let mat_cmp = MatZ::from_str("[[1, 1, 0],[0, 0, 1],[1, 1, 1]]").unwrap();
assert_eq!(mat, MatZq::from((&mat_cmp, &modulus)));§Safety
To use this function safely, make sure that the selected submatrices are part of the matrices, the submatrices are of the same dimensions and the base types are the same. If not, memory leaks, unexpected panics, etc. might occur.
Source§fn set_row(
&mut self,
row_0: impl TryInto<i64> + Display,
other: &Self,
row_1: impl TryInto<i64> + Display,
) -> Result<(), MathError>
fn set_row( &mut self, row_0: impl TryInto<i64> + Display, other: &Self, row_1: impl TryInto<i64> + Display, ) -> Result<(), MathError>
other. Read moreSource§unsafe fn set_row_unchecked(&mut self, row_0: i64, other: &Self, row_1: i64)
unsafe fn set_row_unchecked(&mut self, row_0: i64, other: &Self, row_1: i64)
other. Read moreSource§fn set_column(
&mut self,
col_0: impl TryInto<i64> + Display,
other: &Self,
col_1: impl TryInto<i64> + Display,
) -> Result<(), MathError>
fn set_column( &mut self, col_0: impl TryInto<i64> + Display, other: &Self, col_1: impl TryInto<i64> + Display, ) -> Result<(), MathError>
other. Read moreSource§unsafe fn set_column_unchecked(&mut self, col_0: i64, other: &Self, col_1: i64)
unsafe fn set_column_unchecked(&mut self, col_0: i64, other: &Self, col_1: i64)
other. Read moreSource§fn set_submatrix(
&mut self,
row_self_start: impl TryInto<i64> + Display,
col_self_start: impl TryInto<i64> + Display,
other: &Self,
row_other_start: impl TryInto<i64> + Display,
col_other_start: impl TryInto<i64> + Display,
row_other_end: impl TryInto<i64> + Display,
col_other_end: impl TryInto<i64> + Display,
) -> Result<(), MathError>
fn set_submatrix( &mut self, row_self_start: impl TryInto<i64> + Display, col_self_start: impl TryInto<i64> + Display, other: &Self, row_other_start: impl TryInto<i64> + Display, col_other_start: impl TryInto<i64> + Display, row_other_end: impl TryInto<i64> + Display, col_other_end: impl TryInto<i64> + Display, ) -> Result<(), MathError>
self to entries defined in other.
The entries in self starting from (row_self_start, col_self_start) are set to be
the entries from the submatrix from other defined by (row_other_start, col_other_start)
to (row_other_end, col_other_end) (inclusively).
The original matrix must have sufficiently many entries to contain the defined submatrix. Read moreSource§impl MatrixSwaps for MatZq
impl MatrixSwaps for MatZq
Source§fn swap_entries(
&mut self,
row_0: impl TryInto<i64> + Display,
col_0: impl TryInto<i64> + Display,
row_1: impl TryInto<i64> + Display,
col_1: impl TryInto<i64> + Display,
) -> Result<(), MathError>
fn swap_entries( &mut self, row_0: impl TryInto<i64> + Display, col_0: impl TryInto<i64> + Display, row_1: impl TryInto<i64> + Display, col_1: impl TryInto<i64> + Display, ) -> Result<(), MathError>
Swaps two entries of the specified matrix.
Parameters:
row_0: specifies the row, in which the first entry is locatedcol_0: specifies the column, in which the first entry is locatedrow_1: specifies the row, in which the second entry is locatedcol_1: specifies the column, in which the second entry is located
Negative indices can be used to index from the back, e.g., -1 for
the last element.
Returns an empty Ok if the action could be performed successfully.
Otherwise, a MathError is returned if one of the specified entries is not part of the matrix.
§Examples
use qfall_math::{integer_mod_q::MatZq, traits::MatrixSwaps};
let mut matrix = MatZq::new(4, 3, 5);
matrix.swap_entries(0, 0, 2, 1);§Errors and Failures
- Returns a
MathErrorof typeMathError::OutOfBoundsif row or column are greater than the matrix size.
Source§fn swap_columns(
&mut self,
col_0: impl TryInto<i64> + Display,
col_1: impl TryInto<i64> + Display,
) -> Result<(), MathError>
fn swap_columns( &mut self, col_0: impl TryInto<i64> + Display, col_1: impl TryInto<i64> + Display, ) -> Result<(), MathError>
Swaps two columns of the specified matrix.
Parameters:
col_0: specifies the first column which is swapped with the second onecol_1: specifies the second column which is swapped with the first one
Negative indices can be used to index from the back, e.g., -1 for
the last element.
Returns an empty Ok if the action could be performed successfully.
Otherwise, a MathError is returned if one of the specified columns is not part of the matrix.
§Examples
use qfall_math::{integer_mod_q::MatZq, traits::MatrixSwaps};
let mut matrix = MatZq::new(4, 3, 5);
matrix.swap_columns(0, 2);§Errors and Failures
- Returns a
MathErrorof typeOutOfBoundsif one of the given columns is greater than the matrix.
Source§fn swap_rows(
&mut self,
row_0: impl TryInto<i64> + Display,
row_1: impl TryInto<i64> + Display,
) -> Result<(), MathError>
fn swap_rows( &mut self, row_0: impl TryInto<i64> + Display, row_1: impl TryInto<i64> + Display, ) -> Result<(), MathError>
Swaps two rows of the specified matrix.
Parameters:
row_0: specifies the first row which is swapped with the second onerow_1: specifies the second row which is swapped with the first one
Negative indices can be used to index from the back, e.g., -1 for
the last element.
Returns an empty Ok if the action could be performed successfully.
Otherwise, a MathError is returned if one of the specified rows is not part of the matrix.
§Examples
use qfall_math::{integer_mod_q::MatZq, traits::MatrixSwaps};
let mut matrix = MatZq::new(4, 3, 5);
matrix.swap_rows(0, 2);§Errors and Failures
- Returns a
MathErrorof typeOutOfBoundsif one of the given rows is greater than the matrix.
Source§impl Mul<&MatZ> for &MatZq
impl Mul<&MatZ> for &MatZq
Source§fn mul(self, other: &MatZ) -> Self::Output
fn mul(self, other: &MatZ) -> Self::Output
Implements the Mul trait for MatZq and MatZ.
Mul is implemented for any combination of owned and borrowed values.
Parameters:
other: specifies the value to multiply withself
Returns the product of self and other as a MatZq.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::integer::MatZ;
use std::str::FromStr;
let a = MatZq::from_str("[[2, 1],[1, 2]] mod 3").unwrap();
let b = MatZ::identity(2, 2);
let c = &a * &b;
let d = a * b;
let e = d * &MatZ::identity(2, 2);
let f = &e * MatZ::identity(2, 2);§Panics …
- if the dimensions of
selfandotherdo not match for multiplication.
Source§impl Mul<&MatZq> for &MatZ
impl Mul<&MatZq> for &MatZ
Source§fn mul(self, other: &MatZq) -> Self::Output
fn mul(self, other: &MatZq) -> Self::Output
Implements the Mul trait for MatZ and MatZq.
Mul is implemented for any combination of owned and borrowed values.
Parameters:
other: specifies the value to multiply withself
Returns the product of self and other as a MatZq.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::integer::MatZ;
use std::str::FromStr;
let a = MatZ::identity(2, 2);
let b = MatZq::from_str("[[2, 1],[1, 2]] mod 3").unwrap();
let c = &a * &b;
let d = a * b;
let e = &MatZ::identity(2, 2) * d;
let f = MatZ::identity(2, 2) * &e;§Panics …
- if the dimensions of
selfandotherdo not match for multiplication.
Source§impl Mul<&Z> for &MatZq
impl Mul<&Z> for &MatZq
Source§fn mul(self, scalar: &Z) -> Self::Output
fn mul(self, scalar: &Z) -> Self::Output
Implements the Mul trait for a MatZq matrix with a Z integer.
Mul is implemented for any combination of owned and borrowed values.
Parameters:
scalar: specifies the scalar by which the matrix is multiplied
Returns the product of self and scalar as a MatZq.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::integer::Z;
use std::str::FromStr;
let mat_1 = MatZq::from_str("[[42, 17],[8, 6]] mod 61").unwrap();
let integer = Z::from(3);
let mat_2 = &mat_1 * &integer;Source§impl Mul<&Zq> for &MatZq
impl Mul<&Zq> for &MatZq
Source§fn mul(self, scalar: &Zq) -> Self::Output
fn mul(self, scalar: &Zq) -> Self::Output
Implements the Mul trait for a MatZq matrix with a Zq.
Mul is implemented for any combination of owned and borrowed values.
Parameters:
scalar: specifies the scalar by which the matrix is multiplied
Returns the product of self and scalar as a MatZq.
§Examples
use qfall_math::integer_mod_q::{MatZq, Zq};
use std::str::FromStr;
let mat_1 = MatZq::from_str("[[42, 17],[8, 6]] mod 61").unwrap();
let integer = Zq::from((2, 61));
let mat_2 = &mat_1 * &integer;§Panics …
- if the moduli mismatch.
Source§impl Mul for &MatZq
impl Mul for &MatZq
Source§fn mul(self, other: Self) -> Self::Output
fn mul(self, other: Self) -> Self::Output
Implements the Mul trait for two MatZq values.
Mul is implemented for any combination of MatZq and borrowed MatZq.
Parameters:
other: specifies the value to multiply withself
Returns the product of self and other as a MatZq.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a = MatZq::from_str("[[2, 1],[1, 2]] mod 3").unwrap();
let b = MatZq::from_str("[[1, 0],[0, 1]] mod 3").unwrap();
let c = &a * &b;
let d = a * b;
let e = &c * d;
let f = c * &e;§Panics …
- if the dimensions of
selfandotherdo not match for multiplication. - if the moduli mismatch.
Source§impl MulAssign<&Z> for MatZq
impl MulAssign<&Z> for MatZq
Source§fn mul_assign(&mut self, scalar: &Z)
fn mul_assign(&mut self, scalar: &Z)
Computes the scalar multiplication of self and scalar reusing
the memory of self.
Parameters:
scalar: specifies the value to multiply toself
Returns the scalar of the matrix as a MatZq.
§Examples
use qfall_math::integer::Z;
use qfall_math::integer_mod_q::{MatZq, Zq};
use std::str::FromStr;
let mut a = MatZq::from_str("[[2, 1],[1, 2]] mod 61").unwrap();
let b = Z::from(2);
let c = Zq::from((17, 61));
a *= &b;
a *= b;
a *= 2;
a *= -2;
a *= &c;
a *= c;Source§impl MulAssign<Z> for MatZq
impl MulAssign<Z> for MatZq
Source§fn mul_assign(&mut self, other: Z)
fn mul_assign(&mut self, other: Z)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<Zq> for MatZq
impl MulAssign<Zq> for MatZq
Source§fn mul_assign(&mut self, other: Zq)
fn mul_assign(&mut self, other: Zq)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<i16> for MatZq
impl MulAssign<i16> for MatZq
Source§fn mul_assign(&mut self, other: i16)
fn mul_assign(&mut self, other: i16)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<i32> for MatZq
impl MulAssign<i32> for MatZq
Source§fn mul_assign(&mut self, other: i32)
fn mul_assign(&mut self, other: i32)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<i64> for MatZq
impl MulAssign<i64> for MatZq
Source§fn mul_assign(&mut self, other: i64)
fn mul_assign(&mut self, other: i64)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<i8> for MatZq
impl MulAssign<i8> for MatZq
Source§fn mul_assign(&mut self, other: i8)
fn mul_assign(&mut self, other: i8)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<u16> for MatZq
impl MulAssign<u16> for MatZq
Source§fn mul_assign(&mut self, other: u16)
fn mul_assign(&mut self, other: u16)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<u32> for MatZq
impl MulAssign<u32> for MatZq
Source§fn mul_assign(&mut self, other: u32)
fn mul_assign(&mut self, other: u32)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<u64> for MatZq
impl MulAssign<u64> for MatZq
Source§fn mul_assign(&mut self, other: u64)
fn mul_assign(&mut self, other: u64)
Documentation at MatZq::mul_assign.
Source§impl MulAssign<u8> for MatZq
impl MulAssign<u8> for MatZq
Source§fn mul_assign(&mut self, other: u8)
fn mul_assign(&mut self, other: u8)
Documentation at MatZq::mul_assign.
Source§impl PartialEq for MatZq
impl PartialEq for MatZq
Source§fn eq(&self, other: &Self) -> bool
fn eq(&self, other: &Self) -> bool
Checks if two MatZq instances are equal. Used by the == and != operators.
The values in the matrix as well as the modulus have to be equal.
Parameters:
other: the other value that is compare againstself
Returns true if the elements are equal, otherwise false.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a = MatZq::from_str("[[1, 2],[3, 4]] mod 4").unwrap();
let b = MatZq::from_str("[[1, 2],[2, 4]] mod 4").unwrap();
// These are all equivalent and return false.
let compared: bool = (a == b);
let compared: bool = (&a == &b);
let compared: bool = (a.eq(&b));
let compared: bool = (MatZq::eq(&a, &b));Source§impl Sub<&MatZ> for &MatZq
impl Sub<&MatZ> for &MatZq
Source§fn sub(self, other: &MatZ) -> Self::Output
fn sub(self, other: &MatZ) -> Self::Output
Implements the Sub trait for a MatZq and a MatZ matrix.
Sub is implemented for any combination of owned and borrowed values.
Parameters:
other: specifies the matrix to subtract fromself.
Returns the subtraction of self and other as a MatZq.
§Examples
use qfall_math::{integer::MatZ, integer_mod_q::MatZq};
use std::str::FromStr;
let a = MatZ::from_str("[[1, 2, 3],[3, 4, 5]]").unwrap();
let b = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c = &b - &a;
let d = b.clone() - a.clone();
let e = &b - &a;
let f = b - a;§Panics …
- if the dimensions of both matrices mismatch.
Source§impl Sub<&MatZq> for &MatZ
impl Sub<&MatZq> for &MatZ
Source§fn sub(self, other: &MatZq) -> Self::Output
fn sub(self, other: &MatZq) -> Self::Output
Implements the Sub trait for a MatZ and a MatZq matrix.
Sub is implemented for any combination of owned and borrowed values.
Parameters:
other: specifies the matrix to subtract fromself.
Returns the subtraction of self and other as a MatZq.
§Examples
use qfall_math::{integer::MatZ, integer_mod_q::MatZq};
use std::str::FromStr;
let a = MatZ::from_str("[[1, 2, 3],[3, 4, 5]]").unwrap();
let b = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c = &a - &b;
let d = a.clone() - b.clone();
let e = &a - &b;
let f = a - b;§Panics …
- if the dimensions of both matrices mismatch.
Source§impl Sub for &MatZq
impl Sub for &MatZq
Source§fn sub(self, other: Self) -> Self::Output
fn sub(self, other: Self) -> Self::Output
Implements the Sub trait for two MatZq values.
Sub is implemented for any combination of MatZq and borrowed MatZq.
Parameters:
other: specifies the value to subtract fromself
Returns the result of the subtraction as a MatZq.
§Examples
use qfall_math::integer_mod_q::MatZq;
use std::str::FromStr;
let a: MatZq = MatZq::from_str("[[1, 2, 3],[3, 4, 5]] mod 7").unwrap();
let b: MatZq = MatZq::from_str("[[1, 9, 3],[1, 0, 5]] mod 7").unwrap();
let c: MatZq = &a - &b;
let d: MatZq = a - b;
let e: MatZq = &c - d;
let f: MatZq = c - &e;§Panics …
- if the dimensions of both matrices mismatch.
- if the moduli mismatch.
Source§impl SubAssign<&MatZ> for MatZq
impl SubAssign<&MatZ> for MatZq
Source§fn sub_assign(&mut self, other: &MatZ)
fn sub_assign(&mut self, other: &MatZ)
Documentation at MatZq::sub_assign.
Source§impl SubAssign<&MatZq> for MatZq
impl SubAssign<&MatZq> for MatZq
Source§fn sub_assign(&mut self, other: &Self)
fn sub_assign(&mut self, other: &Self)
Computes the subtraction of self and other reusing
the memory of self.
SubAssign can be used on MatZq in combination with
MatZq and MatZ.
Parameters:
other: specifies the value to subtract fromself
§Examples
use qfall_math::{integer_mod_q::MatZq, integer::MatZ};
let mut a = MatZq::identity(2, 2, 7);
let b = MatZq::new(2, 2, 7);
let c = MatZ::new(2, 2);
a -= &b;
a -= b;
a -= &c;
a -= c;§Panics …
- if the matrix dimensions mismatch.
- if the moduli of the matrices mismatch.
Source§impl SubAssign<MatZ> for MatZq
impl SubAssign<MatZ> for MatZq
Source§fn sub_assign(&mut self, other: MatZ)
fn sub_assign(&mut self, other: MatZ)
Documentation at MatZq::sub_assign.
Source§impl SubAssign for MatZq
impl SubAssign for MatZq
Source§fn sub_assign(&mut self, other: MatZq)
fn sub_assign(&mut self, other: MatZq)
Documentation at MatZq::sub_assign.
Source§impl Tensor for MatZq
impl Tensor for MatZq
Source§fn tensor_product(&self, other: &Self) -> Self
fn tensor_product(&self, other: &Self) -> Self
Computes the tensor product of self with other.
Parameters:
other: the value with which the tensor product is computed.
Returns the tensor product of self with other and panics if the
moduli of the provided matrices mismatch.
§Examples
use qfall_math::integer_mod_q::MatZq;
use qfall_math::traits::Tensor;
use std::str::FromStr;
let mat_1 = MatZq::from_str("[[1, 1],[2, 2]] mod 7").unwrap();
let mat_2 = MatZq::from_str("[[1, 2],[3, 4]] mod 7").unwrap();
let mat_ab = mat_1.tensor_product(&mat_2);
let mat_ba = mat_2.tensor_product(&mat_1);
let res_ab = "[[1, 2, 1, 2],[3, 4, 3, 4],[2, 4, 2, 4],[6, 1, 6, 1]] mod 7";
let res_ba = "[[1, 1, 2, 2],[2, 2, 4, 4],[3, 3, 4, 4],[6, 6, 1, 1]] mod 7";
assert_eq!(mat_ab, MatZq::from_str(res_ab).unwrap());
assert_eq!(mat_ba, MatZq::from_str(res_ba).unwrap());§Panics …
- if the moduli of both matrices mismatch.
Use
tensor_product_safeto get an error instead.