Template Class DCRTPolyInterface

Nested Relationships

Nested Types

Inheritance Relationships

Base Type

Derived Type

Template Parameter Order

  1. typename DerivedType

  2. typename BigVecType

  3. typename LilVecType

  4. template< typename LVT > typename RNSContainerType

Class Documentation

template<typename DerivedType, typename BigVecType, typename LilVecType, template<typename LVT> typename RNSContainerType>
class lbcrypto::DCRTPolyInterface : public lbcrypto::ILElement<DerivedType, BigVecType>

Inheritence diagram for lbcrypto::DCRTPolyInterface:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="lbcrypto::ILElement< DerivedType, BigVecType >" tooltip="lbcrypto::ILElement< DerivedType, BigVecType >"] "3" [label="lbcrypto::Serializable" tooltip="lbcrypto::Serializable"] "4" [label="lbcrypto::DCRTPolyImpl< BigVector >" tooltip="lbcrypto::DCRTPolyImpl< BigVector >"] "1" [label="lbcrypto::DCRTPolyInterface< DerivedType, BigVecType, LilVecType, RNSContainerType >" tooltip="lbcrypto::DCRTPolyInterface< DerivedType, BigVecType, LilVecType, RNSContainerType >" fillcolor="#BFBFBF"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for lbcrypto::DCRTPolyInterface:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="lbcrypto::ILElement< DerivedType, BigVecType >" tooltip="lbcrypto::ILElement< DerivedType, BigVecType >"] "3" [label="lbcrypto::Serializable" tooltip="lbcrypto::Serializable"] "1" [label="lbcrypto::DCRTPolyInterface< DerivedType, BigVecType, LilVecType, RNSContainerType >" tooltip="lbcrypto::DCRTPolyInterface< DerivedType, BigVecType, LilVecType, RNSContainerType >" fillcolor="#BFBFBF"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Ideal lattice for the double-CRT interface representation. The interface contains a methods required for computations on lattices The double-CRT representation of polynomials is a common optimization for lattice encryption operations. Basically, it allows large-modulus polynomials to be represented as multiple smaller-modulus polynomials. The double-CRT representations are discussed theoretically here:

  • Gentry C., Halevi S., Smart N.P. (2012) Homomorphic Evaluation of the AES Circuit. In: Safavi-Naini R., Canetti R. (eds) Advances in Cryptology – CRYPTO 2012. Lecture Notes in Computer Science, vol 7417. Springer, Berlin, Heidelberg

example for the default DerivedType the template types would be… DerivedType - DCRTPolyImpl<BigVector> BigVecType - BigVector LilVecType - NativeVector RNSContainer<LVT> - PolyImpl

Template Parameters
  • DerivedType – Curiously-Recurring-Template-Pattern

  • BigVecType – The Vector type before decomposing the polynomial into CRT

  • LilVecType – The underlaying RNS data structure, a vectors type structure, that will compose the CRT data

  • RNSContainer – The container of LilVecType, a lbcrypto::PolyImpl or vector typically

Subclassed by lbcrypto::DCRTPolyImpl< BigVector >

Public Types

using BigIntType = typename BigVecType::Integer
using Params = ILDCRTParams<BigIntType>
using LilIntType = typename LilVecType::Integer
using TowerType = RNSContainerType<LilVecType>
using PolyLargeType = RNSContainerType<BigVecType>
using DggType = DiscreteGaussianGeneratorImpl<LilVecType>
using DugType = DiscreteUniformGeneratorImpl<LilVecType>
using TugType = TernaryUniformGeneratorImpl<LilVecType>
using BugType = BinaryUniformGeneratorImpl<LilVecType>
typedef struct CRTBasisExtensionPrecomputations CRTBasisExtensionPrecomputations

Public Functions

inline DerivedType &GetDerived()

Get the Derived object, this is apart of the CRTP software design pattern it allows the base class (this one) to implement methods that call the derived objects implementation.

Chapter 21.2 “C++ Templates The Complete Guide” by David Vandevoorde and Nicolai M. Josuttis http://www.informit.com/articles/article.asp?p=31473

Returns

DerivedType&

inline const DerivedType &GetDerived() const
inline DerivedType CloneTowers(uint32_t startTower, uint32_t endTower)

Makes a copy of the DCRTPoly, but it includes only a sequential subset of the towers that the original holds.

Parameters
  • startTower – The index number of the first tower to clone

  • endTower – The index number of the last tower to clone

Returns

new Element

inline virtual DerivedType Clone() const final

Clone the object by making a copy of it and returning the copy.

Returns

new Element

inline virtual DerivedType CloneEmpty() const final

Clone the object, but have it contain nothing.

Returns

new Element

inline virtual DerivedType CloneParametersOnly() const final

Clones the element’s parameters, leaves vector initialized to 0.

Returns

new Element

virtual DerivedType CloneWithNoise(const DiscreteGaussianGeneratorImpl<BigVecType> &dgg, Format format) const override = 0

Clone with noise. This method creates a new DCRTPoly and clones the params. The tower values will be filled up with noise based on the discrete gaussian.

Parameters
  • &dgg – the input discrete Gaussian generator. The dgg will be the seed to populate the towers of the DCRTPoly with random numbers.

  • format – the input format fixed to EVALUATION. Format is a enum type that indicates if the polynomial is in Evaluation representation or Coefficient representation. It is defined in inttypes.h.

virtual Format GetFormat() const override = 0

Get method of the format.

Returns

the format, either COEFFICIENT or EVALUATION

inline const std::shared_ptr<Params> &GetParams() const

returns the parameters of the element.

Returns

the element parameter set.

inline virtual usint GetCyclotomicOrder() const final

returns the element’s cyclotomic order

Returns

returns the cyclotomic order of the element.

inline usint GetRingDimension() const

returns the element’s ring dimension

Returns

returns the ring dimension of the element.

inline virtual const BigIntType &GetModulus() const final

returns the element’s modulus

Returns

returns the modulus of the element.

inline const BigIntType &GetOriginalModulus() const

returns the element’s original modulus, derived from Poly

Returns

returns the modulus of the element.

inline const BigIntType GetRootOfUnity() const

returns the element’s root of unity.

Returns

the element’s root of unity.

inline virtual usint GetLength() const final

Get method for length of each component element. NOTE assumes all components are the same size. (Ring Dimension)

Returns

length of the component element

inline virtual BigIntType &at(usint i) final

Get interpolated value of elements at all tower index i. Note this operation is computationally intense. Does bound checking.

Returns

interpolated value at index i.

inline virtual const BigIntType &at(usint i) const final
inline virtual BigIntType &operator[](usint i) final

Get interpolated value of element at index i. Note this operation is computationally intense. No bound checking.

Returns

interpolated value at index i.

inline virtual const BigIntType &operator[](usint i) const final
inline const std::vector<TowerType> &GetAllElements() const

Get method that returns a vector of all component elements.

Returns

a vector of the component elements.

inline std::vector<TowerType> &GetAllElements()
inline usint GetNumOfElements() const

Get method of the number of component elements, also known as the number of towers.

Returns

the number of component elements.

inline const TowerType &GetElementAtIndex(usint i) const

Get method of individual tower of elements. Note this behavior is different than poly.

Parameters

i – index of tower to be returned.

Returns

a reference to the returned tower

inline void SetElementAtIndex(usint index, const TowerType &element)

Sets element at index.

Parameters
  • index – where the element should be set

  • element – The element to store

inline void SetElementAtIndex(usint index, TowerType &&element)

Sets element at index.

Parameters
  • index – where the element should be set

  • element – The element to store

virtual std::vector<DerivedType> BaseDecompose(usint baseBits, bool evalModeAnswer) const override = 0

Write the element as \sum\limits{i=0}^{\lfloor {\log q/base} \rfloor} {(base^i u_i)} and return the vector of \left\{u_0, u_1,...,u_{\lfloor {\log q/base} \rfloor} \right\} \in R_{{base}^{\lceil {\log q/base} \rceil}}; This is used as a subroutine in the relinearization procedure.

Warning

not efficient and not fast, uses multiprecision arithmetic and will be removed in future. Use

Parameters

baseBits – is the number of bits in the base, i.e., base = 2^{baseBits}.

Returns

is the pointer where the base decomposition vector is stored

virtual std::vector<DerivedType> PowersOfBase(usint baseBits) const override = 0

Generate a vector of PolyImpl’s as \left\{x, {base}*x, {base}^2*x, ..., {base}^{\lfloor {\log q/{base}} \rfloor} \right\}*x, where x is the current PolyImpl object; used as a subroutine in the relinearization procedure to get powers of a certain “base” for the secret key element.

Warning

not efficient and not fast, uses multiprecision arithmetic and will be removed in future. Use

Parameters

baseBits – is the number of bits in the base, i.e., base = 2^{baseBits}.

Returns

is the pointer where the base decomposition vector is stored

inline std::vector<DerivedType> CRTDecompose(uint32_t baseBits) const

CRT basis decomposition of c as [c qi/q]_qi

Parameters

&baseBits – bits in the base for additional digit decomposition if base > 0

Returns

is the pointer where the resulting vector is stored

inline DerivedType &operator=(const TowerType &rhs)
virtual DerivedType &operator=(const DerivedType &rhs) override = 0

Assignment Operator.

Parameters

&rhs – the copied element.

Returns

the resulting element.

virtual DerivedType &operator=(DerivedType &&rhs) override = 0

Move Assignment Operator.

Parameters

&rhs – the copied element.

Returns

the resulting element.

virtual DerivedType &operator=(std::initializer_list<uint64_t> rhs) override = 0

Initalizer list.

Parameters

&rhs – the list to initalized the element.

Returns

the resulting element.

inline DerivedType &operator=(uint64_t val)

Assignment Operator. The usint val will be set at index zero and all other indices will be set to zero.

Parameters

val – is the usint to assign to index zero.

Returns

the resulting vector.

inline DerivedType &operator=(const std::vector<int64_t> &rhs)

Creates a Poly from a vector of signed integers (used for trapdoor sampling)

Parameters

&rhs – the vector to set the PolyImpl to.

Returns

the resulting PolyImpl.

inline DerivedType &operator=(const std::vector<int32_t> &rhs)

Creates a Poly from a vector of signed integers (used for trapdoor sampling)

Parameters

&rhs – the vector to set the PolyImpl to.

Returns

the resulting PolyImpl.

inline DerivedType &operator=(std::initializer_list<std::string> rhs)

Initalizer list.

Parameters

&rhs – the list to set the PolyImpl to.

Returns

the resulting PolyImpl.

virtual DerivedType operator-() const override = 0

Unary minus on a element.

Returns

additive inverse of the an element.

virtual bool operator==(const DerivedType &rhs) const override = 0

Equality operator.

Parameters

&rhs – is the specified element to be compared with this element.

Returns

true if this element represents the same values as the specified element, false otherwise.

virtual DerivedType &operator+=(const DerivedType &rhs) override = 0

Performs an entry-wise addition over all elements of each tower with the towers of the element on the right hand side.

Parameters

&rhs – is the element to add with.

Returns

is the result of the addition.

virtual DerivedType &operator-=(const DerivedType &rhs) override = 0

Performs an entry-wise subtraction over all elements of each tower with the towers of the element on the right hand side.

Parameters

&rhs – is the element to subtract from.

Returns

is the result of the addition.

virtual DerivedType AutomorphismTransform(uint32_t i) const override = 0

Permutes coefficients in a polynomial. Moves the ith index to the first one, it only supports odd indices.

Parameters

&i – is the element to perform the automorphism transform with.

Returns

is the result of the automorphism transform.

virtual DerivedType AutomorphismTransform(uint32_t i, const std::vector<uint32_t> &vec) const override = 0

Performs an automorphism transform operation using precomputed bit reversal indices.

Parameters
  • &i – is the element to perform the automorphism transform with.

  • &vec – a vector with precomputed indices

Returns

is the result of the automorphism transform.

inline virtual DerivedType Transpose() const final

Transpose the ring element using the automorphism operation.

Returns

is the result of the transposition.

virtual DerivedType Plus(const DerivedType &rhs) const override = 0

Performs an addition operation and returns the result.

Parameters

&element – is the element to add with.

Returns

is the result of the addition.

virtual DerivedType Times(const DerivedType &rhs) const override = 0

Performs a multiplication operation and returns the result.

Parameters

&element – is the element to multiply with.

Returns

is the result of the multiplication.

virtual DerivedType Minus(const DerivedType &rhs) const override = 0

Performs a subtraction operation and returns the result.

Parameters

&element – is the element to subtract from.

Returns

is the result of the subtraction.

virtual DerivedType Plus(const BigIntType &rhs) const override = 0

Scalar addition - add an element to the first index of each tower.

Parameters

&element – is the element to add entry-wise.

Returns

is the result of the addition operation.

inline DerivedType Plus(const std::vector<BigIntType> &rhs) const

Scalar addition for elements in CRT format. CRT elements are represented as vector of integer elements which correspond to the represented number modulo the primes in the tower chain (in same order).

Parameters

&element – is the element to add entry-wise.

Returns

is the result of the addition operation.

virtual DerivedType Minus(const BigIntType &rhs) const override = 0

Scalar subtraction - subtract an element to all entries.

Parameters

&element – is the element to subtract entry-wise.

Returns

is the return value of the minus operation.

inline DerivedType Minus(const std::vector<BigIntType> &rhs) const

Scalar subtraction for elements in CRT format. CRT elements are represented as vector of integer elements which correspond to the represented number modulo the primes in the tower chain (in same order).

Parameters

&element – is the element to subtract entry-wise.

Returns

is the result of the subtraction operation.

virtual DerivedType Times(const BigIntType &rhs) const override = 0

Scalar multiplication - multiply all entries.

Parameters

&element – is the element to multiply entry-wise.

Returns

is the return value of the times operation.

virtual DerivedType Times(NativeInteger::SignedNativeInt rhs) const override = 0

Scalar multiplication - multiply by a signed integer.

Parameters

&element – is the element to multiply entry-wise.

Returns

is the return value of the times operation.

inline DerivedType Times(int64_t rhs) const

Scalar multiplication - multiply by a signed integer.

Note

this is need for 128-bit so that the 64-bit inputs can be used.

Parameters

&element – is the element to multiply entry-wise.

Returns

is the return value of the times operation.

inline DerivedType Times(const std::vector<NativeInteger> &rhs) const

Scalar multiplication by an integer represented in CRT Basis.

Parameters

&element – is the element to multiply entry-wise.

Returns

is the return value of the times operation.

inline DerivedType TimesNoCheck(const std::vector<NativeInteger> &rhs) const

Performs a multiplication operation even when the multiplicands have a different number of towers.

Parameters

&element – is the element to multiply with.

Returns

is the result of the multiplication.

inline DerivedType Times(const std::vector<BigIntType> &rhs) const

Scalar modular multiplication by an integer represented in CRT Basis.

Warning

Should remove this, data is truncated to native-word size.

Parameters

&element – is the element to multiply entry-wise.

Returns

is the return value of the times operation.

inline DerivedType MultiplyAndRound(const BigIntType &p, const BigIntType &q) const final

Scalar multiplication followed by division and rounding operation - operation on all entries.

Warning

Will remove, this is only inplace because of BFV

Parameters
  • &p – is the element to multiply entry-wise.

  • &q – is the element to divide entry-wise.

Returns

is the return value of the multiply, divide and followed by rounding operation.

inline DerivedType DivideAndRound(const BigIntType &q) const final

Scalar division followed by rounding operation - operation on all entries.

Warning

Will remove, this is only inplace because of BFV

Parameters

&q – is the element to divide entry-wise.

Returns

is the return value of the divide, followed by rounding operation.

virtual DerivedType Negate() const = 0

Performs a negation operation and returns the result.

Returns

is the result of the negation.

virtual DerivedType &operator+=(const BigIntType &rhs) override = 0
virtual DerivedType &operator+=(const LilIntType &rhs) = 0
virtual DerivedType &operator-=(const BigIntType &rhs) override = 0

Performs a subtraction operation and returns the result.

Parameters

&element – is the element to subtract from.

Returns

is the result of the subtraction.

virtual DerivedType &operator-=(const LilIntType &rhs) = 0
virtual DerivedType &operator*=(const BigIntType &rhs) override = 0

Performs a multiplication operation and returns the result.

Parameters

&element – is the element to multiply by.

Returns

is the result of the multiplication.

virtual DerivedType &operator*=(const LilIntType &rhs) = 0
virtual DerivedType &operator*=(const DerivedType &rhs) override = 0

Performs a multiplication operation and returns the result.

Parameters

&element – is the element to multiply with.

Returns

is the result of the multiplication.

virtual DerivedType MultiplicativeInverse() const override = 0

Performs a multiplicative inverse operation and returns the result.

Returns

is the result of the multiplicative inverse.

inline virtual DerivedType ModByTwo() const final

Perform a modulus by 2 operation. Returns the least significant bit.

Warning

Doesn’t make sense for DCRT

Returns

is the resulting value.

inline DerivedType Mod(const BigIntType &modulus) const final

Modulus - perform a modulus operation. Does proper mapping of [-modulus/2, modulus/2) to [0, modulus)

Warning

Doesn’t make sense for DCRT

Parameters

modulus – is the modulus to use.

Returns

is the return value of the modulus.

inline virtual const BigVecType &GetValues() const final

Get method that should not be used.

Warning

Doesn’t make sense for DCRT

Returns

will throw an error.

inline void SetValues(const BigVecType &values, Format format)

Set method that should not be used, will throw an error.

Warning

Doesn’t make sense for DCRT

Parameters
  • &values

  • format

virtual void SetValuesToZero() = 0

Sets all values of element to zero.

virtual void SetValuesModSwitch(const DerivedType &element, const NativeInteger &modulus) = 0

Sets values with a different modulus.

virtual void AddILElementOne() override = 0

Adds “1” to every entry in every tower.

inline DerivedType AddRandomNoise(const BigIntType &modulus) const

Add uniformly random values to all components except for the first one.

Warning

Doesn’t make sense for DCRT

inline virtual void MakeSparse(uint32_t wFactor) final

Make DCRTPoly Sparse. Sets every index of each tower not equal to zero mod the wFactor to zero.

Warning

Only used by RingSwitching, which is no longer supported. Will be removed in future.

Parameters

&wFactor – ratio between the sparse and none-sparse values.

virtual bool IsEmpty() const override = 0

Returns true if ALL the tower(s) are empty.

Returns

true if all towers are empty

virtual void DropLastElement() = 0

Drops the last element in the double-CRT representation. The resulting DCRTPoly element will have one less tower.

virtual void DropLastElements(size_t i) = 0

Drops the last i elements in the double-CRT representation.

virtual void DropLastElementAndScale(const std::vector<NativeInteger> &QlQlInvModqlDivqlModq, const std::vector<NativeInteger> &qlInvModq) = 0

Drops the last element in the double-CRT representation and scales down by the last CRT modulus. The resulting DCRTPoly element will have one less tower.

Parameters
  • &QlQlInvModqlDivqlModq – precomputed values for [Q^(l)*[Q^(l)^{-1}]_{q_l}/q_l]_{q_i}

  • &QlQlInvModqlDivqlModqPrecon – NTL-specific precomputations

  • &qlInvModq – precomputed values for [q_l^{-1}]_{q_i}

  • &qlInvModqPrecon – NTL-specific precomputations

virtual void ModReduce(const NativeInteger &t, const std::vector<NativeInteger> &tModqPrecon, const NativeInteger &negtInvModq, const NativeInteger &negtInvModqPrecon, const std::vector<NativeInteger> &qlInvModq, const std::vector<NativeInteger> &qlInvModqPrecon) = 0

ModReduces reduces the DCRTPoly element’s composite modulus by dropping the last modulus from the chain of moduli as well as dropping the last tower.

Parameters
  • &t – is the plaintextModulus used for the DCRTPoly

  • &tModqPrecon – NTL-specific precomputations for [t]_{q_i}

  • &negtInvModq – precomputed values for [-t^{-1}]_{q_i}

  • &negtInvModqPrecon – NTL-specific precomputations for [-t^{-1}]_{q_i}

  • &qlInvModq – precomputed values for [q_{l}^{-1}]_{q_i}

  • &qlInvModqPrecon – NTL-specific precomputations for [q_{l}^{-1}]_{q_i}

virtual PolyLargeType CRTInterpolate() const = 0

Interpolates the DCRTPoly to a Poly based on the Chinese Remainder Transform Interpolation. and then returns a Poly with that single element.

Returns

the interpolated ring element as a Poly object.

virtual TowerType DecryptionCRTInterpolate(PlaintextModulus ptm) const = 0
virtual TowerType ToNativePoly() const = 0

If the values are small enough this is used for efficiency.

Warning

This will be replaced with a non-member utility function.

Returns

NativePoly

virtual PolyLargeType CRTInterpolateIndex(usint i) const = 0

Interpolates the DCRTPoly to a Poly based on the Chinese Remainder Transform Interpolation, only at element index i, all other elements are zero. and then returns a Poly with that single element.

Returns

the interpolated ring element as a Poly object.

virtual BigIntType GetWorkingModulus() const = 0

Computes and returns the product of primes in the current moduli chain. Compared to GetModulus, which always returns the product of all primes in the crypto parameters, this method will return a different modulus, based on the towers/moduli that are currently in the chain (some towers are dropped along the way).

Returns

the product of moduli in the current towers.

virtual std::shared_ptr<Params> GetExtendedCRTBasis(const std::shared_ptr<Params> &paramsP) const = 0

Returns the element parameters for DCRTPoly elements in an extended CRT basis, which is the concatenation of the towers currently in “this” DCRTPoly, and the moduli in ParamsP.

Returns

element parameters of the extended basis.

virtual void TimesQovert(const std::shared_ptr<Params> &paramsQ, const std::vector<NativeInteger> &tInvModq, const NativeInteger &t, const NativeInteger &NegQModt, const NativeInteger &NegQModtPrecon) = 0
virtual DerivedType ApproxSwitchCRTBasis(const std::shared_ptr<Params> &paramsQ, const std::shared_ptr<Params> &paramsP, const std::vector<NativeInteger> &QHatInvModq, const std::vector<NativeInteger> &QHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModp, const std::vector<DoubleNativeInt> &modpBarrettMu) const = 0

Performs approximate CRT basis switching: {X}_{Q} -> {X’}_{P} X’ = X + alpha*Q for small alpha {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: [X’]_{p_j} = [\sum_i([x_i*(Q/q_i)^{-1}]_{q_i}*(Q/q_i)]_{p_j}

Source: “A full RNS variant of approximate homomorphic encryption” by Cheon, et. al.

Parameters
  • &paramsQ – parameters for the CRT basis {q_1,…,q_l}

  • &paramsP – parameters for the CRT basis {p_1,…,p_k}

  • &QHatinvModq – precomputed values for [(Q/q_i)^{-1}]_{q_i}

  • &QHatinvModqPrecon – NTL-specific precomputations

  • &QHatModp – precomputed values for [Q/q_i]_{p_j}

  • &modpBarrettMu – 128-bit Barrett reduction precomputed values

Returns

the representation of {X + alpha*Q} in basis {P}.

virtual void ApproxModUp(const std::shared_ptr<Params> &paramsQ, const std::shared_ptr<Params> &paramsP, const std::shared_ptr<Params> &paramsQP, const std::vector<NativeInteger> &QHatInvModq, const std::vector<NativeInteger> &QHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModp, const std::vector<DoubleNativeInt> &modpBarrettMu) = 0

Performs approximate modulus raising: {X}_{Q} -> {X’}_{Q,P}. X’ = X + alpha*Q for small alpha {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: {X}_{Q} -> {X’}_Q : trivial {X}_{Q} -> {X’}_P : use DCRTPoly::ApproxSwitchCRTBasis

Source: “A full RNS variant of approximate homomorphic encryption” by Cheon, et. al.

Parameters
  • &paramsQ – parameters for the CRT basis {q_1,…,q_l}

  • &paramsP – parameters for the CRT basis {p_1,…,p_k}

  • &QHatInvModq – precomputed values for [(Q/q_i)^{-1}]_{q_i}

  • &QHatInvModqPrecon – NTL-specific precomputations

  • &QHatModp – precomputed values for [Q/q_i]_{p_j}

  • &modpBarrettMu – 128-bit Barrett reduction precomputed values for p_j

Returns

the representation of {X + alpha*Q} in basis {Q,P}.

virtual DerivedType ApproxModDown(const std::shared_ptr<Params> &paramsQ, const std::shared_ptr<Params> &paramsP, const std::vector<NativeInteger> &PInvModq, const std::vector<NativeInteger> &PInvModqPrecon, const std::vector<NativeInteger> &PHatInvModp, const std::vector<NativeInteger> &PHatInvModpPrecon, const std::vector<std::vector<NativeInteger>> &PHatModq, const std::vector<DoubleNativeInt> &modqBarrettMu, const std::vector<NativeInteger> &tInvModp, const std::vector<NativeInteger> &tInvModpPrecon, const NativeInteger &t, const std::vector<NativeInteger> &tModqPrecon) const = 0

Performs approximate modulus reduction: {X}_{Q,P} -> {\approx(X/P)}_{Q}. {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: 1) use DCRTPoly::ApproxSwitchCRTBasis : {X}_{P} -> {X’}_{Q} 2) compute : {(X-X’) * P^{-1}}_{Q}

Source: “A full RNS variant of approximate homomorphic encryption” by Cheon, et. al.

Parameters
  • &paramsQ – parameters for the CRT basis {q_1,…,q_l}

  • &paramsP – parameters for the CRT basis {p_1,…,p_k}

  • &PInvModq – precomputed values for (P^{-1} mod q_j)

  • &PInvModqPrecon – NTL-specific precomputations

  • &PHatInvModp – precomputed values for [(P/p_j)^{-1}]_{p_j}

  • &PHatInvModpPrecon – NTL-specific precomputations

  • &PHatModq – precomputed values for [P/p_j]_{q_i}

  • &modqBarrettMu – 128-bit Barrett reduction precomputed values for q_i

  • &tInvModp – precomputed values for [t^{-1}]_{p_j} used in BGVrns

  • t – often corresponds to the plaintext modulus used in BGVrns

Returns

the representation of {\approx(X/P)}_{Q}

virtual DerivedType SwitchCRTBasis(const std::shared_ptr<Params> &paramsP, const std::vector<NativeInteger> &QHatInvModq, const std::vector<NativeInteger> &QHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModp, const std::vector<std::vector<NativeInteger>> &alphaQModp, const std::vector<DoubleNativeInt> &modpBarrettMu, const std::vector<double> &qInv) const = 0

Performs CRT basis switching: {X}_{Q} -> {X}_{P} {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: 1) X=\sum_i[x_i*(Q/q_i)^{-1}]_{q_i}*(Q/q_i)-alpha*Q 2) compute round[[x_i*(Q/q_i)^{-1}]_{q_i} / q_i] to find alpha 3) [X]_{p_j}=[\sum_i[x_i*(Q/q_i)^{-1}]_{q_i}*(Q/q_i)]_{p_j}-[alpha*Q]_{p_j}

Source: Halevi S., Polyakov Y., and Shoup V. An Improved RNS Variant of the BFV Homomorphic Encryption Scheme. Cryptology ePrint Archive, Report 2018/117. (https://eprint.iacr.org/2018/117)

Parameters
  • &paramsP – parameters for the CRT basis {p_1,…,p_k}

  • &QHatInvModq – precomputed values for [(Q/q_i)^{-1}]_{q_i}

  • &QHatInvModqPrecon – NTL-specific precomputations

  • &QHatModp – precomputed values for [Q/q_i]_{p_j}

  • &alphaQModp – precomputed values for [alpha*Q]_{p_j}

  • &modpBarrettMu – 128-bit Barrett reduction precomputed values for p_j @params &qInv precomputed values for 1/q_i

Returns

the representation of {X}_{P}

virtual void ExpandCRTBasis(const std::shared_ptr<Params> &paramsQP, const std::shared_ptr<Params> &paramsP, const std::vector<NativeInteger> &QHatInvModq, const std::vector<NativeInteger> &QHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModp, const std::vector<std::vector<NativeInteger>> &alphaQModp, const std::vector<DoubleNativeInt> &modpBarrettMu, const std::vector<double> &qInv, Format resultFormat) = 0

Performs modulus raising: {X}_{Q} -> {X}_{Q,P} {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: {X}_{Q} -> {X}_P : use DCRTPoly::SwitchCRTBasis combine {X}_{Q} and {X}_{P} Outputs the resulting polynomial in CRT/RNS

Source: Halevi S., Polyakov Y., and Shoup V. An Improved RNS Variant of the BFV Homomorphic Encryption Scheme. Cryptology ePrint Archive, Report 2018/117. (https://eprint.iacr.org/2018/117)

Parameters
  • &paramsQP – parameters for the CRT basis {q_1,…,q_l,p_1,…,p_k}

  • &paramsP – parameters for the CRT basis {p_1,…,p_k}

  • &QHatInvModq – precomputed values for [QInv_i]_{q_i}

  • &QHatInvModqPrecon – NTL-specific precomputations

  • &QHatModp – precomputed values for [QHat_i]_{p_j}

  • &alphaQModp – precomputed values for [alpha*Q]_{p_j}

  • &modpBarrettMu – 128-bit Barrett reduction precomputed values for p_j @params &qInv precomputed values for 1/q_i

  • resultFormat – Specifies the format we want the result to be in

virtual void ExpandCRTBasisReverseOrder(const std::shared_ptr<Params> &paramsQP, const std::shared_ptr<Params> &paramsP, const std::vector<NativeInteger> &QHatInvModq, const std::vector<NativeInteger> &QHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModp, const std::vector<std::vector<NativeInteger>> &alphaQModp, const std::vector<DoubleNativeInt> &modpBarrettMu, const std::vector<double> &qInv, Format resultFormat) = 0

Performs modulus raising in reverse order: {X}_{Q} -> {X}_{P,Q}.

virtual void FastExpandCRTBasisPloverQ(const CRTBasisExtensionPrecomputations &precomputed) = 0
virtual void ExpandCRTBasisQlHat(const std::shared_ptr<Params> &paramsQ, const std::vector<NativeInteger> &QlHatModq, const std::vector<NativeInteger> &QlHatModqPrecon, const usint sizeQ) = 0
virtual TowerType ScaleAndRound(const NativeInteger &t, const std::vector<NativeInteger> &tQHatInvModqDivqModt, const std::vector<NativeInteger> &tQHatInvModqDivqModtPrecon, const std::vector<NativeInteger> &tQHatInvModqBDivqModt, const std::vector<NativeInteger> &tQHatInvModqBDivqModtPrecon, const std::vector<double> &tQHatInvModqDivqFrac, const std::vector<double> &tQHatInvModqBDivqFrac) const = 0

Performs scale and round: {X}_{Q} -> {\round(t/Q*X)}_t {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: [\sum_i x_i*[t*QHatInv_i/q_i]_t + Round(\sum_i x_i*{t*QHatInv_i/q_i})]_t

Source: Halevi S., Polyakov Y., and Shoup V. An Improved RNS Variant of the BFV Homomorphic Encryption Scheme. Cryptology ePrint Archive, Report 2018/117. (https://eprint.iacr.org/2018/117)

Parameters
  • &t – often corresponds to the plaintext modulus

  • &tQHatInvModqDivqModt – precomputed values for [Floor{t*QHatInv_i/q_i}]_t

  • &tQHatInvModqDivqModtPrecon – NTL-specific precomputations

  • &tQHatInvModqBDivqModt – precomputed values for [Floor{t*QHatInv_i*B/q_i}]_t used when CRT moduli are 45..60 bits long

  • &tQHatInvBDivqModtPrecon – NTL-specific precomputations used when CRT moduli are 45..60 bits long

  • &tQHatInvModqDivqFrac – precomputed values for Frac{t*QHatInv_i/q_i}

  • &tQHatInvBDivqFrac – precomputed values for Frac{t*QHatInv_i*B/q_i} used when CRT moduli are 45..60 bits long

Returns

the result of computation as a polynomial with native 64-bit coefficients

virtual DerivedType ApproxScaleAndRound(const std::shared_ptr<Params> &paramsP, const std::vector<std::vector<NativeInteger>> &tPSHatInvModsDivsModp, const std::vector<DoubleNativeInt> &modpBarretMu) const = 0

Computes approximate scale and round: {X}_{Q,P} -> {\approx{t/Q * X}}_{P} {Q} = {q_1,…,q_l} {P} = {p_1,…,p_k}.

Brief algorithm: Let S = {Q,P} 1) [\sum_k x_k * alpha_k]_{p_j} 2) alpha_k = [Floor[t*P*[[SHatInv_k]_{s_k}/s_k]]_{p_j}

Source: Halevi S., Polyakov Y., and Shoup V. An Improved RNS Variant of the BFV Homomorphic Encryption Scheme. Cryptology ePrint Archive, Report 2018/117. (https://eprint.iacr.org/2018/117)

Parameters
  • &paramsP – parameters for the CRT basis {p_1,…,p_k}

  • &tPSHatInvModsDivsModp – precomputed values for [\floor[t*P*[[SHatInv_k]_{s_k}/s_k]]_{p_j}

  • &modpBarretMu – 128-bit Barrett reduction precomputed values for p_j

Returns

the result {\approx{t/Q * X}}_{P}

virtual DerivedType ScaleAndRound(const std::shared_ptr<Params> &paramsOutput, const std::vector<std::vector<NativeInteger>> &tOSHatInvModsDivsModo, const std::vector<double> &tOSHatInvModsDivsFrac, const std::vector<DoubleNativeInt> &modoBarretMu) const = 0

Computes scale and round: {X}_{I,O} -> {t/I * X}_{O} {I} = {i_1,…,i_l} {O} = {o_1,…,o_k} O, the output modulus can be either P or Q, and I is the other one.

Brief algorithm: Let S = {I,O} 1) [\sum_k x_k * alpha_k + Round(\sum_k beta_k * x_k)]_{o_j} 2) alpha_k = [Floor[t*O*[[SHatInv_k]_{s_k}/s_k]]_{o_j} 3) beta_k = {t*O*[[SHatInv_k]_{s_k}/s_k}

Source: Halevi S., Polyakov Y., and Shoup V. An Improved RNS Variant of the BFV Homomorphic Encryption Scheme. Cryptology ePrint Archive, Report 2018/117. (https://eprint.iacr.org/2018/117)

Parameters
  • &paramsOutput – parameters for the CRT basis {o_1,…,o_k}.

  • &tOSHatInvModsDivsModo – precomputed values for [\floor[t*O*[[SHatInv_k]_{s_k}/s_k]]_{o_j}

  • &tPSHatInvModsDivsFrac – precomputed values for {t*O*[[SHatInv_k]_{s_k}/s_k}

  • &modoBarretMu – 128-bit Barrett reduction precomputed values for o_j

Returns

the result {t/I * X}_{O}

virtual TowerType ScaleAndRound(const std::vector<NativeInteger> &moduliQ, const NativeInteger &t, const NativeInteger &tgamma, const std::vector<NativeInteger> &tgammaQHatModq, const std::vector<NativeInteger> &tgammaQHatModqPrecon, const std::vector<NativeInteger> &negInvqModtgamma, const std::vector<NativeInteger> &negInvqModtgammaPrecon) const = 0

Computes scale and round for fast rounding: {X}_{Q} -> {\round(t/Q * X)}_t {Q} = {q_1,…,q_l}.

Brief algorithm:

Source: Jean-Claude Bajard and Julien Eynard and Anwar Hasan and Vincent Zucca. A Full RNS Variant of FV like Somewhat Homomorphic Encryption Schemes. Cryptology ePrint Archive: Report 2016/510. (https://eprint.iacr.org/2016/510)

Parameters
  • &moduliQ – moduli {q_1,…,q_l}

  • &t – often corresponds to the plaintext modulus

  • &tgamma – t * gamma : t * 2^26 reduction

  • &tgammaQHatModq – [t*gamma*(Q/q_i)]_{q_i}

  • &tgammaQHatModqPrecon – NTL-specific precomputations

  • &negInvqModtgamma – [-q^{-1}]_{t*gamma}

  • &negInvqModtgammaPrecon – NTL-specific precomputations

Returns

virtual void ScaleAndRoundPOverQ(const std::shared_ptr<Params> &paramsQ, const std::vector<NativeInteger> &pInvModq) = 0

Computes scale and round for BFV encryption mode EXTENDED: {X}_{Qp} -> {\round(1/p * X)}_Q {Q} = {q_1,…,q_l}.

Source: Andrey Kim and Yuriy Polyakov and Vincent Zucca. Revisiting Homomorphic Encryption Schemes for Finite Fields. Cryptology ePrint Archive: Report 2021/204. (https://eprint.iacr.org/2021/204.pdf)

Parameters
  • &paramsQ – Parameters for moduli {q_1,…,q_l}

  • &pInvModq – p^{-1}_{q_i}

Returns

virtual void FastBaseConvqToBskMontgomery(const std::shared_ptr<Params> &paramsQBsk, const std::vector<NativeInteger> &moduliQ, const std::vector<NativeInteger> &moduliBsk, const std::vector<DoubleNativeInt> &modbskBarrettMu, const std::vector<NativeInteger> &mtildeQHatInvModq, const std::vector<NativeInteger> &mtildeQHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModbsk, const std::vector<uint64_t> &QHatModmtilde, const std::vector<NativeInteger> &QModbsk, const std::vector<NativeInteger> &QModbskPrecon, const uint64_t &negQInvModmtilde, const std::vector<NativeInteger> &mtildeInvModbsk, const std::vector<NativeInteger> &mtildeInvModbskPrecon) = 0

Expands basis: {X}_{Q} -> {X}_{Q,Bsk,mtilde} mtilde is a redundant modulus used to remove q overflows generated from fast conversion. Outputs the resulting polynomial in CRT/RNS {Q} = {q_1,…,q_l} {Bsk} = {bsk_1,…,bsk_k}.

Source: Jean-Claude Bajard and Julien Eynard and Anwar Hasan and Vincent Zucca. A Full RNS Variant of FV like Somewhat Homomorphic Encryption Schemes. Cryptology ePrint Archive: Report 2016/510. (https://eprint.iacr.org/2016/510)

Parameters
  • paramsQBsk – container of QBsk moduli and roots on unity

  • &moduliQ – basis {Q} = {q_1,q_2,…,q_l}

  • &moduliBsk – basis {Bsk U mtilde} …

  • &modbskBarrettMu – 128-bit Barrett reduction precomputed values for bsk_j

  • &mtildeQHatInvModq – [mtilde*(Q/q_i)^{-1}]_{q_i}

  • &mtildeQHatInvModqPrecon – NTL-specific precomputations

  • &QHatModbsk – [Q/q_i]_{bsk_j}

  • &QHatModmtilde – [Q/q_i]_{mtilde}

  • &QModbsk – [Q]_{bsk_j}

  • &QModbskPrecon – NTL-specific precomputations

  • &negQInvModmtilde – [-Q^{-1}]_{mtilde}

  • &mtildeInvModbsk – [mtilde^{-1}]_{bsk_j}

  • &mtildeInvModbskPrecon – NTL-specific precomputations

virtual void FastRNSFloorq(const NativeInteger &t, const std::vector<NativeInteger> &moduliQ, const std::vector<NativeInteger> &moduliBsk, const std::vector<DoubleNativeInt> &modbskBarrettMu, const std::vector<NativeInteger> &tQHatInvModq, const std::vector<NativeInteger> &tQHatInvModqPrecon, const std::vector<std::vector<NativeInteger>> &QHatModbsk, const std::vector<std::vector<NativeInteger>> &qInvModbsk, const std::vector<NativeInteger> &tQInvModbsk, const std::vector<NativeInteger> &tQInvModbskPrecon) = 0

Computes scale and floor: {X}_{Q,Bsk} -> {\floor{t/Q * X}}_{Bsk} {Q} = {q_1,…,q_l} {Bsk} = {bsk_1,…,bsk_k} Outputs the resulting polynomial in CRT/RNS.

Source: Jean-Claude Bajard and Julien Eynard and Anwar Hasan and Vincent Zucca. A Full RNS Variant of FV like Somewhat Homomorphic Encryption Schemes. Cryptology ePrint Archive: Report 2016/510. (https://eprint.iacr.org/2016/510)

Parameters
  • &t – plaintext modulus

  • &moduliQ – {Q} = {q_1,…,q_l}

  • &moduliBsk – {Bsk} = {bsk_1,…,bsk_k}

  • &modbskBarrettMu – 128-bit Barrett reduction precomputed values for bsk_j

  • &tQHatInvModq – [(Q/q_i)^{-1}]_{q_i}

  • &tQHatInvModqPrecon – NTL-specific precomputations

  • &QHatModbsk – [Q/q_i]_{bsk_i}

  • &qInvModbsk – [(q_i)^{-1}]_{bsk_j}

  • &tQInvModbsk – [t*Q^{-1}]_{bsk_j}

  • &tQInvModbskPrecon – NTL-specific precomputations

virtual void FastBaseConvSK(const std::shared_ptr<Params> &paramsQ, const std::vector<DoubleNativeInt> &modqBarrettMu, const std::vector<NativeInteger> &moduliBsk, const std::vector<DoubleNativeInt> &modbskBarrettMu, const std::vector<NativeInteger> &BHatInvModb, const std::vector<NativeInteger> &BHatInvModbPrecon, const std::vector<NativeInteger> &BHatModmsk, const NativeInteger &BInvModmsk, const NativeInteger &BInvModmskPrecon, const std::vector<std::vector<NativeInteger>> &BHatModq, const std::vector<NativeInteger> &BModq, const std::vector<NativeInteger> &BModqPrecon) = 0

Converts basis: {X}_{Q,Bsk} -> {X}_{Bsk} {Q} = {q_1,…,q_l} {Bsk} = {bsk_1,…,bsk_k} using Shenoy Kumaresan method. Outputs the resulting polynomial in CRT/RNS.

Source: Jean-Claude Bajard and Julien Eynard and Anwar Hasan and Vincent Zucca. A Full RNS Variant of FV like Somewhat Homomorphic Encryption Schemes. Cryptology ePrint Archive: Report 2016/510. (https://eprint.iacr.org/2016/510)

Note in the source paper, B is referred to by M.

Parameters
  • &paramsQParams for Q

  • &modqBarrettMu – precomputed Barrett Mu for q_i

  • &moduliBsk – basis {Bsk} = {bsk_1,…,bsk_k}

  • &modbskBarrettMu – precomputed Barrett Mu for bsk_j

  • &BHatInvModb – [(B/b_j)^{-1}]_{b_j}

  • &BHatInvModbPrecon – NTL precomptations for [(B/b_j)^{-1}]_{b_j}

  • &BHatModmsk – [B/b_j]_{msk}

  • &BInvModmsk – [B^{-1}]_{msk}

  • &BInvModmskPrecon – NTL precomptation for [B^{-1}]_{msk}

  • &BHatModq – [B/b_j]_{q_i}

  • &BModq – [B]_{q_i}

  • &BModqPrecon – NTL precomptations for [B]_{q_i}

virtual void SwitchFormat() override = 0

Convert from Coefficient to CRT or vice versa; calls FFT and inverse FFT.

See also

SetFormat(format) instead

Warning

use

virtual void OverrideFormat(const Format f) = 0

Sets format to value without performing NTT. Only use if you know what you’re doing.

inline void SwitchModulus(const BigIntType &modulus, const BigIntType &rootOfUnity, const BigIntType &modulusArb, const BigIntType &rootOfUnityArb) final

Switch modulus and adjust the values.

Parameters
  • &modulus – is the modulus to be set

  • &rootOfUnity – is the corresponding root of unity for the modulus

  • &modulusArb – is the modulus used for arbitrary cyclotomics CRT

  • &rootOfUnityArb – is the corresponding root of unity for the modulus ASSUMPTION: This method assumes that the caller provides the correct rootOfUnity for the modulus

virtual void SwitchModulusAtIndex(size_t index, const BigIntType &modulus, const BigIntType &rootOfUnity) = 0

Switch modulus at tower i and adjust the values.

Parameters
  • index – is the index for the tower

  • &modulus – is the modulus to be set

  • &rootOfUnity – is the corresponding root of unity for the modulus ASSUMPTION: This method assumes that the caller provides the correct rootOfUnity for the modulus

virtual bool InverseExists() const override = 0

Determines if inverse exists.

Returns

is the Boolean representation of the existence of multiplicative inverse.

inline virtual double Norm() const final

Returns the infinity norm, basically the largest value in the ring element.

Returns

is the largest value in the ring element.

inline const std::string GetElementName() const

Public Static Functions

static inline std::function<DerivedType()> Allocator(const std::shared_ptr<Params> &params, Format format)

Create lambda that allocates a zeroed element for the case when it is called from a templated class.

Parameters
  • params – the params to use.

  • format – - EVALUATION or COEFFICIENT

static inline std::function<DerivedType()> MakeDiscreteGaussianCoefficientAllocator(const std::shared_ptr<Params> &params, Format resultFormat, double stddev)

Allocator for discrete uniform distribution.

Parameters
  • paramsParams instance that is is passed.

  • resultFormat – resultFormat for the polynomials generated.

  • stddev – standard deviation for the discrete gaussian generator.

Returns

the resulting vector.

static inline std::function<DerivedType()> MakeDiscreteUniformAllocator(const std::shared_ptr<Params> &params, Format format)

Allocator for discrete uniform distribution.

Parameters
  • paramsParams instance that is is passed.

  • format – format for the polynomials generated.

Returns

the resulting vector.

Friends

inline friend std::ostream &operator<<(std::ostream &os, const DerivedType &vec)

ostream operator

Parameters
  • os – the input preceding output stream

  • vec – the element to add to the output stream.

Returns

a resulting concatenated output stream

inline friend DerivedType operator+(const DerivedType &a, const DerivedType &b)

Element-element addition operator.

Parameters
  • a – first element to add.

  • b – second element to add.

Returns

the result of the addition operation.

inline friend DerivedType operator+(const DerivedType &a, const BigIntType &b)

Element-integer addition operator.

Parameters
  • a – first element to add.

  • b – integer to add.

Returns

the result of the addition operation.

inline friend DerivedType operator+(const BigIntType &a, const DerivedType &b)

BigIntType-element addition operator.

Parameters
  • a – integer to add.

  • b – element to add.

Returns

the result of the addition operation.

inline friend DerivedType operator+(const DerivedType &a, const std::vector<BigIntType> &b)

Element-integer addition operator with CRT integer.

Parameters
  • a – first element to add.

  • b – integer to add.

Returns

the result of the addition operation.

inline friend DerivedType operator+(const std::vector<BigIntType> &a, const DerivedType &b)

BigIntType-element addition operator with CRT integer.

Parameters
  • a – integer to add.

  • b – element to add.

Returns

the result of the addition operation.

inline friend DerivedType operator-(const DerivedType &a, const DerivedType &b)

Element-element subtraction operator.

Parameters
  • a – element to subtract from.

  • b – element to subtract.

Returns

the result of the subtraction operation.

inline friend DerivedType operator-(const DerivedType &a, const std::vector<BigIntType> &b)

Element-integer subtraction operator with CRT integer.

Parameters
  • a – first element to subtract.

  • b – integer to subtract.

Returns

the result of the subtraction operation.

inline friend DerivedType operator-(const std::vector<BigIntType> &a, const DerivedType &b)

BigIntType-element subtraction operator with CRT integer.

Parameters
  • a – integer to subtract.

  • b – element to subtract.

Returns

the result of the subtraction operation.

inline friend DerivedType operator-(const DerivedType &a, const BigIntType &b)

Element-integer subtraction operator.

Parameters
  • a – element to subtract from.

  • b – integer to subtract.

Returns

the result of the subtraction operation.

inline friend DerivedType operator*(const DerivedType &a, const DerivedType &b)

Element-element multiplication operator.

Parameters
  • a – element to multiply.

  • b – element to multiply.

Returns

the result of the multiplication operation.

inline friend DerivedType operator*(const DerivedType &a, const BigIntType &b)

Element-integer multiplication operator.

Parameters
  • a – element to multiply.

  • b – integer to multiply.

Returns

the result of the multiplication operation.

inline friend DerivedType operator*(const DerivedType &a, const std::vector<BigIntType> &b)

Element-CRT number multiplication operator.

Parameters
  • a – element to multiply.

  • b – integer to multiply, in CRT format.

Returns

the result of the multiplication operation.

inline friend DerivedType operator*(const BigIntType &a, const DerivedType &b)

BigIntType-element multiplication operator.

Parameters
  • a – integer to multiply.

  • b – element to multiply.

Returns

the result of the multiplication operation.

inline friend DerivedType operator*(const DerivedType &a, int64_t b)

Element-signed-integer multiplication operator.

Parameters
  • a – element to multiply.

  • b – integer to multiply.

Returns

the result of the multiplication operation.

inline friend DerivedType operator*(int64_t a, const DerivedType &b)

signed-BigIntType-element multiplication operator.

Parameters
  • a – integer to multiply.

  • b – element to multiply.

Returns

the result of the multiplication operation.

struct CRTBasisExtensionPrecomputations

Public Functions

inline CRTBasisExtensionPrecomputations(const std::shared_ptr<Params> &paramsQlPl0, const std::shared_ptr<Params> &paramsPl0, const std::shared_ptr<Params> &paramsQl0, const std::vector<NativeInteger> &mPlQHatInvModq0, const std::vector<NativeInteger> &mPlQHatInvModqPrecon0, const std::vector<std::vector<NativeInteger>> &qInvModp0, const std::vector<DoubleNativeInt> &modpBarrettMu0, const std::vector<NativeInteger> &PlHatInvModp0, const std::vector<NativeInteger> &PlHatInvModpPrecon0, const std::vector<std::vector<NativeInteger>> &PlHatModq0, const std::vector<std::vector<NativeInteger>> &alphaPlModq0, const std::vector<DoubleNativeInt> &modqBarrettMu0, const std::vector<double> &pInv0)

Public Members

std::shared_ptr<Params> paramsQlPl
std::shared_ptr<Params> paramsPl
std::shared_ptr<Params> paramsQl
std::vector<NativeInteger> mPlQHatInvModq
std::vector<NativeInteger> mPlQHatInvModqPrecon
std::vector<std::vector<NativeInteger>> qInvModp
std::vector<DoubleNativeInt> modpBarrettMu
std::vector<NativeInteger> PlHatInvModp
std::vector<NativeInteger> PlHatInvModpPrecon
std::vector<std::vector<NativeInteger>> PlHatModq
std::vector<std::vector<NativeInteger>> alphaPlModq
std::vector<DoubleNativeInt> modqBarrettMu
std::vector<double> pInv