Plda Class Reference

#include <plda.h>

Collaboration diagram for Plda:

Public Member Functions

 Plda ()
 
 Plda (const Plda &other)
 
double TransformIvector (const PldaConfig &config, const VectorBase< double > &ivector, int32 num_enroll_examples, VectorBase< double > *transformed_ivector) const
 Transforms an iVector into a space where the within-class variance is unit and between-class variance is diagonalized. More...
 
float TransformIvector (const PldaConfig &config, const VectorBase< float > &ivector, int32 num_enroll_examples, VectorBase< float > *transformed_ivector) const
 float version of the above (not BaseFloat because we'd be implementing it twice for the same type if BaseFloat == double). More...
 
double LogLikelihoodRatio (const VectorBase< double > &transformed_enroll_ivector, int32 num_enroll_utts, const VectorBase< double > &transformed_test_ivector) const
 Returns the log-likelihood ratio log (p(test_ivector | same) / p(test_ivector | different)). More...
 
void SmoothWithinClassCovariance (double smoothing_factor)
 This function smooths the within-class covariance by adding to it, smoothing_factor (e.g. More...
 
void ApplyTransform (const Matrix< double > &in_transform)
 Apply a transform to the PLDA model. More...
 
int32 Dim () const
 
void Write (std::ostream &os, bool binary) const
 
void Read (std::istream &is, bool binary)
 

Protected Member Functions

void ComputeDerivedVars ()
 

Protected Attributes

Vector< double > mean_
 
Matrix< double > transform_
 
Vector< double > psi_
 
Vector< double > offset_
 

Private Member Functions

Pldaoperator= (const Plda &other)
 
double GetNormalizationFactor (const VectorBase< double > &transformed_ivector, int32 num_examples) const
 This returns a normalization factor, which is a quantity we must multiply "transformed_ivector" by so that it has the length that it "should" have. More...
 

Friends

class PldaEstimator
 
class PldaUnsupervisedAdaptor
 

Detailed Description

Definition at line 74 of file plda.h.

Constructor & Destructor Documentation

◆ Plda() [1/2]

Plda ( )
inline

Definition at line 76 of file plda.h.

76 { }

◆ Plda() [2/2]

Plda ( const Plda other)
inlineexplicit

Definition at line 78 of file plda.h.

78  :
79  mean_(other.mean_),
80  transform_(other.transform_),
81  psi_(other.psi_),
82  offset_(other.offset_) {
83  };
Vector< double > offset_
Definition: plda.h:155
Matrix< double > transform_
Definition: plda.h:149
Vector< double > psi_
Definition: plda.h:152
Vector< double > mean_
Definition: plda.h:148

Member Function Documentation

◆ ApplyTransform()

void ApplyTransform ( const Matrix< double > &  in_transform)

Apply a transform to the PLDA model.

This is mostly used for projecting the parameters of the model into a lower dimensional space, i.e. in_transform.NumRows() <= in_transform.NumCols(), typically for speaker diarization with a PCA transform.

Definition at line 220 of file plda.cc.

References SpMatrix< Real >::AddMat2Sp(), MatrixBase< Real >::AddMatMat(), VectorBase< Real >::AddMatVec(), VectorBase< Real >::ApplyFloor(), Plda::ComputeDerivedVars(), kaldi::ComputeNormalizingTransform(), VectorBase< Real >::CopyFromVec(), Plda::Dim(), SpMatrix< Real >::Eig(), SpMatrix< Real >::Invert(), KALDI_ASSERT, KALDI_WARN, kaldi::kNoTrans, kaldi::kTrans, Plda::mean_, VectorBase< Real >::Min(), rnnlm::n, MatrixBase< Real >::NumCols(), MatrixBase< Real >::NumRows(), Plda::psi_, Vector< Real >::Resize(), Matrix< Real >::Resize(), kaldi::SortSvd(), and Plda::transform_.

Referenced by main().

220  {
221  KALDI_ASSERT(in_transform.NumRows() <= Dim()
222  && in_transform.NumCols() == Dim());
223 
224  // Apply in_transform to mean_.
225  Vector<double> mean_new(in_transform.NumRows());
226  mean_new.AddMatVec(1.0, in_transform, kNoTrans, mean_, 0.0);
227  mean_.Resize(in_transform.NumRows());
228  mean_.CopyFromVec(mean_new);
229 
230  SpMatrix<double> between_var(in_transform.NumCols()),
231  within_var(in_transform.NumCols()),
232  psi_mat(in_transform.NumCols()),
233  between_var_new(Dim()),
234  within_var_new(Dim());
235  Matrix<double> transform_invert(transform_);
236 
237  // Next, compute the between_var and within_var that existed
238  // prior to diagonalization.
239  psi_mat.AddDiagVec(1.0, psi_);
240  transform_invert.Invert();
241  within_var.AddMat2(1.0, transform_invert, kNoTrans, 0.0);
242  between_var.AddMat2Sp(1.0, transform_invert, kNoTrans, psi_mat, 0.0);
243 
244  // Next, transform the variances using the input transformation.
245  between_var_new.AddMat2Sp(1.0, in_transform, kNoTrans, between_var, 0.0);
246  within_var_new.AddMat2Sp(1.0, in_transform, kNoTrans, within_var, 0.0);
247 
248  // Finally, we need to recompute psi_ and transform_. The remainder of
249  // the code in this function is a lightly modified copy of
250  // PldaEstimator::GetOutput().
251  Matrix<double> transform1(Dim(), Dim());
252  ComputeNormalizingTransform(within_var_new, &transform1);
253  // Now transform is a matrix that if we project with it,
254  // within_var becomes unit.
255  // between_var_proj is between_var after projecting with transform1.
256  SpMatrix<double> between_var_proj(Dim());
257  between_var_proj.AddMat2Sp(1.0, transform1, kNoTrans, between_var_new, 0.0);
258 
259  Matrix<double> U(Dim(), Dim());
260  Vector<double> s(Dim());
261  // Do symmetric eigenvalue decomposition between_var_proj = U diag(s) U^T,
262  // where U is orthogonal.
263  between_var_proj.Eig(&s, &U);
264 
265  KALDI_ASSERT(s.Min() >= 0.0);
266  int32 n;
267  s.ApplyFloor(0.0, &n);
268  if (n > 0) {
269  KALDI_WARN << "Floored " << n << " eigenvalues of between-class "
270  << "variance to zero.";
271  }
272  // Sort from greatest to smallest eigenvalue.
273  SortSvd(&s, &U);
274 
275  // The transform U^T will make between_var_proj diagonal with value s
276  // (i.e. U^T U diag(s) U U^T = diag(s)). The final transform that
277  // makes within_var unit and between_var diagonal is U^T transform1,
278  // i.e. first transform1 and then U^T.
279  transform_.Resize(Dim(), Dim());
280  transform_.AddMatMat(1.0, U, kTrans, transform1, kNoTrans, 0.0);
281  psi_.Resize(Dim());
282  psi_.CopyFromVec(s);
284 }
Matrix< double > transform_
Definition: plda.h:149
kaldi::int32 int32
void Resize(MatrixIndexT length, MatrixResizeType resize_type=kSetZero)
Set vector to a specified size (can be zero).
static void ComputeNormalizingTransform(const SpMatrix< Real > &covar, MatrixBase< Real > *proj)
This function computes a projection matrix that when applied makes the covariance unit (i...
Definition: plda.cc:46
void CopyFromVec(const VectorBase< Real > &v)
Copy data from another vector (must match own size).
Vector< double > psi_
Definition: plda.h:152
struct rnnlm::@11::@12 n
void AddMatMat(const Real alpha, const MatrixBase< Real > &A, MatrixTransposeType transA, const MatrixBase< Real > &B, MatrixTransposeType transB, const Real beta)
#define KALDI_WARN
Definition: kaldi-error.h:150
void ComputeDerivedVars()
Definition: plda.cc:57
Vector< double > mean_
Definition: plda.h:148
#define KALDI_ASSERT(cond)
Definition: kaldi-error.h:185
void Resize(const MatrixIndexT r, const MatrixIndexT c, MatrixResizeType resize_type=kSetZero, MatrixStrideType stride_type=kDefaultStride)
Sets matrix to a specified size (zero is OK as long as both r and c are zero).
int32 Dim() const
Definition: plda.h:140
void SortSvd(VectorBase< Real > *s, MatrixBase< Real > *U, MatrixBase< Real > *Vt, bool sort_on_absolute_value)
Function to ensure that SVD is sorted.

◆ ComputeDerivedVars()

void ComputeDerivedVars ( )
protected

Definition at line 57 of file plda.cc.

References VectorBase< Real >::AddMatVec(), Plda::Dim(), KALDI_ASSERT, kaldi::kNoTrans, Plda::mean_, Plda::offset_, Vector< Real >::Resize(), and Plda::transform_.

Referenced by Plda::ApplyTransform(), PldaEstimator::GetOutput(), Plda::Read(), and Plda::SmoothWithinClassCovariance().

57  {
58  KALDI_ASSERT(Dim() > 0);
59  offset_.Resize(Dim());
61 }
Vector< double > offset_
Definition: plda.h:155
Matrix< double > transform_
Definition: plda.h:149
void Resize(MatrixIndexT length, MatrixResizeType resize_type=kSetZero)
Set vector to a specified size (can be zero).
void AddMatVec(const Real alpha, const MatrixBase< Real > &M, const MatrixTransposeType trans, const VectorBase< Real > &v, const Real beta)
Add matrix times vector : this <– beta*this + alpha*M*v.
Definition: kaldi-vector.cc:92
Vector< double > mean_
Definition: plda.h:148
#define KALDI_ASSERT(cond)
Definition: kaldi-error.h:185
int32 Dim() const
Definition: plda.h:140

◆ Dim()

int32 Dim ( ) const
inline

◆ GetNormalizationFactor()

double GetNormalizationFactor ( const VectorBase< double > &  transformed_ivector,
int32  num_examples 
) const
private

This returns a normalization factor, which is a quantity we must multiply "transformed_ivector" by so that it has the length that it "should" have.

This comment explains the thinking behind the function LogLikelihoodRatio.

We assume "transformed_ivector" is an iVector in the transformed space (i.e., mean-subtracted, and multiplied by transform_). The covariance it "should" have in this space is + I/num_examples.

The reference is "Probabilistic Linear Discriminant Analysis" by Sergey Ioffe, ECCV 2006.

I'm looking at the un-numbered equation between eqs. (4) and (5), that says P(u^p | u^g_{1...n}) = N (u^p | {n }{n + I} {u}^g, I + {}{n + I})

Here, the superscript ^p refers to the "probe" example (e.g. the example to be classified), and u^g_1 is the first "gallery" example, i.e. the first training example of that class. is the between-class covariance matrix, assumed to be diagonalized, and I can be interpreted as the within-class covariance matrix which we have made unit.

We want the likelihood ratio P(u^p | u^g_{1..n}) / P(u^p), where the numerator is the probability of u^p given that it's in that class, and the denominator is the probability of u^p with no class assumption at all (e.g. in its own class).

The expression above even works for n = 0 (e.g. the denominator of the likelihood ratio), where it gives us P(u^p) = N(u^p | 0, I + ) i.e. it's distributed with zero mean and covarance (within + between). The likelihood ratio we want is: N(u^p | {n }{n + I} {u}^g, I + {}{n + I}) / N(u^p | 0, I + ) where {u}^g is the mean of the "gallery examples"; and we can expand the log likelihood ratio as

  • 0.5 [ (u^p - m) (I + /(n + I))^{-1} (u^p - m) + logdet(I + /(n + I)) ] + 0.5 [u^p (I + ) u^p + logdet(I + ) ] where m = (n )/(n + I) {u}^g.

Definition at line 99 of file plda.cc.

References VectorBase< Real >::Add(), VectorBase< Real >::ApplyPow(), Plda::Dim(), VectorBase< Real >::InvertElements(), KALDI_ASSERT, Plda::psi_, and kaldi::VecVec().

Referenced by Plda::TransformIvector().

101  {
102  KALDI_ASSERT(num_examples > 0);
103  // Work out the normalization factor. The covariance for an average over
104  // "num_examples" training iVectors equals \Psi + I/num_examples.
105  Vector<double> transformed_ivector_sq(transformed_ivector);
106  transformed_ivector_sq.ApplyPow(2.0);
107  // inv_covar will equal 1.0 / (\Psi + I/num_examples).
108  Vector<double> inv_covar(psi_);
109  inv_covar.Add(1.0 / num_examples);
110  inv_covar.InvertElements();
111  // "transformed_ivector" should have covariance (\Psi + I/num_examples), i.e.
112  // within-class/num_examples plus between-class covariance. So
113  // transformed_ivector_sq . (I/num_examples + \Psi)^{-1} should be equal to
114  // the dimension.
115  double dot_prod = VecVec(inv_covar, transformed_ivector_sq);
116  return sqrt(Dim() / dot_prod);
117 }
Vector< double > psi_
Definition: plda.h:152
#define KALDI_ASSERT(cond)
Definition: kaldi-error.h:185
int32 Dim() const
Definition: plda.h:140
Real VecVec(const VectorBase< Real > &a, const VectorBase< Real > &b)
Returns dot product between v1 and v2.
Definition: kaldi-vector.cc:37

◆ LogLikelihoodRatio()

double LogLikelihoodRatio ( const VectorBase< double > &  transformed_enroll_ivector,
int32  num_enroll_utts,
const VectorBase< double > &  transformed_test_ivector 
) const

Returns the log-likelihood ratio log (p(test_ivector | same) / p(test_ivector | different)).

transformed_enroll_ivector is an average over utterances for that speaker. Both transformed_enroll_vector and transformed_test_ivector are assumed to have been transformed by the function TransformIvector(). Note: any length normalization will have been done while computing the transformed iVectors.

Definition at line 153 of file plda.cc.

References VectorBase< Real >::Add(), VectorBase< Real >::AddVec(), VectorBase< Real >::ApplyPow(), Plda::Dim(), rnnlm::i, VectorBase< Real >::InvertElements(), kaldi::kUndefined, M_LOG_2PI, Plda::psi_, VectorBase< Real >::SumLog(), and kaldi::VecVec().

Referenced by main().

156  {
157  int32 dim = Dim();
158  double loglike_given_class, loglike_without_class;
159  { // work out loglike_given_class.
160  // "mean" will be the mean of the distribution if it comes from the
161  // training example. The mean is \frac{n \Psi}{n \Psi + I} \bar{u}^g
162  // "variance" will be the variance of that distribution, equal to
163  // I + \frac{\Psi}{n\Psi + I}.
164  Vector<double> mean(dim, kUndefined);
165  Vector<double> variance(dim, kUndefined);
166  for (int32 i = 0; i < dim; i++) {
167  mean(i) = n * psi_(i) / (n * psi_(i) + 1.0)
168  * transformed_train_ivector(i);
169  variance(i) = 1.0 + psi_(i) / (n * psi_(i) + 1.0);
170  }
171  double logdet = variance.SumLog();
172  Vector<double> sqdiff(transformed_test_ivector);
173  sqdiff.AddVec(-1.0, mean);
174  sqdiff.ApplyPow(2.0);
175  variance.InvertElements();
176  loglike_given_class = -0.5 * (logdet + M_LOG_2PI * dim +
177  VecVec(sqdiff, variance));
178  }
179  { // work out loglike_without_class. Here the mean is zero and the variance
180  // is I + \Psi.
181  Vector<double> sqdiff(transformed_test_ivector); // there is no offset.
182  sqdiff.ApplyPow(2.0);
183  Vector<double> variance(psi_);
184  variance.Add(1.0); // I + \Psi.
185  double logdet = variance.SumLog();
186  variance.InvertElements();
187  loglike_without_class = -0.5 * (logdet + M_LOG_2PI * dim +
188  VecVec(sqdiff, variance));
189  }
190  double loglike_ratio = loglike_given_class - loglike_without_class;
191  return loglike_ratio;
192 }
#define M_LOG_2PI
Definition: kaldi-math.h:60
Real SumLog() const
Returns sum of the logs of the elements.
kaldi::int32 int32
Vector< double > psi_
Definition: plda.h:152
struct rnnlm::@11::@12 n
int32 Dim() const
Definition: plda.h:140
Real VecVec(const VectorBase< Real > &a, const VectorBase< Real > &b)
Returns dot product between v1 and v2.
Definition: kaldi-vector.cc:37

◆ operator=()

Plda& operator= ( const Plda other)
private

◆ Read()

void Read ( std::istream &  is,
bool  binary 
)

Definition at line 34 of file plda.cc.

References Plda::ComputeDerivedVars(), kaldi::ExpectToken(), Plda::mean_, Plda::psi_, Vector< Real >::Read(), Matrix< Real >::Read(), and Plda::transform_.

34  {
35  ExpectToken(is, binary, "<Plda>");
36  mean_.Read(is, binary);
37  transform_.Read(is, binary);
38  psi_.Read(is, binary);
39  ExpectToken(is, binary, "</Plda>");
41 }
Matrix< double > transform_
Definition: plda.h:149
void Read(std::istream &in, bool binary, bool add=false)
read from stream.
Vector< double > psi_
Definition: plda.h:152
void ExpectToken(std::istream &is, bool binary, const char *token)
ExpectToken tries to read in the given token, and throws an exception on failure. ...
Definition: io-funcs.cc:191
void ComputeDerivedVars()
Definition: plda.cc:57
Vector< double > mean_
Definition: plda.h:148
void Read(std::istream &in, bool binary, bool add=false)
Read function using C++ streams.

◆ SmoothWithinClassCovariance()

void SmoothWithinClassCovariance ( double  smoothing_factor)

This function smooths the within-class covariance by adding to it, smoothing_factor (e.g.

0.1) times the between-class covariance (it's implemented by modifying transform_). This is to compensate for situations where there were too few utterances per speaker get a good estimate of the within-class covariance, and where the leading elements of psi_ were as a result very large.

We now revise our estimate of the within-class covariance to this larger value. This means that the transform has to change to as to make this new, larger covariance unit. And our between-class covariance in this space is now less.

Definition at line 195 of file plda.cc.

References VectorBase< Real >::AddVec(), VectorBase< Real >::ApplyPow(), Plda::ComputeDerivedVars(), Plda::Dim(), KALDI_ASSERT, KALDI_LOG, MatrixBase< Real >::MulRowsVec(), Plda::psi_, VectorBase< Real >::Set(), and Plda::transform_.

Referenced by main().

195  {
196  KALDI_ASSERT(smoothing_factor >= 0.0 && smoothing_factor <= 1.0);
197  // smoothing_factor > 1.0 is possible but wouldn't really make sense.
198 
199  KALDI_LOG << "Smoothing within-class covariance by " << smoothing_factor
200  << ", Psi is initially: " << psi_;
201  Vector<double> within_class_covar(Dim());
202  within_class_covar.Set(1.0); // It's now the current within-class covariance
203  // (a diagonal matrix) in the space transformed
204  // by transform_.
205  within_class_covar.AddVec(smoothing_factor, psi_);
210 
211  psi_.DivElements(within_class_covar);
212  KALDI_LOG << "New value of Psi is " << psi_;
213 
214  within_class_covar.ApplyPow(-0.5);
215  transform_.MulRowsVec(within_class_covar);
216 
218 }
Matrix< double > transform_
Definition: plda.h:149
Vector< double > psi_
Definition: plda.h:152
void ComputeDerivedVars()
Definition: plda.cc:57
void MulRowsVec(const VectorBase< Real > &scale)
Equivalent to (*this) = diag(scale) * (*this).
#define KALDI_ASSERT(cond)
Definition: kaldi-error.h:185
void ApplyPow(Real power)
Take all elements of vector to a power.
Definition: kaldi-vector.h:179
int32 Dim() const
Definition: plda.h:140
#define KALDI_LOG
Definition: kaldi-error.h:153

◆ TransformIvector() [1/2]

double TransformIvector ( const PldaConfig config,
const VectorBase< double > &  ivector,
int32  num_enroll_examples,
VectorBase< double > *  transformed_ivector 
) const

Transforms an iVector into a space where the within-class variance is unit and between-class variance is diagonalized.

The only anticipated use of this function is to pre-transform iVectors before giving them to the function LogLikelihoodRatio (it's done this way for efficiency because a given iVector may be used multiple times in LogLikelihoodRatio and we don't want to repeat the matrix multiplication

If config.normalize_length == true, it will also normalize the iVector's length by multiplying by a scalar that ensures that ivector^T inv_var ivector = dim. In this case, "num_enroll_examples" comes into play because it affects the expected covariance matrix of the iVector. The normalization factor is returned, even if config.normalize_length == false, in which case the normalization factor is computed but not applied. If config.simple_length_normalization == true, then an alternative normalization factor is computed that causes the iVector length to be equal to the square root of the iVector dimension.

Definition at line 120 of file plda.cc.

References VectorBase< Real >::AddMatVec(), VectorBase< Real >::CopyFromVec(), VectorBase< Real >::Dim(), Plda::Dim(), Plda::GetNormalizationFactor(), KALDI_ASSERT, kaldi::kNoTrans, VectorBase< Real >::Norm(), PldaConfig::normalize_length, Plda::offset_, VectorBase< Real >::Scale(), PldaConfig::simple_length_norm, and Plda::transform_.

Referenced by main(), Plda::TransformIvector(), and kaldi::TransformIvectors().

123  {
124  KALDI_ASSERT(ivector.Dim() == Dim() && transformed_ivector->Dim() == Dim());
125  double normalization_factor;
126  transformed_ivector->CopyFromVec(offset_);
127  transformed_ivector->AddMatVec(1.0, transform_, kNoTrans, ivector, 1.0);
128  if (config.simple_length_norm)
129  normalization_factor = sqrt(transformed_ivector->Dim())
130  / transformed_ivector->Norm(2.0);
131  else
132  normalization_factor = GetNormalizationFactor(*transformed_ivector,
133  num_examples);
134  if (config.normalize_length)
135  transformed_ivector->Scale(normalization_factor);
136  return normalization_factor;
137 }
Vector< double > offset_
Definition: plda.h:155
Matrix< double > transform_
Definition: plda.h:149
double GetNormalizationFactor(const VectorBase< double > &transformed_ivector, int32 num_examples) const
This returns a normalization factor, which is a quantity we must multiply "transformed_ivector" by so...
Definition: plda.cc:99
#define KALDI_ASSERT(cond)
Definition: kaldi-error.h:185
int32 Dim() const
Definition: plda.h:140

◆ TransformIvector() [2/2]

float TransformIvector ( const PldaConfig config,
const VectorBase< float > &  ivector,
int32  num_enroll_examples,
VectorBase< float > *  transformed_ivector 
) const

float version of the above (not BaseFloat because we'd be implementing it twice for the same type if BaseFloat == double).

Definition at line 140 of file plda.cc.

References VectorBase< Real >::CopyFromVec(), VectorBase< Real >::Dim(), and Plda::TransformIvector().

143  {
144  Vector<double> tmp(ivector), tmp_out(ivector.Dim());
145  float ans = TransformIvector(config, tmp, num_examples, &tmp_out);
146  transformed_ivector->CopyFromVec(tmp_out);
147  return ans;
148 }
double TransformIvector(const PldaConfig &config, const VectorBase< double > &ivector, int32 num_enroll_examples, VectorBase< double > *transformed_ivector) const
Transforms an iVector into a space where the within-class variance is unit and between-class variance...
Definition: plda.cc:120
void CopyFromVec(const VectorBase< Real > &v)
Copy data from another vector (must match own size).
MatrixIndexT Dim() const
Returns the dimension of the vector.
Definition: kaldi-vector.h:64

◆ Write()

void Write ( std::ostream &  os,
bool  binary 
) const

Definition at line 26 of file plda.cc.

References Plda::mean_, Plda::psi_, Plda::transform_, VectorBase< Real >::Write(), MatrixBase< Real >::Write(), and kaldi::WriteToken().

26  {
27  WriteToken(os, binary, "<Plda>");
28  mean_.Write(os, binary);
29  transform_.Write(os, binary);
30  psi_.Write(os, binary);
31  WriteToken(os, binary, "</Plda>");
32 }
void Write(std::ostream &out, bool binary) const
write to stream.
void Write(std::ostream &Out, bool binary) const
Writes to C++ stream (option to write in binary).
Matrix< double > transform_
Definition: plda.h:149
Vector< double > psi_
Definition: plda.h:152
void WriteToken(std::ostream &os, bool binary, const char *token)
The WriteToken functions are for writing nonempty sequences of non-space characters.
Definition: io-funcs.cc:134
Vector< double > mean_
Definition: plda.h:148

Friends And Related Function Documentation

◆ PldaEstimator

friend class PldaEstimator
friend

Definition at line 145 of file plda.h.

◆ PldaUnsupervisedAdaptor

friend class PldaUnsupervisedAdaptor
friend

Definition at line 146 of file plda.h.

Member Data Documentation

◆ mean_

◆ offset_

Vector<double> offset_
protected

Definition at line 155 of file plda.h.

Referenced by Plda::ComputeDerivedVars(), and Plda::TransformIvector().

◆ psi_

◆ transform_


The documentation for this class was generated from the following files: