ATLAS Offline Software
AlSymMat.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #include "GaudiKernel/StatusCode.h"
6 
10 #include "TrkAlgebraUtils/AlMat.h"
11 
12 #include <algorithm>
13 #include <cfloat> //for DBL_EPSILON
14 #include <cmath>
15 #include <cstdint>
16 #include <fstream>
17 #include <iomanip>
18 #include <iostream>
19 #include <stdexcept>
20 
21 #include <TMatrixDSparse.h>
22 
23 
24 extern "C" {
25  void dsptrf_(char*, const int*, double*, int*, int* );
26 }
27 
28 extern "C" {
29  void dsptri_(char*, const int*, double*, int*, double*, int* );
30 }
31 
32 extern "C" {
33  void dspev_(char*, char*, const int*, double*,
34  double*, double*, const int*, double*, int* );
35 }
36 
37 namespace Trk {
38 
39 //______________________________________________________________________________
41  : m_ptr_data(nullptr)
42  , m_pathbin("./")
43  , m_pathtxt("./")
44 {
45  m_matrix_type = 1;
46 }
47 
48 
49 //______________________________________________________________________________
51  : AlSymMatBase(N)
52  , m_ptr_data(new double[m_nele])
53  , m_pathbin("./")
54  , m_pathtxt("./")
55 {
56  m_matrix_type = 1;
57 
58  double* p = m_ptr_data + m_nele;
59  while (p > m_ptr_data) *(--p) = 0.;
60 
61 }
62 
63 //______________________________________________________________________________
65  : AlSymMatBase(m)
66  , m_ptr_data(new double[m_nele])
67  , m_pathbin(m.m_pathbin)
68  , m_pathtxt(m.m_pathtxt)
69 {
70  m_matrix_type = 1;
71  m_size = m.size();
72  m_nele = m.m_nele;
73 
74  copy(m);
75 }
76 
77 //______________________________________________________________________________
79  : AlSymMatBase (m.size()),
80  m_ptr_data (new double[m_nele]),
81  m_pathbin (m.pathBin()),
82  m_pathtxt (m.pathTxt())
83 {
85 }
86 
87 //______________________________________________________________________________
89 {
90  m_matrix_type = 1;
91  m_size = m.size();
92  m_nele = m_size*(m_size+1)/2;
93  m_pathbin = m.pathBin();
94  m_pathtxt = m.pathTxt();
95 
96  m_ptr_data = new double[m_nele];
97  copy(m);
98  return *this;
99 }
100 
101 //______________________________________________________________________________
103 {
104  if( m_ptr_data != nullptr ) delete [] m_ptr_data;
105  //ptr_map.clear();
106 }
107 
108 //______________________________________________________________________________
110 {
111  // this one implements the fast copy alg.
112  if( size() != m.size()) {
113  throw std::range_error( "AlSymMat::copy: size do not match!" );
114  }
115 
116  double * p = m_ptr_data + m_nele;
117  double* r = m.m_ptr_data + m_nele;
118  while (p > m_ptr_data) *(--p) = *(--r);
119 }
120 
121 //______________________________________________________________________________
123 {
124  //
125  if( size() != m.size()) {
126  throw std::range_error( "AlSymMat::copy: size do not match!" );
127  }
128 
129  (*this) = 0.;
130  //Convert Storage System
131  long int i, j;
133  for (pos = m.ptrMap()->begin(); pos!=m.ptrMap()->end(); pos++) {
134  m.elem(pos->first, i, j);
135  AlSymMat::elemr(i,j) = pos->second;
136  }
137 
138  }
139 
140 //______________________________________________________________________________
141 void AlSymMat::copy(const AlMat& m)
142 {
143  // copy the lower triangle only!
144  int si = size();
145  if( si != m.nrow() || si != m.ncol() ) {
146  throw std::range_error("AlSymMat::copy: sizes do not match!" );
147  }
148 
149  for( int i=0; i<si; i++ )
150  for( int j=0; j<=i; j++ )
151  elemr(i, j) = m.elemc(i, j);
152 }
153 
154 //______________________________________________________________________________
156 {
157  if (this==&m)
158  return *this;
159 
160  if ( m_ptr_data != m.m_ptr_data ) {
161  reSize(m.size());
162  copy(m);
163  }
164  return *this;
165 }
166 
167 //______________________________________________________________________________
169 {
170  if( m.nrow() != m.ncol() ) {
171  throw std::range_error( "AlSymMat::operator=: a squared matrix is required!" );
172  }
173 
174  reSize(m.nrow());
175  copy(m);
176 
177  return *this;
178 }
179 
180 //______________________________________________________________________________
182 {
183  double * p = m_ptr_data + m_nele;
184  while (p > m_ptr_data) *(--p) = d;
185 
186  return *this;
187 }
188 
189 //______________________________________________________________________________
191 {
192  if( size() != m.size()) {
193  throw std::range_error( "AlSymMat::operator+: size do not match!" );
194  }
195 
196  AlSymMat b(size());
197  double * p = m_ptr_data + m_nele;
198  double * q = m.m_ptr_data + m_nele;
199  double * r = b.m_ptr_data + m_nele;
200  while (p > m_ptr_data) *(--r) = (*(--p))+(*(--q));
201 
202  return b;
203 }
204 
205 //______________________________________________________________________________
207 {
208  if( size() != m.size()){
209  throw std::range_error( "AlSymMat::operator+=: size do not match!" );
210  }
211 
212  double * p = m_ptr_data + m_nele;
213  double * q = m.m_ptr_data + m_nele;
214  while (p > m_ptr_data) *(--p) += *(--q);
215 
216  return *this;
217 }
218 
219 //______________________________________________________________________________
221 {
222  if( size() != m.size()) {
223  throw std::range_error( "AlSymMat::operator-: size do not match!" );
224  }
225 
226  AlSymMat b(size());
227  double * p = m_ptr_data + m_nele;
228  double * q = m.m_ptr_data + m_nele;
229  double * r = b.m_ptr_data + m_nele;
230  while (p > m_ptr_data) *(--r) = (*(--p))-(*(--q));
231 
232  return b;
233 }
234 
235 //______________________________________________________________________________
237 {
238  if( size() != m.size()) {
239  throw std::range_error( "AlSymMat::operator-=: size do not match!" );
240  }
241 
242  double * p = m_ptr_data + m_nele;
243  double * q = m.m_ptr_data + m_nele;
244  while (p > m_ptr_data) *(--p) -= *(--q);
245 
246  return *this;
247 }
248 
249 //______________________________________________________________________________
251 {
252  if( size() != m.size() ) {
253  throw std::range_error( "AlSymMat::operator*: size do not match!" );
254  }
255 
256  AlMat b( size(), size());
257  double x;
258  long int k, l, n, ii, jj;
259  long int siz(size());
260 
261  // optimised version of symmetric matrix product:
262  for( long int i=0; i<siz; i++ )
263  for( long int j=0; j<siz; j++ ) {
264  x = 0.;
265  l = i < j ? i : j;
266  n = i < j ? j : i;
267  ii = (i+1)*i/2;
268  jj = (j+1)*j/2;
269  for( k=0; k<l; k++)
270  x += (*(m_ptr_data+ii+k))*(*(m.m_ptr_data+jj+k));
271  if( i<j ) {
272  for( k=l; k<n; k++)
273  x += (*(m_ptr_data+(k+1)*k/2+i))*(*(m.m_ptr_data+jj+k));
274  }
275  else {
276  for( k=l; k<n; k++)
277  x += (*(m_ptr_data+ii+k))*(*(m.m_ptr_data+(k+1)*k/2+j));
278  }
279  for( k=n; k<siz; k++)
280  x += (*(m_ptr_data+(k+1)*k/2+i))*(*(m.m_ptr_data+(k+1)*k/2+j));
281  b.elemr(i,j) = x;
282  }
283 
284  return b;
285 }
286 
287 //______________________________________________________________________________
289 {
290  if( size() != m.nrow() ) {
291  throw std::range_error("AlSymMat::operator*: size do not match!" );
292  }
293 
294  AlMat b( size(), m.ncol());
295 
296  double x;
297  long int k, ii;
298  long int siz(size()), nco(m.ncol());
299 
300  // optimised version of matrix product:
301  for( long int i=0; i<siz; i++ )
302  for( long int j=0; j<nco; j++ ) {
303  x = 0.;
304  ii = (i+1)*i/2;
305  for( k=0; k<i; k++)
306  x += (*(m_ptr_data+ii+k))*m.elemc(k,j);
307  for( k=i; k<siz; k++)
308  x += (*(m_ptr_data+(k+1)*k/2+i))*m.elemc(k,j);
309  b.elemr(i,j) = x;
310  }
311 
312  return b;
313 }
314 
315 //______________________________________________________________________________
317 {
318  if( size() != v.size() ) {
319  throw std::range_error( "AlSymMat::operator*: size do not match! " );
320  }
321 
322  AlVec b(size());
323 
324  long int ii;
325 
326  for( long int i=0; i<size(); i++ ) {
327  ii = (i+1)*i/2;
328  for( long int j=0; j<i; j++ )
329  b[i] += (*(m_ptr_data+ii+j))*v[j];
330  for( long int j=i; j<size(); j++ )
331  b[i] += (*(m_ptr_data+(j+1)*j/2+i))*v[j];
332  }
333 
334  return b;
335 }
336 
337 //______________________________________________________________________________
339 {
340  double * p = m_ptr_data + m_nele;
341  while (p > m_ptr_data) *(--p) *= d;
342 
343  return *this;
344 }
345 
346 //______________________________________________________________________________
347 AlSymMat AlSymMat::operator*(const double& d) const
348 {
349  AlSymMat a(size());
350  double * p = m_ptr_data + m_nele;
351  double * q = a.m_ptr_data + m_nele;
352  while (p > m_ptr_data) *(--q) = (*(--p))*d;
353 
354  return a;
355 }
356 
357 //______________________________________________________________________________
359 {
360  // add a protection to detect singular matrices:
361  // double det = (*this).determinant();
362  // double detcut = 1.E-20;
363 
364  // if (fabs(det)>detcut) {
365 
366  int ierr = 0;
367  char uplo = 'U';
368  int N = size();
369  int * ipiv = new int[N];
370  double * work = new double[N];
371  double * mx = m_ptr_data;
372 
373  dsptrf_(&uplo, &N, mx, ipiv, &ierr);
374 
375  dsptri_(&uplo, &N, mx, ipiv, work, &ierr);
376 
377  delete [] ipiv;
378  delete [] work;
379 
380  // } else {
381  // ierr=999;
382  // }
383  return ierr;
384 }
385 
386 //______________________________________________________________________________
388 {
389  double deter = 1.;
390 
391  // get determinant using LU decomposition + Crout algorithm
392  AlMat A(size(),size());
393  A = (*this);
394 
395  double AMAX,DUM,SUM;
396  //int CODE;
397  double TINY = 1.E-20;
398  double D;
399  int N=size();
400 
401  D=1.; //CODE=0;
402 
403  for(int I=0;I<N;I++) {
404  AMAX=0.;
405  for(int J=0;J<N;J++)
406  if (fabs(A[I][J])>AMAX)
407  AMAX=fabs(A[I][J]);
408  if(AMAX<TINY) {
409 // CODE = 1; // matrix is singular
410  return 0;
411  }
412  }
413 
414  for(int J=0;J<N;J++) {
415  for(int I=0;I<N;I++) {
416  SUM = A[I][J];
417  int KMAX = (I<J)?I:J;
418  for(int K=0;K<KMAX;K++)
419  SUM = SUM - A[I][K]*A[K][J];
420  A[I][J]=SUM;
421  }
422 
423  // find pivot and exchange if necessary
424  int IMAX=J;
425  for (int I=J+1;I<N;I++) {
426  if(fabs(A[I][J])>fabs(A[IMAX][J])) IMAX=I;
427  }
428 
429  if(IMAX!=J) {
430  for(int K=0;K<N;K++) {
431  DUM=A[IMAX][K];
432  A[IMAX][K]=A[J][K];
433  A[J][K]=DUM;
434  }
435  D=-D;
436  }
437 
438  if(J<N && (fabs(A[J][J])>TINY)) {
439  DUM = 1./ A[J][J];
440  for(int I=J+1;I<N;I++) A[I][J] = A[I][J]*DUM;
441  }
442  }
443 
444  deter=deter*D;
445  for(int I=0;I<N;I++)
446  deter = deter*A[I][I];
447 
448  return deter;
449 }
450 
451 //______________________________________________________________________________
452 int AlSymMat::diagonalize(char jobz, AlVec &w, AlMat &z)
453 {
454  int info = 0;
455  char uplo = 'U';
456  int N = size();
457  double * work = new double[3*N];
458 
459  double * ap = m_ptr_data;
460 
461  dspev_(&jobz, &uplo, &N, ap, w.ptrData(), z.ptrData(), &N, work, &info);
462 
463  delete [] work;
464 
465  return info;
466 }
467 
468 //______________________________________________________________________________
469 // index here is module number within matrix (size=6*max_index)
471 {
472  const int shift=6;
473 
474  if((shift*index)<size() && index>=0) {
475  for(int row=(shift*index); row<(size()-shift); row+=shift) {
476  for(int k=0; k<shift; k++) {
477  for(int col=0; col<row+k+1; col++) {
478  // cout << "(" << row+k << "," << col << ") -> " << *(m_ptr_data+elem(row+k,col)) << endl;
479 
480  if(col<shift*index)
481  *(m_ptr_data+elem(row+k,col))=*(m_ptr_data+elem(row+k+shift,col));
482  else
483  *(m_ptr_data+elem(row+k,col))=*(m_ptr_data+elem(row+k+shift,col+shift));
484  }
485  }
486  }
487  }
488 }
489 
490 
491 //______________________________________________________________________________
492 int AlSymMat::RemoveCollsRows(std::vector<int> indices)
493 {
494  int n = indices.size();
495  if (n==0) {
496  return m_size;
497  }
498  if (n>m_size) {
499  throw std::range_error("Vector of indices larger than matrix size.");
500  }
501 
502  // first sort the list of indices descending
503  // maybe not the fastest way to do that but it works
504  for (int i=0;i<n-1;i++)
505  for (int j=i+1; j<n; j++)
506  if (indices[j] > indices[i]) {
507  int itmp = indices[i];
508  indices[i] = indices[j];
509  indices[j] = itmp;
510  }
511 
512  // remove rows and columns starting from largest indices
513  for (int i=0;i<n;i++) {
514  int index = indices[i];
515  if (index > m_size-1) {
516  throw std::out_of_range("AlSymMat::RemoveCollsRows: Index goes beyond matrix.");
517  }
518 
519  for (int j=index; j<m_size-1; j++)
520  for (int col=0; col<=j; col++)
521  if(col<index)
522  *(m_ptr_data+elem(j,col)) = *(m_ptr_data+elem(j+1,col));
523  else
524  *(m_ptr_data+elem(j,col)) = *(m_ptr_data+elem(j+1,col+1));
525 
526  m_size--;
527  }
528 
529  return m_size;
530 }
531 
532 //______________________________________________________________________________
533 // index here is alignment parameter to remove
534 void AlSymMat::RemoveAlignPar(int index, int control)
535 {
536  // Dynamic shiftRow
537  int shiftRow = 1;
538  int counterRow = 0;
539 
540  // Dynamic shiftCol
541  int shiftCol;
542  int counterCol = 0;
543 
544  index = index-control;
545 
546  for(int row=index; row<(size()-shiftRow); row++) {
547 
548  shiftCol = 0;
549  counterCol = 0;
550 
551  for(int col=0; col<row+1; col++) {
552 
553  if (col==index) {
554  shiftCol = 1;
555  counterCol=0;
556  }
557 
558  *(m_ptr_data+elem(row,col))=*(m_ptr_data+elem(row+shiftRow,col+shiftCol));
559 
560  counterCol++;
561  if (counterCol==5-control) {
562  counterCol=0;
563  shiftCol++;
564  }
565  }
566 
567  counterRow++;
568  if (counterRow==5-control) {
569  counterRow=0;
570  shiftRow++;
571  }
572  }
573 }
574 
575 //______________________________________________________________________________
576 void AlSymMat::SetPathBin(const std::string &path)
577 {
578  m_pathbin.assign(path);
579 }
580 
581 //______________________________________________________________________________
582 void AlSymMat::SetPathTxt(const std::string &path)
583 {
584  m_pathtxt.assign(path);
585 }
586 
587 //______________________________________________________________________________
588 StatusCode AlSymMat::Write(const std::string &filename, bool binary,
589  bool square, double, float version)
590 {
591  std::ofstream outmat;
592 
593  int32_t msizz = square ? size() : (-1)*size();
594 
595  if(binary) {
596  outmat.open((m_pathbin+filename).c_str(), std::ios::binary);
597  if(outmat.fail())
598  return StatusCode::FAILURE;
599  // if triangular, change size sign to distinguish
600  // from square format
601  outmat.write((char*)&msizz, sizeof (msizz));
602  outmat.write((char*)&version, sizeof (version));
603  }
604  else {
605  outmat.open((m_pathtxt+filename).c_str());
606  if(outmat.fail())
607  return StatusCode::FAILURE;
608  outmat.setf(std::ios::fixed);
609  outmat.setf(std::ios::showpoint);
610  outmat.precision(3);
611  outmat << "msizz: " << std::setw(6) << msizz << std::endl;
612  outmat << "AlSymMat version: " << std::setw(6) << version << std::endl;
613  }
614 
615  double melem=0;
616 
617  if (!square) {
618  for( int i=0; i<size(); i++) {
619  for( int j=0; j<=i; j++) {
620  melem = *(m_ptr_data+(i+1)*i/2+j);
621  if(binary)
622  outmat.write((char*)&(melem), sizeof (melem));
623  else
624  outmat << std::setw(14) << melem;
625  }
626  if(!binary)
627  outmat << std::endl;
628  }
629  }
630  else {
631  for( int i=0; i<size(); i++) {
632  for( int j=0; j<size(); j++) {
633  if(i>=j) melem = *(m_ptr_data+(i+1)*i/2+j);
634  else melem = *(m_ptr_data+(j+1)*j/2+i);
635  if(binary)
636  outmat.write((char*)&(melem), sizeof (melem));
637  else
638  outmat << std::setw(14) << melem;
639  }
640  if(!binary)
641  outmat << std::endl;
642  }
643  }
644  outmat.close();
645  return StatusCode::SUCCESS;
646 }
647 
648 //______________________________________________________________________________
649 StatusCode AlSymMat::CheckMatVersion(const std::string& filename, bool &StdUnits)
650 {
651  std::ifstream inmat((filename).c_str(), std::ios::binary);
652  if(inmat.fail())
653  return StatusCode::FAILURE;
654 
655  int32_t msiz=0;
656  inmat.read((char*)&msiz, sizeof (msiz));
657 
658  float version=0.0;
659  inmat.read((char*)&version, sizeof (version));
660 
661  StdUnits = version>=2.0;
662 
663  inmat.close();
664 
665  // std::cout << "AlSymMat::StdUnits: " << StdUnits << std::endl;
666 
667  return StatusCode::SUCCESS;
668 }
669 
670 //______________________________________________________________________________
671 StatusCode AlSymMat::Read(const std::string &filename, int &dofs,
672  bool &triang, float &version)
673 {
674  bool stdUnits = true;
675  if (StatusCode::SUCCESS != CheckMatVersion(m_pathbin+filename, stdUnits))
676  return StatusCode::FAILURE;
677 
678  // std::cout << "AlSymMat::StdUnits: " << stdUnits << std::endl;
679 
680  std::ifstream inmat((m_pathbin+filename).c_str(), std::ios::binary);
681  if(inmat.fail())
682  return StatusCode::FAILURE;
683 
684  int32_t msiz=0;
685  inmat.read((char*)&msiz, sizeof (msiz));
686  dofs = abs(msiz);
687  m_size = abs(msiz);
688 
689  if (stdUnits)
690  inmat.read((char*)&version, sizeof (version));
691 
692  double melem=0;
693 
694  if(msiz>0) { // square format
695  triang=false;
696  for(int i=0; i<msiz; i++) {
697  for(int j=0; j<msiz; j++) {
698  inmat.read((char*)&melem, sizeof (melem));
699  if( i>=j )
700  *(m_ptr_data+(i+1)*i/2+j) = melem;
701  }
702  }
703  }
704  else { // triangular format
705  triang=true;
706  msiz = (-1)*msiz;
707  // std::cout << "msiz=" << msiz << std::endl;
708  for( int i=0; i<msiz; i++) {
709  for( int j=0; j<=i; j++) {
710  inmat.read((char*)&melem, sizeof (melem));
711  *(m_ptr_data+(i+1)*i/2+j) = melem;
712  }
713  }
714  }
715 
716  inmat.close();
717  return StatusCode::SUCCESS;
718 }
719 
720 //______________________________________________________________________________
721 StatusCode AlSymMat::ReadProjected(const std::string &filename, int &dofs,
722  bool &triang, float &version)
723 {
724  std::ifstream inmat((m_pathbin+filename).c_str(), std::ios::binary);
725  if(inmat.fail())
726  return StatusCode::FAILURE;
727 
728  int32_t msiz=0;
729  inmat.read((char*)&msiz, sizeof (msiz));
730  dofs = abs(msiz);
731  m_size = abs(msiz);
732 
733  inmat.read((char*)&version, sizeof (version));
734 
735  double melem=0;
736 
737  if(msiz>0) { // square format
738  triang=false;
739  for(int i=0; i<msiz; i++) {
740  for(int j=0; j<msiz; j++) {
741  inmat.read((char*)&melem, sizeof (melem));
742  if( i>=j )
743  *(m_ptr_data+(i+1)*i/2+j) = melem;
744  }
745  }
746  }
747  else { // triangular format
748  triang=true;
749  msiz = (-1)*msiz;
750  // std::cout << "msiz=" << msiz << std::endl;
751  for( int i=0; i<msiz; i++) {
752  for( int j=0; j<=i; j++) {
753  inmat.read((char*)&melem, sizeof (melem));
754  *(m_ptr_data+(i+1)*i/2+j) = melem;
755  }
756  }
757  }
758 
759  inmat.close();
760  return StatusCode::SUCCESS;
761 }
762 
763 //______________________________________________________________________________
764 void AlSymMat::reSize(long int Nnew)
765 {
766  if ( Nnew != size() ) {
767  double * p = m_ptr_data;
768  long int ii;
769 
770  m_nele = Nnew*(Nnew+1)/2;
771  m_ptr_data = new double[m_nele];
772  long int Size_old = size();
773  m_size = Nnew;
774  long int k = size() <= Size_old ? size() : Size_old;
775  for( long int i=0; i<k; i++ ) {
776  ii = (i+1)*i/2;
777  for( long int j=0; j<=i; j++ )
778  *(m_ptr_data+ii+j) = *(p+ii+j);
779  }
780  delete [] p;
781  }
782 }
783 
784 //______________________________________________________________________________
785 double& AlSymMat::elemr(long int i,long int j)
786 {
787 #ifdef _DEBUG
788  if( i<0 ) {
789  throw std::underflow_error( "AlSymMat::elemr: Index 1 < zero!" );
790  }
791  if( i>=size() ) {
792  throw std::overflow_error( "AlSymMat::elemr: Index 1 too large!" );
793  }
794  if( j<0 ) {
795  throw std::underflow_error( "AlSymMat::elemr: Index 2 < zero!" );
796  }
797  if( j>=size() ) {
798  throw std::overflow_error( "AlSymMat::elemr: Index 2 too large!" );
799  }
800 #endif
801  if( j<=i )
802  return *(m_ptr_data+(i+1)*i/2+j);
803 
804  return *(m_ptr_data+(j+1)*j/2+i);
805 }
806 
807 //______________________________________________________________________________
808 double AlSymMat::elemc(long int i,long int j) const
809 {
810 #ifdef _DEBUG
811  if( i<0 ) {
812  throw std::underflow_error( "AlSymMat::elemc: Index 1 < zero!" );
813  }
814  if( i>=size() ) {
815  throw std::overflow_error( "AlSymMat::elemc: Index 1 too large!" );
816  }
817  if( j<0 ) {
818  throw std::underflow_error( "AlSymMat::elemc: Index 2 < zero!" );
819  }
820  if( j>=size() ) {
821  throw std::overflow_error( "AlSymMat::elemc: Index 2 too large!" );
822  }
823 #endif
824  if( j<=i )
825  return *(m_ptr_data+(i+1)*i/2+j);
826 
827  return *(m_ptr_data+(j+1)*j/2+i);
828 }
829 
830 
832 
833  int si = size();
834 
835  int nnz(0);
836 
837  for( int i=0; i<si; i++ )
838  for( int j=0; j<=i; j++ )
839  if( *(m_ptr_data+(i+1)*i/2+j) != 0 )
840  ++nnz;
841 
842 
843  int *irow = new int[nnz];
844  int *icol = new int[nnz];
845  double *val = new double[nnz];
846  int n(0);
847 
848  for (int i=0;i<si;i++) {
849  for (int j=0;j<=i;j++) {
850  if ( *(m_ptr_data+(i+1)*i/2+j) != 0.) {
851  irow[n] = i;
852  icol[n] = j;
853  val[n] = *(m_ptr_data+(i+1)*i/2+j);
854  ++n;
855  }
856  }
857  }
858 
859 
860  TMatrixDSparse* myTMatrix = new TMatrixDSparse(0,si-1,0,si-1);
861  myTMatrix->SetMatrixArray(nnz,irow,icol,val);
862 
863  delete [] irow;
864  delete [] icol;
865  delete [] val;
866 
867 
868  return myTMatrix;
869 
870 }
871 
872 
873 } // end namespace Trk
grepfile.info
info
Definition: grepfile.py:38
Trk::AlSymMat::RemoveAlignPar
virtual void RemoveAlignPar(int, int) override final
Definition: AlSymMat.cxx:534
query_example.row
row
Definition: query_example.py:24
beamspotman.r
def r
Definition: beamspotman.py:676
Trk::const_mapiterator
datamap::const_iterator const_mapiterator
Definition: AlSymMatBase.h:30
python.SystemOfUnits.m
int m
Definition: SystemOfUnits.py:91
athena.path
path
python interpreter configuration --------------------------------------—
Definition: athena.py:128
Trk::z
@ z
global position (cartesian)
Definition: ParamDefs.h:57
Trk::AlSymMat::operator+
virtual AlSymMat operator+(const AlSymMat &) const
Definition: AlSymMat.cxx:190
Trk::AlVec
Definition: AlVec.h:23
AlSpaMat.h
index
Definition: index.py:1
Trk::AlSymMatBase
Definition: AlSymMatBase.h:38
hist_file_dump.d
d
Definition: hist_file_dump.py:137
Trk::indices
std::pair< long int, long int > indices
Definition: AlSymMatBase.h:24
Trk::AlSymMat::Write
virtual StatusCode Write(const std::string &, bool, bool, double, float) override final
Definition: AlSymMat.cxx:588
CSV_InDetExporter.new
new
Definition: CSV_InDetExporter.py:145
Trk::AlSymMat::operator+=
virtual AlSymMat & operator+=(const AlSymMat &)
Definition: AlSymMat.cxx:206
AlVec.h
Trk::AlSymMat::m_pathbin
std::string m_pathbin
Definition: AlSymMat.h:89
Trk::AlSymMat::elem
long int elem(long int, long int) const
Definition: AlSymMat.h:96
fitman.mx
mx
Definition: fitman.py:520
UploadAMITag.l
list l
Definition: UploadAMITag.larcaf.py:158
JetTiledMap::N
@ N
Definition: TiledEtaPhiMap.h:44
AlSymMat.h
Trk::AlSymMat::invert
virtual int invert() override final
Definition: AlSymMat.cxx:358
Trk::AlSymMat::SetPathBin
virtual void SetPathBin(const std::string &) override final
Definition: AlSymMat.cxx:576
python.AtlRunQueryParser.ap
ap
Definition: AtlRunQueryParser.py:826
Trk::AlSymMat::operator*=
virtual AlSymMat & operator*=(const double &)
Definition: AlSymMat.cxx:338
atlasStyleMacro.icol
int icol
Definition: atlasStyleMacro.py:13
Trk::AlSymMat
Definition: AlSymMat.h:26
TMatrixDSparse
class TMatrixTSparse< double > TMatrixDSparse
Definition: AlSpaMat.h:14
Trk::AlSymMatBase::size
long int size() const
Definition: AlSymMatBase.h:145
Trk::AlSymMat::SetPathTxt
virtual void SetPathTxt(const std::string &) override final
Definition: AlSymMat.cxx:582
Trk::AlSymMat::operator=
AlSymMat & operator=(const AlSpaMat &)
Definition: AlSymMat.cxx:88
python.setupRTTAlg.size
int size
Definition: setupRTTAlg.py:39
A
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:210
lumiFormat.i
int i
Definition: lumiFormat.py:85
dsptrf_
void dsptrf_(char *, const int *, double *, int *, int *)
beamspotman.n
n
Definition: beamspotman.py:731
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
Trk::AlSymMat::RemoveModule
virtual void RemoveModule(int) override final
Definition: AlSymMat.cxx:470
Trk::AlSymMat::Read
virtual StatusCode Read(const std::string &, int &, bool &, float &) override final
Definition: AlSymMat.cxx:671
AlMat.h
Trk::AlMat
Definition: AlMat.h:27
Trk::AlSymMatBase::m_nele
long int m_nele
Definition: AlSymMatBase.h:104
Trk::AlSymMat::elemc
virtual double elemc(long int, long int) const override final
Definition: AlSymMat.cxx:808
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
Trk::AlSymMat::m_pathtxt
std::string m_pathtxt
Definition: AlSymMat.h:90
Trk::AlSymMat::makeTMatrix
virtual TMatrixDSparse * makeTMatrix() override final
Definition: AlSymMat.cxx:831
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Trk::AlSymMat::copy
void copy(const AlSymMat &m)
Definition: AlSymMat.cxx:109
plotBeamSpotMon.b
b
Definition: plotBeamSpotMon.py:77
query_example.col
col
Definition: query_example.py:7
dspev_
void dspev_(char *, char *, const int *, double *, double *, double *, const int *, double *, int *)
Trk::AlSymMat::ReadProjected
virtual StatusCode ReadProjected(const std::string &, int &, bool &, float &) override final
Definition: AlSymMat.cxx:721
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
get_generator_info.version
version
Definition: get_generator_info.py:33
Trk::AlSymMat::diagonalize
virtual int diagonalize(char jobz, AlVec &w, AlMat &z) override final
Definition: AlSymMat.cxx:452
a
TList * a
Definition: liststreamerinfos.cxx:10
Trk::AlSymMat::AlSymMat
AlSymMat()
Definition: AlSymMat.cxx:40
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
dsptri_
void dsptri_(char *, const int *, double *, int *, double *, int *)
Trk::AlSymMat::RemoveCollsRows
virtual int RemoveCollsRows(std::vector< int >) override final
Definition: AlSymMat.cxx:492
Trk::AlSymMat::reSize
virtual void reSize(long int Nnew) override final
Definition: AlSymMat.cxx:764
CaloCellTimeCorrFiller.filename
filename
Definition: CaloCellTimeCorrFiller.py:24
extractSporadic.q
list q
Definition: extractSporadic.py:98
Trk::AlSymMat::operator*
virtual AlMat operator*(const AlSymMat &) const
Definition: AlSymMat.cxx:250
Trk::AlSymMatBase::m_matrix_type
int m_matrix_type
Definition: AlSymMatBase.h:100
I
#define I(x, y, z)
Definition: MD5.cxx:116
Trk::AlSymMat::operator-
virtual AlSymMat operator-(const AlSymMat &) const
Definition: AlSymMat.cxx:220
python.IoTestsLib.w
def w
Definition: IoTestsLib.py:200
Trk::x
@ x
Definition: ParamDefs.h:55
Trk::AlSymMatBase::m_size
long int m_size
Definition: AlSymMatBase.h:103
Trk::AlSymMat::m_ptr_data
double * m_ptr_data
Definition: AlSymMat.h:87
Trk::AlSymMat::~AlSymMat
~AlSymMat()
Definition: AlSymMat.cxx:102
Trk::AlSymMat::determinant
virtual double determinant() override final
Definition: AlSymMat.cxx:387
Trk::AlSymMat::CheckMatVersion
virtual StatusCode CheckMatVersion(const std::string &, bool &) override final
Definition: AlSymMat.cxx:649
Trk::AlSpaMat
Definition: AlSpaMat.h:27
Trk::v
@ v
Definition: ParamDefs.h:78
Trk::AlSymMat::elemr
virtual double & elemr(long int, long int) override final
Definition: AlSymMat.cxx:785
fitman.k
k
Definition: fitman.py:528
Trk::AlSymMat::operator-=
virtual AlSymMat & operator-=(const AlSymMat &)
Definition: AlSymMat.cxx:236