-
Notifications
You must be signed in to change notification settings - Fork 0
/
SparseSolver.cpp
1103 lines (928 loc) · 25.8 KB
/
SparseSolver.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
* This file is part of qpOASES.
*
* qpOASES -- An Implementation of the Online Active Set Strategy.
*
* Copyright (C) 2012 by Andreas Waechter. All rights reserved.
*
* qpOASES is free software; you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public
* License as published by the Free Software Foundation; either
* version 2.1 of the License, or (at your option) any later version.
*
* qpOASES is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
* See the GNU Lesser General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with qpOASES; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
*
*/
/**
* \file src/SparseSolver.cpp
* \author Andreas Waechter, Dennis Janka
* \version 3.2
* \date 2012-2015
*
* Interfaces to sparse linear solvers that are used in a Schur-complement
* implementation in qpOASES.
*/
#include "SparseSolver.hpp"
#ifndef __MATLAB__
# include <cstdarg>
void MyPrintf(const char* pformat, ... );
#else
# include <mex.h>
# define MyPrintf mexPrintf
#endif
//#define __DEBUG_ITER__
BEGIN_NAMESPACE_QPOASES
/*****************************************************************************
* P U B L I C *
****************************************************************************/
/*
* S p a r s e S o l v e r B a s e
*/
SparseSolver::SparseSolver( )
{
}
/*
* S p a r s e S o l v e r B a s e
*/
SparseSolver::SparseSolver( const SparseSolver& rhs )
{
copy( rhs );
}
/*
* ~ S p a r s e S o l v e r B a s e
*/
SparseSolver::~SparseSolver( )
{
clear( );
}
/*
* o p e r a t o r =
*/
SparseSolver& SparseSolver::operator=( const SparseSolver& rhs )
{
if ( this != &rhs )
{
clear( );
copy( rhs );
}
return *this;
}
/*
* r e s e t
*/
returnValue SparseSolver::reset( )
{
return SUCCESSFUL_RETURN;
}
/*
* g e t N e g a t i v e E i g e n v a l u e s
*/
int_t SparseSolver::getNegativeEigenvalues( )
{
return -1;
}
/*
* g e t R a n k
*/
int_t SparseSolver::getRank( )
{
return -1;
}
/*
* g e t Z e r o P i v o t s
*/
returnValue SparseSolver::getZeroPivots( int_t *&zeroPivots )
{
if ( zeroPivots ) delete[] zeroPivots;
zeroPivots = 0;
return SUCCESSFUL_RETURN;
}
/*****************************************************************************
* P R O T E C T E D *
*****************************************************************************/
/*
* c l e a r
*/
returnValue SparseSolver::clear( )
{
return SUCCESSFUL_RETURN;
}
/*
* c o p y
*/
returnValue SparseSolver::copy( const SparseSolver& rhs
)
{
return SUCCESSFUL_RETURN;
}
#ifdef SOLVER_MA27
/*****************************************************************************
*****************************************************************************
*****************************************************************************
* M A 2 7 S P A R E S E S O L V E R *
*****************************************************************************
*****************************************************************************
*****************************************************************************/
#define MA27ID ma27id_
#define MA27AD ma27ad_
#define MA27BD ma27bd_
#define MA27CD ma27cd_
extern "C" {
void MA27ID(fint* ICNTL, double* CNTL);
void MA27AD(fint *N, fint *NZ, const fint *IRN, const fint* ICN,
fint *IW, fint* LIW, fint* IKEEP, fint *IW1,
fint* NSTEPS, fint* IFLAG, fint* ICNTL,
double* CNTL, fint *INFO, double* OPS);
void MA27BD(fint *N, fint *NZ, const fint *IRN, const fint* ICN,
double* A, fint* LA, fint* IW, fint* LIW,
fint* IKEEP, fint* NSTEPS, fint* MAXFRT,
fint* IW1, fint* ICNTL, double* CNTL,
fint* INFO);
void MA27CD(fint *N, double* A, fint* LA, fint* IW,
fint* LIW, double* W, fint* MAXFRT,
double* RHS, fint* IW1, fint* NSTEPS,
fint* ICNTL, double* CNTL);
}
/*****************************************************************************
* P U B L I C *
****************************************************************************/
/*
* S p a r s e S o l v e r B a s e
*/
Ma27SparseSolver::Ma27SparseSolver( ) : SparseSolver()
{
a_ma27 = 0;
irn_ma27 = 0;
jcn_ma27 = 0;
iw_ma27 = 0;
ikeep_ma27 = 0;
clear( );
/* Set default options for MA27 */
MA27ID(icntl_ma27, cntl_ma27);
icntl_ma27[0] = 0; // Suppress error messages
icntl_ma27[1] = 0; // Suppress diagnostic messages
cntl_ma27[0] = 1e-8; // Set pivot tolerance
}
/*
* S p a r s e S o l v e r B a s e
*/
Ma27SparseSolver::Ma27SparseSolver( const Ma27SparseSolver& rhs )
{
copy( rhs );
}
/*
* ~ S p a r s e S o l v e r B a s e
*/
Ma27SparseSolver::~Ma27SparseSolver( )
{
clear( );
}
/*
* o p e r a t o r =
*/
Ma27SparseSolver& Ma27SparseSolver::operator=( const SparseSolver& rhs )
{
const Ma27SparseSolver* ma27_rhs = dynamic_cast<const Ma27SparseSolver*>(&rhs);
if (!ma27_rhs)
{
fprintf(getGlobalMessageHandler()->getOutputFile(),"Error in Ma27SparseSolver& Ma27SparseSolver::operator=( const SparseSolver& rhs )\n");
throw; // TODO: More elegant exit?
}
if ( this != ma27_rhs )
{
clear( );
SparseSolver::operator=( rhs );
copy( *ma27_rhs );
}
return *this;
}
/*
* s e t M a t r i x D a t a
*/
returnValue Ma27SparseSolver::setMatrixData( int_t dim_,
int_t numNonzeros_,
const int_t* const irn,
const int_t* const jcn,
const real_t* const avals
)
{
reset( );
dim = dim_;
numNonzeros = numNonzeros_;
if ( numNonzeros_ > 0 )
{
a_ma27 = new double[numNonzeros_];
irn_ma27 = new fint[numNonzeros_];
jcn_ma27 = new fint[numNonzeros_];
numNonzeros=0;
for (int_t i=0; i<numNonzeros_; ++i)
if ( avals[i] != 0 )
{
a_ma27[numNonzeros] = avals[i];
irn_ma27[numNonzeros] = irn[i];
jcn_ma27[numNonzeros] = jcn[i];
numNonzeros++;
}
}
else
{
numNonzeros = 0;
a_ma27 = 0;
irn_ma27 = 0;
jcn_ma27 = 0;
}
return SUCCESSFUL_RETURN;
}
/*
* f a c t o r i z e
*/
returnValue Ma27SparseSolver::factorize( )
{
if ( dim == 0 )
{
have_factorization = true;
neig = 0;
rank = 0;
return SUCCESSFUL_RETURN;
}
/******************************************
* Call MA27AD for symbolic factorization *
******************************************/
// Overstimation factor for LIW (20% recommended in MA27 documentation)
const double LiwFact = 2.0; // This is 200% overestimation
liw_ma27 = (fint)(LiwFact*(double(2*numNonzeros+3*dim+1)));
iw_ma27 = new fint[liw_ma27];
ikeep_ma27 = new fint[3*dim];
fint iflag_ma27 = 0;
double ops_ma27;
fint info_ma27[20];
fint* iw1_ma27 = new fint[2*dim];
MA27AD(&dim, &numNonzeros, irn_ma27, jcn_ma27, iw_ma27, &liw_ma27, ikeep_ma27,
iw1_ma27, &nsteps_ma27, &iflag_ma27, icntl_ma27, cntl_ma27,
info_ma27, &ops_ma27);
/* Receive some information from MA27AD */
fint iflag = info_ma27[0]; // Information flag
fint ierror = info_ma27[1]; // Error flag
fint nrlnec = info_ma27[4]; // recommended value for la
fint nirnec = info_ma27[5]; // recommended value for liw
if (iflag != 0)
{
MyPrintf("MA27AD returns iflag = %d with ierror = %d\n", iflag, ierror);
delete [] iw1_ma27;
clear( );
return THROWERROR(RET_MATRIX_FACTORISATION_FAILED);
}
/* Allocate memory for actual factorization */
delete [] iw_ma27;
double liw_init_factor = 5.0; // This could be an option.
liw_ma27 = (fint)(liw_init_factor * (double)(nirnec));
iw_ma27 = new fint[liw_ma27];
double la_init_factor = 20.0; // This could be an option.
la_ma27 = getMax(numNonzeros,(fint)(la_init_factor * (double)(nrlnec)));
double* a_new = new double[la_ma27];
for (int_t i=0; i<numNonzeros; ++i)
a_new[i] = a_ma27[i];
delete [] a_ma27;
a_ma27 = a_new;
/*******************************************
* Call MA27BD for numerical factorization *
*******************************************/
MA27BD(&dim, &numNonzeros, irn_ma27, jcn_ma27, a_ma27,
&la_ma27, iw_ma27, &liw_ma27, ikeep_ma27, &nsteps_ma27,
&maxfrt_ma27, iw1_ma27, icntl_ma27, cntl_ma27, info_ma27);
delete [] iw1_ma27;
/* Receive some information from MA27BD */
iflag = info_ma27[0]; // Information flag
ierror = info_ma27[1]; // Error flag
neig = info_ma27[14]; // Number of negative eigenvalues
if (iflag == 3)
{
rank = info_ma27[1];
return RET_KKT_MATRIX_SINGULAR;
}
else if (iflag == -5)
{ //DJ: I think this is more severe. Can this actually happen?
rank = -1;
return RET_KKT_MATRIX_SINGULAR;
}
else if (iflag != 0)
{
MyPrintf("MA27BD returns iflag = %d with ierror = %d\n", iflag, ierror);
clear( );
return THROWERROR(RET_MATRIX_FACTORISATION_FAILED);
}
else
rank = dim;
have_factorization = true;
return SUCCESSFUL_RETURN;
}
/*
* s o l v e
*/
returnValue Ma27SparseSolver::solve( int_t dim_,
const real_t* const rhs,
real_t* const sol
)
{
/* consistency check */
if ( dim_ != dim )
return THROWERROR( RET_INVALID_ARGUMENTS );
if ( !have_factorization )
{
MyPrintf("Factorization not called before solve in Ma27SparseSolver::solve.\n");
return THROWERROR( RET_INVALID_ARGUMENTS );
}
if ( dim == 0 )
return SUCCESSFUL_RETURN;
/* Call MA27CD to solve the system */
double* w_ma27 = new double[maxfrt_ma27];
fint* iw1_ma27 = new fint[nsteps_ma27];
/* MA27CD overwrites rhs */
for (int_t i=0; i<dim; ++i) sol[i] = rhs[i];
MA27CD(&dim, a_ma27, &la_ma27, iw_ma27, &liw_ma27, w_ma27, &maxfrt_ma27,
sol, iw1_ma27, &nsteps_ma27, icntl_ma27, cntl_ma27);
delete [] iw1_ma27;
delete [] w_ma27;
return SUCCESSFUL_RETURN;
}
/*
* r e s e t
*/
returnValue Ma27SparseSolver::reset( )
{
/* AW: We probably want to avoid resetting factorization in QProblem */
if ( SparseSolver::reset( ) != SUCCESSFUL_RETURN )
return THROWERROR( RET_RESET_FAILED );
clear( );
return SUCCESSFUL_RETURN;
}
/*
* g e t N e g a t i v e E i g e n v a l u e s
*/
int_t Ma27SparseSolver::getNegativeEigenvalues( )
{
if( !have_factorization )
return -1;
else
return neig;
}
/*
* g e t R a n k
*/
int_t Ma27SparseSolver::getRank( )
{
return rank;
}
/*****************************************************************************
* P R O T E C T E D *
*****************************************************************************/
/*
* c l e a r
*/
returnValue Ma27SparseSolver::clear( )
{
delete [] a_ma27;
delete [] irn_ma27;
delete [] jcn_ma27;
delete [] iw_ma27;
delete [] ikeep_ma27;
dim = -1;
numNonzeros = -1;
neig = -1;
rank = -1;
la_ma27 = -1;
a_ma27 = 0;
irn_ma27 = 0;
jcn_ma27 = 0;
liw_ma27 = -1;
iw_ma27 = 0;
ikeep_ma27 = 0;
nsteps_ma27 = -1;
maxfrt_ma27 = -1;
have_factorization = false;
return SUCCESSFUL_RETURN;
}
/*
* c o p y
*/
returnValue Ma27SparseSolver::copy( const Ma27SparseSolver& rhs
)
{
dim = rhs.dim;
numNonzeros = rhs.numNonzeros;
la_ma27 = rhs.la_ma27;
if ( rhs.a_ma27 != 0 )
{
if (rhs.have_factorization)
{
a_ma27 = new double[la_ma27];
memcpy( a_ma27,rhs.a_ma27,la_ma27*sizeof(double) );
}
else
{
a_ma27 = new double[numNonzeros];
memcpy( a_ma27,rhs.a_ma27,numNonzeros*sizeof(double) );
}
}
else
a_ma27 = 0;
if ( rhs.irn_ma27 != 0 )
{
irn_ma27 = new fint[numNonzeros];
memcpy( irn_ma27,rhs.irn_ma27,numNonzeros*sizeof(fint) );
}
else
irn_ma27 = 0;
if ( rhs.jcn_ma27 != 0 )
{
jcn_ma27 = new fint[numNonzeros];
memcpy( jcn_ma27,rhs.jcn_ma27,numNonzeros*sizeof(fint) );
}
else
jcn_ma27 = 0;
for ( int_t i=0; i<30; ++i)
icntl_ma27[i] = rhs.icntl_ma27[i];
for ( int_t i=0; i<5; ++i)
cntl_ma27[i] = rhs.cntl_ma27[i];
liw_ma27 = rhs.liw_ma27;
if ( rhs.iw_ma27 != 0 )
{
iw_ma27 = new fint[liw_ma27];
memcpy( iw_ma27,rhs.iw_ma27,liw_ma27*sizeof(fint) );
}
else
iw_ma27 = 0;
if ( rhs.ikeep_ma27 != 0 )
{
ikeep_ma27 = new fint[3*dim];
memcpy( ikeep_ma27,rhs.ikeep_ma27,3*dim*sizeof(fint) );
}
else
ikeep_ma27 = 0;
nsteps_ma27 = rhs.nsteps_ma27;
maxfrt_ma27 = rhs.maxfrt_ma27;
have_factorization = rhs.have_factorization;
neig = rhs.neig;
rank = rhs.rank;
return SUCCESSFUL_RETURN;
}
#endif // SOLVER_MA27
#ifdef SOLVER_MA57
/*****************************************************************************
*****************************************************************************
*****************************************************************************
* M A 5 7 S P A R E S E S O L V E R *
*****************************************************************************
*****************************************************************************
*****************************************************************************/
#define MA57ID ma57id_
#define MA57AD ma57ad_
#define MA57BD ma57bd_
#define MA57CD ma57cd_
extern "C"
{
/*
* MA57ID -- Initialize solver.
*/
extern void MA57ID (
double *cntl,
fint *icntl);
/*
* MA57AD -- Symbolic Factorization.
*/
extern void MA57AD (
fint *n, /* Order of matrix. */
fint *ne, /* Number of entries. */
const fint *irn, /* Matrix nonzero row structure */
const fint *jcn, /* Matrix nonzero column structure */
fint *lkeep, /* Workspace for the pivot order of lenght 3*n */
fint *keep, /* Workspace for the pivot order of lenght 3*n */
/* Automatically iflag = 0; ikeep pivot order iflag = 1 */
fint *iwork, /* Integer work space. */
fint *icntl, /* Integer Control parameter of length 30*/
fint *info, /* Statistical Information; Integer array of length 20 */
double *rinfo);/* Double Control parameter of length 5 */
/*
* MA57BD -- Numerical Factorization.
*/
extern void MA57BD (
fint *n, /* Order of matrix. */
fint *ne, /* Number of entries. */
double *a, /* Numerical values. */
double *fact, /* Entries of factors. */
fint *lfact, /* Length of array `fact'. */
fint *ifact, /* Indexing info for factors. */
fint *lifact, /* Length of array `ifact'. */
fint *lkeep, /* Length of array `keep'. */
fint *keep, /* Integer array. */
fint *iwork, /* Workspace of length `n'. */
fint *icntl, /* Integer Control parameter of length 20. */
double *cntl, /* Double Control parameter of length 5. */
fint *info, /* Statistical Information; Integer array of length 40. */
double *rinfo); /* Statistical Information; Real array of length 20. */
/*
* MA57CD -- Solution.
*/
extern void MA57CD (
fint *job, /* Solution job. Solve for... */
/* JOB <= 1: A */
/* JOB == 2: PLP^t */
/* JOB == 3: PDP^t */
/* JOB >= 4: PL^t P^t */
fint *n, /* Order of matrix. */
double *fact, /* Entries of factors. */
fint *lfact, /* Length of array `fact'. */
fint *ifact, /* Indexing info for factors. */
fint *lifact, /* Length of array `ifact'. */
fint *nrhs, /* Number of right hand sides. */
double *rhs, /* Numerical Values. */
fint *lrhs, /* Leading dimensions of `rhs'. */
double *work, /* Real workspace. */
fint *lwork, /* Length of `work', >= N*NRHS. */
fint *iwork, /* Integer array of length `n'. */
fint *icntl, /* Integer Control parameter array of length 20. */
fint *info); /* Statistical Information; Integer array of length 40. */
}
/*****************************************************************************
* P U B L I C *
****************************************************************************/
/*
* S p a r s e S o l v e r B a s e
*/
Ma57SparseSolver::Ma57SparseSolver( ) : SparseSolver()
{
a_ma57 = 0;
irn_ma57 = 0;
jcn_ma57 = 0;
fact_ma57 = 0;
ifact_ma57 = 0;
pivots = 0;
clear( );
/* Set default options for MA57 */
MA57ID( cntl_ma57, icntl_ma57 );
icntl_ma57[0] = -1; // Suppress error messages
icntl_ma57[1] = -1; // Suppress warning messages
icntl_ma57[2] = -1; // Suppress monitoring messages
//icntl_ma57[4] = 4; // Print everything (for debugging)
icntl_ma57[15] = 1; // Place small pivots at the end of the factorization (default: 0)
/// \todo good default values?
//cntl_ma57[1] = 5.0e-16; // Pivots smaller than this are treated as zero and are placed at the end of the factorization (default: 1e-20)
//cntl_ma57[0] = 0.5; // Set pivot tolerance: Higher values = more stable but slower/less sparse (default: 0.01, max 0.5)
}
/*
* S p a r s e S o l v e r B a s e
*/
Ma57SparseSolver::Ma57SparseSolver( const Ma57SparseSolver& rhs )
{
copy( rhs );
}
/*
* ~ S p a r s e S o l v e r B a s e
*/
Ma57SparseSolver::~Ma57SparseSolver( )
{
clear( );
}
/*
* o p e r a t o r =
*/
Ma57SparseSolver& Ma57SparseSolver::operator=( const SparseSolver& rhs )
{
const Ma57SparseSolver* ma57_rhs = dynamic_cast<const Ma57SparseSolver*>(&rhs);
if (!ma57_rhs)
{
fprintf(getGlobalMessageHandler()->getOutputFile(),"Error in Ma57SparseSolver& Ma57SparseSolver::operator=( const SparseSolver& rhs )\n");
throw; // TODO: More elegant exit?
}
if ( this != ma57_rhs )
{
clear( );
SparseSolver::operator=( rhs );
copy( *ma57_rhs );
}
return *this;
}
/*
* s e t M a t r i x D a t a
*/
returnValue Ma57SparseSolver::setMatrixData( int_t dim_,
int_t numNonzeros_,
const int_t* const irn,
const int_t* const jcn,
const real_t* const avals
)
{
reset( );
dim = dim_;
numNonzeros = numNonzeros_;
if ( numNonzeros_ > 0 )
{
a_ma57 = new double[numNonzeros_];
irn_ma57 = new fint[numNonzeros_];
jcn_ma57 = new fint[numNonzeros_];
numNonzeros=0;
for (int_t i=0; i<numNonzeros_; ++i)
if ( isZero(avals[i]) == BT_FALSE )
{
a_ma57[numNonzeros] = avals[i];
irn_ma57[numNonzeros] = irn[i];
jcn_ma57[numNonzeros] = jcn[i];
numNonzeros++;
}
}
else
{
numNonzeros = 0;
a_ma57 = 0;
irn_ma57 = 0;
jcn_ma57 = 0;
}
return SUCCESSFUL_RETURN;
}
/*
* f a c t o r i z e
*/
returnValue Ma57SparseSolver::factorize( )
{
if ( dim == 0 )
{
have_factorization = true;
neig = 0;
rank = 0;
return SUCCESSFUL_RETURN;
}
/******************************************
* Call MA57AD for symbolic factorization *
******************************************/
fint lkeep_ma57 = 5*dim + numNonzeros + getMax(numNonzeros,dim) + 42;
fint *keep_ma57 = new fint[lkeep_ma57];
fint *iwork_ma57 = new fint[5*dim];
fint info_ma57[40];
double rinfo_ma57[20];
MA57AD(&dim, &numNonzeros, irn_ma57, jcn_ma57, &lkeep_ma57, keep_ma57,
iwork_ma57, icntl_ma57, info_ma57, rinfo_ma57);
/* Receive some information from MA57AD */
fint iflag = info_ma57[0]; // Information flag
fint ierror = info_ma57[1]; // Error flag
if (iflag != 0)
{
MyPrintf("MA57AD returns iflag = %d with ierror = %d\n", iflag, ierror);
delete [] keep_ma57;
delete [] iwork_ma57;
clear( );
return THROWERROR(RET_MATRIX_FACTORISATION_FAILED);
}
/* Allocate memory for actual factorization */
double lfact_factor = 10.0; // This could be an option
lfact_ma57 = (fint)(lfact_factor * (double)(info_ma57[8]));
fact_ma57 = new double[lfact_ma57];
lifact_ma57 = (fint)(lfact_factor * (double)(info_ma57[9]));
ifact_ma57 = new int_t[lifact_ma57];
/*******************************************
* Call MA57BD for numerical factorization *
*******************************************/
MA57BD( &dim, &numNonzeros, a_ma57, fact_ma57, &lfact_ma57,
ifact_ma57, &lifact_ma57, &lkeep_ma57, keep_ma57,
iwork_ma57, icntl_ma57, cntl_ma57, info_ma57, rinfo_ma57 );
delete [] iwork_ma57;
delete [] keep_ma57;
/* Receive some information from MA57BD */
iflag = info_ma57[0]; // Information flag
ierror = info_ma57[1]; // Error flag
neig = info_ma57[23]; // Number of negative eigenvalues
rank = info_ma57[24]; // Rank of matrix
/* Read pivot sequence (see MA57UD source code) */
pivots = new fint[dim];
fint nrows, ncols;
fint nblk = ifact_ma57[2];
int_t iwpos = 3; // = 4-1
int_t k, kk, count = 0;
for ( k=0; k<nblk; k++ )
{
ncols = ifact_ma57[iwpos];
nrows = ifact_ma57[iwpos+1];
for ( kk=0; kk<nrows; kk++ )
pivots[count++] = ifact_ma57[iwpos+2+kk]-1; //convert Fortran to C indices!
iwpos = iwpos+ncols+2;
}
if (iflag == 4)
{
//MyPrintf("dim = %i, rank = %i. Pivots: ", dim, rank);
//for( k=rank; k<dim; k++ )
//MyPrintf("%i ", pivots[k]);
//MyPrintf("\n");
return RET_KKT_MATRIX_SINGULAR;
}
else if (iflag != 0)
{
MyPrintf("MA57BD returns iflag = %d with ierror = %d\n", iflag, ierror);
clear( );
return THROWERROR(RET_MATRIX_FACTORISATION_FAILED);
}
have_factorization = true;
return SUCCESSFUL_RETURN;
}
/*
* s o l v e
*/
returnValue Ma57SparseSolver::solve( int_t dim_,
const real_t* const rhs,
real_t* const sol
)
{
/* consistency check */
if ( dim_ != dim )
return THROWERROR( RET_INVALID_ARGUMENTS );
if ( !have_factorization )
{
MyPrintf("Factorization not called before solve in Ma57SparseSolver::solve.\n");
return THROWERROR( RET_INVALID_ARGUMENTS );
}
if ( dim == 0 )
return SUCCESSFUL_RETURN;
/* Call MA57CD to solve the system */
fint job_ma57 = 1;
fint nrhs_ma57 = 1;
fint lrhs_ma57 = dim;
fint info_ma57[40];
fint lwork_ma57 = dim*nrhs_ma57;
double* work_ma57 = new double[lwork_ma57];
fint* iwork_ma57 = new fint[dim];
/* MA57CD overwrites rhs */
for (int_t i=0; i<dim; ++i) sol[i] = rhs[i];
MA57CD(&job_ma57, &dim, fact_ma57, &lfact_ma57, ifact_ma57, &lifact_ma57,
&nrhs_ma57, sol, &lrhs_ma57, work_ma57, &lwork_ma57, iwork_ma57,
icntl_ma57, info_ma57);
delete [] work_ma57;
delete [] iwork_ma57;
fint iflag = info_ma57[0]; // Information flag
fint ierror = info_ma57[1]; // Error flag
if (iflag != 0)
{
MyPrintf("MA57CD returns iflag = %d with ierror = %d\n", iflag, ierror);
clear( );
return THROWERROR(RET_MATRIX_FACTORISATION_FAILED);
}
return SUCCESSFUL_RETURN;
}
/*
* r e s e t
*/
returnValue Ma57SparseSolver::reset( )
{
/* AW: We probably want to avoid resetting factorization in QProblem */
if ( SparseSolver::reset( ) != SUCCESSFUL_RETURN )
return THROWERROR( RET_RESET_FAILED );
clear( );
return SUCCESSFUL_RETURN;
}
/*
* g e t N e g a t i v e E i g e n v a l u e s
*/
int_t Ma57SparseSolver::getNegativeEigenvalues( )
{
if( !have_factorization )
return -1;
else
return neig;
}
/*
* g e t R a n k
*/
int_t Ma57SparseSolver::getRank( )
{
return rank;
}
/*
* g e t Z e r o P i v o t s
*/
returnValue Ma57SparseSolver::getZeroPivots( int_t *&zeroPivots )
{
for ( int_t k=0; k<dim-rank; k++ )
zeroPivots[k] = pivots[rank+k];
return SUCCESSFUL_RETURN;
}
/*****************************************************************************
* P R O T E C T E D *
*****************************************************************************/
/*
* c l e a r
*/
returnValue Ma57SparseSolver::clear( )
{
delete [] a_ma57;
delete [] irn_ma57;
delete [] jcn_ma57;
delete [] fact_ma57;
delete [] ifact_ma57;
delete [] pivots;
dim = -1;
numNonzeros = -1;
neig = -1;
rank = -1;
pivots = 0;
a_ma57 = 0;
irn_ma57 = 0;
jcn_ma57 = 0;
fact_ma57 = 0;
lfact_ma57 = -1;
ifact_ma57 = 0;
lifact_ma57 = -1;
have_factorization = false;
return SUCCESSFUL_RETURN;
}
/*
* c o p y
*/
returnValue Ma57SparseSolver::copy( const Ma57SparseSolver& rhs