-
Notifications
You must be signed in to change notification settings - Fork 40
/
main_iter.c
3407 lines (2901 loc) · 155 KB
/
main_iter.c
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
/*******************************************************************************
* Copyright (c) 2018, College of William & Mary
* All rights reserved.
*
* Redistribution and use in source and binary forms, with or without
* modification, are permitted provided that the following conditions are met:
* * Redistributions of source code must retain the above copyright
* notice, this list of conditions and the following disclaimer.
* * Redistributions in binary form must reproduce the above copyright
* notice, this list of conditions and the following disclaimer in the
* documentation and/or other materials provided with the distribution.
* * Neither the name of the College of William & Mary nor the
* names of its contributors may be used to endorse or promote products
* derived from this software without specific prior written permission.
*
* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
* ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
* WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
* DISCLAIMED. IN NO EVENT SHALL THE COLLEGE OF WILLIAM & MARY BE LIABLE FOR ANY
* DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
* (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
* LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
* ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
* (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
* SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*
* PRIMME: https://github.com/primme/primme
* Contact: Andreas Stathopoulos, a n d r e a s _at_ c s . w m . e d u
*******************************************************************************
* File: main_iter.c
*
* Purpose - This is the main Davidson-type iteration
*
******************************************************************************/
#ifndef THIS_FILE
#define THIS_FILE "../eigs/main_iter.c"
#endif
#include "numerical.h"
#include "template_normal.h"
#include "common_eigs.h"
/* Keep automatically generated headers under this section */
#ifndef CHECK_TEMPLATE
#include "main_iter.h"
#include "convergence.h"
#include "correction.h"
#include "init.h"
#include "ortho.h"
#include "restart.h"
#include "solve_projection.h"
#include "update_projection.h"
#include "update_W.h"
#include "auxiliary_eigs.h"
#include "auxiliary_eigs_normal.h"
#include "sketch.h"
#endif
#ifdef SUPPORTED_TYPE
/*----------------------------------------------------------------------------*
* The following are needed for the Dynamic Method Switching
*----------------------------------------------------------------------------*/
#ifndef MAIN_ITER_PRIVATE_H
#define MAIN_ITER_PRIVATE_H
typedef struct {
/* Time measurements for various components of the solver */
double MV_PR; /* OPp operator MV+PR time */
double MV; /* MV time only */
double PR; /* PRecond only */
double qmr_only; /* c_q only the QMR iteration */
double qmr_plus_MV_PR; /* a QMR plus operators */
double gdk_plus_MV_PR; /* b GD plus operators */
double gdk_plus_MV; /* b_m GD plus MV (i.e., GD without correction) */
/* The following two are not currently updated: */
double project_locked; /* c_p projection time per locked eigenvector in QMR*/
double reortho_locked; /* c_g (usu 2*c_p) ortho per locked vector in outer */
/* Average convergence estimates. Updated at restart/switch/convergence */
double gdk_conv_rate; /* convergence rate of all (|r|/|r0|) seen for GD+k */
double jdq_conv_rate; /* convergence rate of all (|r|/|r0|) seen for JDQMR*/
double JDQMR_slowdown; /* log(gdRate)/log(jdRate) restricted in [1.1, 2.5] */
double ratio_MV_outer; /* TotalMV/outerIts (average of last two updates) */
/* To maintain a convergence rate of all residual reductions seen we need: */
int nextReset; /* When to reset the following averaging sums */
double gdk_sum_logResReductions;/* Sumof all log(residual reductions) in GD*/
double jdq_sum_logResReductions;/* Sumof all log(residual reductions) in JD*/
double gdk_sum_MV; /* Total num of MV resulting in these GD reductions */
double jdq_sum_MV; /* Total num of MV resulting in these JD reductions */
int nevals_by_gdk; /* Number of evals found by GD+k since last reset */
int nevals_by_jdq; /* Number of evals found by JDQMR since last reset */
/* Variables to remember MV/its/time/resNorm/etc since last update/switch */
int numIt_0; /*Remembers starting outer its/MVs since switched */
int numMV_0; /* to current method, or since an epair converged */
double timer_0; /*Remembers starting time since switched to a method*/
/* or since an epair converged with that method */
double time_in_inner; /*Accumulative time spent in inner iterations */
/* since last switch or since an epair converged */
double resid_0; /*First residual norm of the convergence of a method*/
/* since last switch or since an epair converged */
/* Weighted ratio of expected times, for final method recommendation */
double accum_jdq_gdk; /*Expected ratio of accumulative times of JDQMR/GD+k*/
double accum_jdq; /* Accumulates jdq_times += ratio*(gdk+MV+PR) */
double accum_gdk; /* Accumulates gdk_times += gdk+MV+PR */
} primme_CostModel;
#endif /* MAIN_ITER_PRIVATE_H */
#if 0
STATIC void displayModel(primme_CostModel *model);
#endif
TEMPLATE_PLEASE
int print_timings_Sprimme(PRIMME_INT basisSize, primme_context ctx) {
primme_params *primme = ctx.primme;
printf("Basis size: %-" PRIMME_INT_P "\n", basisSize);
printf("Iterations: %-" PRIMME_INT_P "\n", primme->stats.numOuterIterations);
printf("Restarts : %-" PRIMME_INT_P "\n", primme->stats.numRestarts);
printf("Matvecs : %-" PRIMME_INT_P "\n", primme->stats.numMatvecs);
printf("Preconds : %-" PRIMME_INT_P "\n", primme->stats.numPreconds);
printf("Elapsed Time : %-22.10E\n", primme->stats.elapsedTime);
printf("MatVec Time : %-22.10E\n", primme->stats.timeMatvec);
printf("Precond Time : %-22.10E\n", primme->stats.timePrecond);
printf("Ortho Time : %-22.10E\n", primme->stats.timeOrtho);
printf("GlobalSum Time : %-22.10E\n", primme->stats.timeGlobalSum);
printf("Broadcast Time : %-22.10E\n", primme->stats.timeBroadcast);
printf("SketchedMatvec Time : %-22.10E\n", primme->stats.timeSketchMatvec);
printf("Rayleigh-Ritz Time : %-22.10E\n", primme->stats.timeRR);
printf("Stabilization Time : %-22.10E\n", primme->stats.timeStabilization);
printf("Residual Time : %-22.10E\n", primme->stats.timeResiduals);
return 0;
}
/******************************************************************************
* Subroutine main_iter - This routine implements a more general, parallel,
* block (Jacobi)-Davidson outer iteration with a variety of options.
*
* A coarse outline of the algorithm performed is as follows:
*
* 1. Initialize basis V
* 2. while (converged Ritz vectors have become unconverged) do
* 3. Update W=A*V, and H=V'*A*V, and solve the H eigenproblem
* 4. while (not all Ritz pairs have converged) do
* 5. while (maxBasisSize has not been reached)
* 6. Adjust the block size if necessary
* 7. Compute the residual norms of each block vector and
* check for convergence. If a Ritz vector has converged,
* target an unconverged Ritz vector to replace it.
* 8. Solve the correction equation for each residual vector
* 9. Insert corrections in V, orthogonalize, update W and H
* 10. Solve the new H eigenproblem to obtain the next Ritz pairs
* 11. endwhile
* 12. Restart V with appropriate vectors
* 13. If (locking) Lock out Ritz pairs that converged since last restart
* 14. endwhile
* 15. if Ritz vectors have become unconverged, reset convergence flags
* 16. endwhile
*
*
* OUTPUT arrays and parameters
* ----------------------------
* evals The converged Ritz values. It is maintained in sorted order.
* If locking is not engaged, values are copied to this array only
* upon return. If locking is engaged, then values are copied to this
* array as they converge.
*
* resNorms The residual norms corresponding to the converged Ritz vectors
*
* INPUT/OUTPUT arrays and parameters
* ----------------------------------
* evecs Stores initial guesses. Upon return, it contains the converged Ritz
* vectors. If locking is not engaged, then converged Ritz vectors
* are copied to this array just before return.
*
* primme.initSize: On output, it stores the number of converged eigenvectors.
* If smaller than numEvals and locking is used, there are
* only primme.initSize vectors in evecs.
* Without locking all numEvals approximations are in evecs
* but only the initSize ones are converged.
* During the execution, access to primme.initSize gives
* the number of converged pairs up to that point. The pairs
* are available in evals and evecs, but only when locking is used
* ret successful state
* numRet number of returned pairs in evals, evecs, and resNorms
*
* Return Value
* ------------
* error code
******************************************************************************/
TEMPLATE_PLEASE
int main_iter_Sprimme(HEVAL *evals, SCALAR *evecs, PRIMME_INT ldevecs,
HREAL *resNorms, double startTime, int *ret, int *numRet,
primme_context ctx) {
/* Default error is something is wrong in this function */
*ret = PRIMME_MAIN_ITER_FAILURE;
primme_params *primme = ctx.primme;
double elapsed_time = primme_wTimer();
double solve_timer;
primme->stats.elapsedTime = 0.0; // XXX: Added this for debugging purposes - Heather
/* primme parameters */
int i; /* Loop variable */
int blockSize; /* Current block size */
int availableBlockSize; /* There is enough space in basis for this block */
int basisSize; /* Current size of the basis V */
int numLocked; /* Number of locked Ritz vectors */
int numGuesses; /* Current number of initial guesses in evecs */
int nextGuess; /* Index of the next initial guess in evecs */
int numConverged; /* Number of converged Ritz pairs */
int targetShiftIndex; /* Target shift used in the QR factorization */
int recentlyConverged; /* Number of target Ritz pairs that have */
/* converged during the current iteration. */
int maxRecentlyConverged;/* Maximum value allowed for recentlyConverged */
int numConvergedStored; /* Numb of Ritzvecs temporarily stored in evecs */
/* to allow for skew projectors w/o locking */
int LockingProblem; /* Flag==1 if practically converged pairs locked */
int restartLimitReached; /* True when maximum restarts performed */
int nprevhVecs; /* Number of vectors stored in prevhVecs */
int numArbitraryVecs; /* Columns in hVecs computed with RR instead of */
/* the current extraction method. */
int maxEvecsSize; /* Maximum capacity of evecs array */
int numPrevRitzVals = 0; /* Size of the prevRitzVals updated in correction*/
int touch=0; /* param used in inner solver stopping criteria */
int *perm=NULL; /* Permutation on the locked pairs */
int *flags; /* Indicates which Ritz values have converged */
int *map; /* Indices of the closest vector from prevhVecs */
int *lockedFlags=NULL; /* Flags for the locked pairs */
int *ipivot; /* The pivot for the Mfact factorization of M */
int *iev; /* Evalue index each block vector corresponds to */
SCALAR *V; /* Basis vectors */
PRIMME_INT ldV; /* The leading dimension of V */
SCALAR *W; /* Work space storing A*V */
PRIMME_INT ldW; /* The leading dimension of W */
SCALAR *BV=NULL; /* Work space storing B*V */
PRIMME_INT ldBV; /* The leading dimension of BV */
HSCALAR *H; /* Upper triangular portion of V'*A*V */
HSCALAR *VtBV = NULL; /* Upper triangular portion of V'*B*V */
HSCALAR *fVtBV = NULL; /* Cholesky factor of VtBV */
HSCALAR *QtQ = NULL; /* Upper triangular portion of Q'*Q */
HSCALAR *fQtQ = NULL; /* Cholesky factor of QtQ */
HSCALAR *M = NULL; /* The projection Q'*K*B*Q, where Q = [evecs, x] */
/* x is the current Ritz vector and K is a */
/* hermitian preconditioner. */
HSCALAR *Mfact = NULL; /* The factorization of M=Q'KBQ */
SCALAR *Bevecs = NULL; /* B*evecs */
PRIMME_INT ldBevecs; /* Leading dimension of Bevecs */
SCALAR *evecsHat = NULL; /* K^{-1}*B*evecs */
PRIMME_INT ldevecsHat=0; /* The leading dimension of evecsHat */
HSCALAR *hVecs; /* Eigenvectors of H */
HSCALAR *hU = NULL; /* Left singular vectors of R */
HSCALAR *prevhVecs=NULL; /* hVecs from previous iteration */
int numQR; /* Maximum number of QR factorizations */
SCALAR *Q = NULL; /* QR decompositions for harmonic or refined */
PRIMME_INT ldQ; /* The leading dimension of Q */
HSCALAR *R = NULL; /* projection: (A-target[i])*V = QR */
HSCALAR *QtV = NULL; /* Q'*V */
HSCALAR *hVecsRot = NULL; /* transformation of hVecs in arbitrary vectors */
HEVAL *hVals; /* Eigenvalues of H */
HREAL *hSVals = NULL; /* Singular values of R */
HEVAL *prevRitzVals; /* Eigenvalues of H at previous outer iteration */
/* by robust shifting algorithm in correction.c */
HREAL *basisNorms; /* Residual norms of basis at pairs */
HREAL *blockNorms; /* Residual norms corresponding to current block */
/* vectors. */
double smallestResNorm; /* the smallest residual norm in the block */
int reset=0; /* Flag to reset V and W */
int restartsSinceReset=0;/* Restart since last reset of V and W */
int wholeSpace=0; /* search subspace reach max size */
/* Runtime measurement variables for dynamic method switching */
primme_CostModel CostModel; /* Structure holding the runtime estimates of */
/* the parameters of the model.Only visible here */
double tstart=0.0; /* Timing variable for accumulative time spent */
int maxNumRandoms = 10; /* We do not allow more than 10 randomizations */
/* -------------------------------------------------------------- */
/* Allocate objects */
/* -------------------------------------------------------------- */
maxEvecsSize = primme->numOrthoConst + primme->numEvals;
if (primme->projectionParams.projection != primme_proj_RR) {
numQR = 1;
}
else {
numQR = 0;
}
/* Use leading dimension ldOPs for the large dimension mats: V, W and Q */
ldV = ldW = ldBV = ldQ = primme->ldOPs;
ldBevecs = primme->massMatrixMatvec ? primme->ldOPs : ldevecs;
if (primme->massMatrixMatvec) {
CHKERR(Num_malloc_Sprimme(ldBV*primme->maxBasisSize, &BV, ctx));
}
CHKERR(Num_malloc_Sprimme(primme->ldOPs*primme->maxBasisSize, &V, ctx));
CHKERR(Num_malloc_Sprimme(primme->ldOPs*primme->maxBasisSize, &W, ctx));
if (numQR > 0) {
CHKERR(Num_malloc_Sprimme(primme->ldOPs*primme->maxBasisSize*numQR, &Q, ctx));
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize * numQR, &R, ctx));
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize * numQR, &hU, ctx));
}
if (primme->projectionParams.projection == primme_proj_harmonic) {
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize * numQR, &QtV, ctx));
}
CHKERR(Num_malloc_SHprimme(primme->maxBasisSize * primme->maxBasisSize, &H, ctx));
CHKERR(Num_malloc_SHprimme(primme->maxBasisSize * primme->maxBasisSize,
&hVecs, ctx));
int maxRank; /* maximum size of the main space being orthonormalize */
if (primme->locking) {
maxRank = primme->numOrthoConst + primme->numEvals + primme->maxBasisSize;
} else {
maxRank = primme->numOrthoConst + primme->maxBasisSize;
}
if (primme->orth == primme_orth_explicit_I) {
CHKERR(Num_malloc_SHprimme(maxRank * maxRank, &VtBV, ctx));
CHKERR(Num_malloc_SHprimme(maxRank * maxRank, &fVtBV, ctx));
}
int ldVtBV = maxRank;
int ldfVtBV = maxRank;
if (primme->orth == primme_orth_explicit_I && numQR) {
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize * numQR, &QtQ, ctx));
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize * numQR, &fQtQ, ctx));
}
int ldQtQ = primme->maxBasisSize;
int ldfQtQ = primme->maxBasisSize;
if (primme->projectionParams.projection == primme_proj_refined
|| primme->projectionParams.projection == primme_proj_harmonic) {
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize * numQR, &hVecsRot, ctx));
}
if (primme->correctionParams.precondition &&
primme->correctionParams.maxInnerIterations != 0 &&
primme->correctionParams.projectors.RightQ &&
primme->correctionParams.projectors.SkewQ ) {
ldevecsHat = primme->nLocal;
CHKERR(Num_malloc_Sprimme(ldevecsHat * maxEvecsSize, &evecsHat, ctx));
CHKERR(Num_malloc_SHprimme(maxEvecsSize * maxEvecsSize, &M, ctx));
CHKERR(Num_malloc_SHprimme(maxEvecsSize * maxEvecsSize, &Mfact, ctx));
}
if (primme->massMatrixMatvec && primme->locking) {
CHKERR(Num_malloc_Sprimme(ldBevecs * maxEvecsSize, &Bevecs, ctx));
}
CHKERR(Num_malloc_SHprimme(
primme->maxBasisSize * primme->maxBasisSize, &prevhVecs, ctx));
CHKERR(KIND(Num_malloc_RHprimme, Num_malloc_SHprimme)(primme->maxBasisSize,
&hVals, ctx));
if (numQR > 0) {
CHKERR(Num_malloc_RHprimme(primme->maxBasisSize, &hSVals, ctx));
}
CHKERR(KIND(Num_malloc_RHprimme, Num_malloc_SHprimme)(
primme->maxBasisSize + primme->numEvals, &prevRitzVals, ctx));
CHKERR(Num_malloc_RHprimme(primme->maxBlockSize, &blockNorms, ctx));
CHKERR(Num_malloc_RHprimme(primme->maxBasisSize, &basisNorms, ctx));
CHKERR(
Num_zero_matrix_RHprimme(basisNorms, 1, primme->maxBasisSize, 1, ctx));
/* Integer objects */
CHKERR(Num_malloc_iprimme(primme->numEvals, &perm, ctx));
if (primme->locking) {
CHKERR(Num_malloc_iprimme(primme->numEvals, &lockedFlags, ctx));
}
CHKERR(Num_malloc_iprimme(primme->maxBasisSize, &flags, ctx));
CHKERR(Num_malloc_iprimme(primme->maxBasisSize, &map, ctx));
CHKERR(Num_malloc_iprimme(primme->maxBlockSize, &iev, ctx));
CHKERR(Num_malloc_iprimme(maxEvecsSize, &ipivot, ctx));
/* -------------------------------------------------------------- */
/* Initialize counters and flags */
/* -------------------------------------------------------------- */
primme->stats.numOuterIterations = 0;
primme->stats.numRestarts = 0;
primme->stats.numMatvecs = 0;
primme->stats.numPreconds = 0;
primme->stats.numGlobalSum = 0;
primme->stats.numBroadcast = 0;
primme->stats.volumeGlobalSum = 0;
primme->stats.volumeBroadcast = 0;
primme->stats.flopsDense = 0;
primme->stats.numOrthoInnerProds = 0.0;
primme->stats.elapsedTime = 0.0;
primme->stats.timeMatvec = 0.0;
primme->stats.timePrecond = 0.0;
primme->stats.timeOrtho = 0.0;
primme->stats.timeGlobalSum = 0.0;
primme->stats.timeBroadcast = 0.0;
primme->stats.timeRR = 0.0;
primme->stats.timeSketchMatvec = 0.0;
primme->stats.timeStabilization = 0.0;
primme->stats.timeResiduals = 0.0;
primme->stats.timeDense = 0.0;
primme->stats.estimateMinEVal = HUGE_VAL;
primme->stats.estimateMaxEVal = -HUGE_VAL;
primme->stats.estimateLargestSVal = -HUGE_VAL;
primme->stats.estimateBNorm = primme->massMatrixMatvec ? -HUGE_VAL : 1.0;
primme->stats.estimateInvBNorm = primme->massMatrixMatvec ? -HUGE_VAL : 1.0;
primme->stats.maxConvTol = 0.0;
primme->stats.estimateResidualError = 0.0;
primme->stats.lockingIssue = 0;
numLocked = 0;
LockingProblem = 0;
blockSize = 0;
*numRet = 0;
for (i=0; i<primme->numEvals; i++) perm[i] = i;
/* -------------------------------------- */
/* Quick return for matrix of dimension 1 */
/* -------------------------------------- */
if (primme->numEvals == 0) {
primme->initSize = 0;
*ret = 0;
goto clean;
}
/* -------------------- */
/* Initialize the basis */
/* -------------------- */
CHKERR(init_basis_Sprimme(V, primme->nLocal, ldV, W, ldW, BV, ldBV, evecs,
ldevecs, Bevecs, ldBevecs, evecsHat, primme->nLocal, M, maxEvecsSize,
Mfact, 0, ipivot, VtBV, ldVtBV, fVtBV, ldfVtBV, maxRank, &basisSize,
&nextGuess, &numGuesses, ctx));
/* Now initSize will store the number of converged pairs */
primme->initSize = 0;
/* ----------------------------------------------------------- */
/* Dynamic method switch means we need to decide whether to */
/* allow inner iterations based on runtime timing measurements */
/* ----------------------------------------------------------- */
if (primme->dynamicMethodSwitch > 0) {
initializeModel(&CostModel, primme);
CostModel.MV = primme->stats.timeMatvec/primme->stats.numMatvecs;
if (primme->numEvals < 5 ||
primme->maxBasisSize + (primme->locking ? primme->numEvals : 0) >=
primme->n)
primme->dynamicMethodSwitch = 1; /* Start tentatively GD+k */
else
primme->dynamicMethodSwitch = 3; /* Start GD+k for 1st pair */
primme->correctionParams.maxInnerIterations = 0;
}
/* ---------------------------------------------------------------------- */
/* Outer most loop */
/* Without locking, restarting can cause converged Ritz values to become */
/* unconverged. Keep performing JD iterations until they remain converged */
/* ---------------------------------------------------------------------- */
while (
primme->stats.numMatvecs < primme->maxMatvecs &&
(primme->maxOuterIterations == 0 ||
primme->stats.numOuterIterations < primme->maxOuterIterations)) {
if (reset > 0) PRINTF(5, "Resetting V, W and QR");
/* Reset convergence flags. This may only reoccur without locking */
primme->initSize = numConverged = numConvergedStored = numLocked;
reset = 0;
for (i=0; i<primme->maxBasisSize; i++)
flags[i] = UNCONVERGED;
/* Compute the initial H and solve for its eigenpairs */
targetShiftIndex = 0;
if (numQR) {
int nQ = 0;
CHKERR(update_Q_Sprimme(BV ? BV : V, primme->nLocal, ldBV, W, ldW, Q,
ldQ, R, primme->maxBasisSize, QtQ, ldQtQ, fQtQ, ldfQtQ,
primme->targetShifts[targetShiftIndex], 0, basisSize, &nQ, ctx));
CHKERRM(numQR && basisSize != nQ, -1, "Not supported deficient QR");
}
if (H)
CHKERR(update_projection_Sprimme(V, ldV, W, ldW, H,
primme->maxBasisSize, primme->nLocal, 0, basisSize,
KIND(1 /*symmetric*/, 0 /* unsymmetric */), ctx));
if (QtV) {
CHKERR(update_projection_Sprimme(
Q, ldQ, V, ldV, QtV, primme->maxBasisSize, primme->nLocal, 0,
basisSize, 0 /*unsymmetric*/, ctx));
}
solve_timer = primme_wTimer(); // XXX: Debugging
CHKERR(solve_H_SHprimme(H, basisSize, primme->maxBasisSize,
VtBV
? &VtBV[(primme->numOrthoConst + numLocked) * ldVtBV +
primme->numOrthoConst + numLocked]
: NULL,
ldVtBV, R, primme->maxBasisSize, QtV, primme->maxBasisSize, QtQ,
ldQtQ, hU, basisSize, hVecs, basisSize, hVals, hSVals,
numConverged, ctx));
primme->stats.timeRR += primme_wTimer() - solve_timer;
numArbitraryVecs = 0;
maxRecentlyConverged = availableBlockSize = blockSize = 0;
smallestResNorm = HUGE_VAL;
primme->stats.estimateResidualError = 0.0;
if (!primme->locking) primme->stats.maxConvTol = 0.0;
blockSize = 0;
restartsSinceReset = 0;
/* -------------------------------------------------------------- */
/* Begin the iterative process. Keep restarting until all of the */
/* required eigenpairs have been found (no verification) */
/* -------------------------------------------------------------- */
while (numConverged < primme->numEvals &&
primme->stats.numMatvecs < primme->maxMatvecs &&
(primme->maxOuterIterations == 0 ||
primme->stats.numOuterIterations <
primme->maxOuterIterations) &&
!wholeSpace) {
nprevhVecs = 0;
int candidates_prepared = 0; /* Has prepare_candidates been called? */
/* ----------------------------------------------------------------- */
/* Main block Davidson loop. */
/* Keep adding vectors to the basis V until the basis has reached */
/* maximum size or the basis plus the locked vectors span the entire */
/* space. Once this happens, restart with a smaller basis. */
/* ----------------------------------------------------------------- */
while (basisSize < primme->maxBasisSize &&
primme->stats.numMatvecs < primme->maxMatvecs &&
( primme->maxOuterIterations == 0 ||
primme->stats.numOuterIterations < primme->maxOuterIterations) ) {
primme->stats.numOuterIterations++;
// XXX: FOR TESTING
if(primme->procID == 0 && primme->stats.numOuterIterations % 100 == 0){
primme->stats.elapsedTime = primme_wTimer() - elapsed_time;
CHKERR(print_timings_Sprimme(basisSize, ctx));
}
/* When QR are computed and there are more than one target shift, */
/* limit blockSize and the converged values to one. */
if (primme->numTargetShifts > numConverged+1 && numQR) {
availableBlockSize = 1;
maxRecentlyConverged = numConverged-numLocked+1;
}
else {
availableBlockSize = primme->maxBlockSize;
maxRecentlyConverged = max(0, primme->numEvals-numConverged);
}
/* Limit blockSize to vacant vectors in the basis */
availableBlockSize = min(availableBlockSize, primme->maxBasisSize-basisSize);
/* Limit blockSize to remaining values to converge plus one */
availableBlockSize = min(availableBlockSize, maxRecentlyConverged+1);
/* Set the block with the first unconverged pairs */
if (availableBlockSize > 0) {
/* If locking with GD and no preconditioning is running, don't */
/* check practically convergence now, unless whole space. */
int practConvCheck = 0;
if (primme->n <= basisSize + numLocked + primme->numOrthoConst) {
practConvCheck = 1;
} else if (primme->locking && !primme->massMatrixMatvec &&
!primme->correctionParams.precondition &&
primme->correctionParams.maxInnerIterations == 0) {
practConvCheck = -1;
}
CHKERR(prepare_candidates(V, ldV, W, ldW, BV, ldBV,
primme->nLocal, H, primme->maxBasisSize, basisSize,
&V[basisSize * ldV], &W[basisSize * ldW],
BV ? &BV[basisSize * ldBV] : NULL,
1 /* compute approx vectors and residuals */, hVecs,
basisSize, hVals, hSVals, flags, maxRecentlyConverged,
blockNorms, blockSize, availableBlockSize, evecs,
numLocked, ldevecs, Bevecs, ldBevecs, evals, resNorms,
targetShiftIndex, iev, &blockSize, &recentlyConverged,
&numArbitraryVecs, &smallestResNorm, hVecsRot,
primme->maxBasisSize, numConverged, basisNorms, &reset,
VtBV, ldVtBV, prevhVecs, nprevhVecs, primme->maxBasisSize,
practConvCheck, map, startTime, ctx));
assert(recentlyConverged >= 0);
candidates_prepared = 1;
}
else {
blockSize = recentlyConverged = 0;
}
/* If the total number of converged pairs, including the */
/* recentlyConverged ones, are greater than or equal to the */
/* target number of eigenvalues, attempt to restart, verify */
/* their convergence, lock them if necessary, and return. */
/* For locking interior, restart and lock now any converged. */
/* If Q, restart after an eigenpair converged to recompute */
/* QR with a different shift. */
/* Also if it has been converged as many pairs as initial */
/* guesses has introduced in V, then restart and introduce */
/* new guesses. */
numConverged += recentlyConverged;
/* Report iteration */
CHKERR(monitorFun_Sprimme(hVals, basisSize, flags, iev, blockSize,
basisNorms, numConverged, evals, numLocked, lockedFlags,
resNorms, -1, -1.0, NULL, 0.0, primme_event_outer_iteration,
startTime, ctx));
/* Reset touch every time an eigenpair converges */
if (recentlyConverged > 0) touch = 0;
if (primme->dynamicMethodSwitch > 0) {
if (CostModel.resid_0 == -1.0L) /* remember the very */
CostModel.resid_0 = blockNorms[0]; /* first residual */
/* If some pairs converged OR we evaluate JDQMR at every step, */
/* update convergence statistics and consider switching */
if (recentlyConverged > 0 || primme->dynamicMethodSwitch == 2)
{
CostModel.MV =
primme->stats.timeMatvec/primme->stats.numMatvecs;
int ret0 = update_statistics(&CostModel, primme, tstart,
recentlyConverged, 0, numConverged, blockNorms[0]);
if (ret0) switch (primme->dynamicMethodSwitch) {
/* for few evals (dyn=1) evaluate GD+k only at restart*/
case 3:
CHKERR(switch_from_GDpk(&CostModel,ctx));
break;
case 2: case 4:
CHKERR(switch_from_JDQMR(&CostModel,ctx));
}
}
}
if (numConverged >= primme->numEvals ||
(primme->locking && numConverged > numLocked &&
primme->target != primme_smallest &&
primme->target != primme_largest &&
(primme->projectionParams.projection ==
primme_proj_RR ||
primme->target == primme_closest_geq ||
primme->target == primme_closest_leq)) ||
targetShiftIndex < 0 ||
(blockSize == 0 && recentlyConverged > 0) ||
/* NOTE: use the same condition as in restart_refined */
(numQR && fabs(primme->targetShifts[targetShiftIndex] -
primme->targetShifts[min(
primme->numTargetShifts - 1,
numConverged)]) >=
max(primme->aNorm,
primme->stats.estimateLargestSVal)) ||
(numConverged >= nextGuess - primme->numOrthoConst &&
numGuesses > 0)) {
break;
}
if (blockSize > 0) {
/* Solve the correction equations with the new blockSize Ritz */
/* vectors and residuals. */
/* If dynamic method switching, time the inner method */
if (primme->dynamicMethodSwitch > 0) {
tstart = primme_wTimer(); /* accumulate correction time */
}
CHKERR(solve_correction_Sprimme(V, ldV, W, ldW, BV ? BV : V,
ldBV, evecs, ldevecs, Bevecs ? Bevecs : evecs,
Bevecs ? ldBevecs : ldevecs, evecsHat, ldevecsHat, Mfact,
ipivot, evals, numLocked, numConvergedStored, hVals,
prevRitzVals, &numPrevRitzVals, flags, basisSize,
blockNorms, iev, blockSize, &touch, startTime, ctx));
/* ------------------------------------------------------ */
/* If dynamic method switch, accumulate inner method time */
/* ------------------------------------------------------ */
if (primme->dynamicMethodSwitch > 0)
CostModel.time_in_inner += primme_wTimer() - tstart;
} /* end of else blocksize=0 */
/* If locking with GD and no preconditioning is running, check */
/* practically convergence here. The heuristic needs Vlocked' */
/* times r, which is going to be store in Rlocked. */
HSCALAR *Rlocked = NULL;
int ldRlocked = primme->numOrthoConst + numLocked;
int blockSize0 = blockSize;
if (primme->locking && !primme->massMatrixMatvec &&
!primme->correctionParams.precondition &&
primme->correctionParams.maxInnerIterations == 0) {
CHKERR(
Num_malloc_SHprimme(ldRlocked * blockSize, &Rlocked, ctx));
}
/* Orthogonalize the corrections with respect to each other and */
/* the current basis. If basis hasn't expanded with any new */
/* vector, try a random vector */
for (i=0; i < maxNumRandoms; i++) {
int basisSizeOut;
CHKERR(Bortho_block_Sprimme(V, ldV, VtBV, ldVtBV, fVtBV, ldfVtBV,
NULL, 0, basisSize, basisSize + blockSize - 1, evecs,
ldevecs, primme->numOrthoConst + numLocked, BV, ldBV,
i == 0 ? Rlocked : NULL, ldRlocked, primme->nLocal,
maxRank, &basisSizeOut, ctx));
blockSize = basisSizeOut - basisSize;
if (blockSize > 0 || availableBlockSize <= 0) break;
PRINTF(5, "Warning: adding random vector "
"to the search subspace");
CHKERR(Num_larnv_Sprimme(2, primme->iseed, primme->nLocal,
&V[ldV * basisSize], ctx));
blockSize = 1;
}
/* The method failing in expanding the search subspace may be an */
/* indication that a) the orthogonalization error is to large for */
/* distinguishing the new components on the computed correction, */
/* or b) there is no more pairs that satisfy the convergence */
/* criterion. If prepare_candidates gave a candidate, we guess */
/* we are in case a), and otherwise we are in case b). */
if (i >= maxNumRandoms) {
if (availableBlockSize > 0 && blockSize0 <= 0 && reset == 0) {
wholeSpace = 1;
} else {
reset = 2;
}
blockSize = 0;
if (primme->locking && !primme->massMatrixMatvec &&
!primme->correctionParams.precondition &&
primme->correctionParams.maxInnerIterations == 0) {
CHKERR(Num_free_SHprimme(Rlocked, ctx));
}
break;
}
// If locking with GD and no preconditioning is running, check
// practically convergence here. A pair in the block is marked as
// practically converged if the next expression is small enough to
// pass the convergence criterion:
// sqrt(|r|^2 - |Rlocked|^2*(1 + |Xx|^2)),
// where Rlocked is V_locked'*r and Xx is V_locked'*x.
// The term Xx is not totally clear. It appears in the perturbation
// analysis of the eigenproblem with V_locked and x. And it is
// necessary to pass the next test for complex:
// tests/testi-100-LOBPCG_OrthoBasis_Window-100-primme_closest_abs-primme_proj_refined.F
if (primme->locking && !primme->massMatrixMatvec &&
!primme->correctionParams.precondition &&
primme->correctionParams.maxInnerIterations == 0) {
if (numLocked > 0) {
for (i = 0; i < blockSize0 && numConverged < primme->numEvals;
i++) {
HREAL normXx = 0.0;
if (primme->orth == primme_orth_explicit_I) {
HSCALAR *Xx = NULL;
CHKERR(Num_malloc_SHprimme(numLocked, &Xx, ctx));
CHKERR(Num_zero_matrix_SHprimme(
Xx, numLocked, 1, numLocked, ctx));
CHKERR(Num_gemv_SHprimme("N", numLocked, basisSize, 1.0,
&VtBV[ldVtBV * numLocked], ldVtBV,
&hVecs[iev[i] * basisSize], 1, 0.0, Xx, 1, ctx));
normXx = ABS(
Num_dot_SHprimme(numLocked, Xx, 1, Xx, 1, ctx));
CHKERR(Num_free_SHprimme(Xx, ctx));
}
HREAL normRlockedi = ABS(
Num_dot_SHprimme(ldRlocked, &Rlocked[ldRlocked * i],
1, &Rlocked[ldRlocked * i], 1, ctx));
HREAL newBlockNorm =
sqrt(max(blockNorms[i] * blockNorms[i] -
normRlockedi * (1. + normXx),
0.0));
CHKERR(check_convergence_Sprimme(&V[(basisSize + i) * ldV],
ldV, 1 /* given X */, NULL, 0, 0 /* not given R */,
evecs, numLocked, ldevecs, Bevecs, ldBevecs, VtBV,
ldVtBV, 0, 1, &flags[iev[i]], &newBlockNorm,
&hVals[iev[i]], &reset,
-1 /* don't check practically convergence */, ctx));
basisNorms[iev[i]] = newBlockNorm;
if (flags[iev[i]] == CONVERGED) {
flags[iev[i]] = PRACTICALLY_CONVERGED;
numConverged++;
/* Report a pair was soft converged */
CHKERR(monitorFun_Sprimme(hVals, basisSize, flags,
&iev[i], 1, basisNorms, numConverged, NULL, 0,
NULL, NULL, -1, -1.0, NULL, 0.0,
primme_event_converged, startTime, ctx));
}
}
}
CHKERR(Num_free_SHprimme(Rlocked, ctx));
if (numConverged > numLocked &&
primme->target != primme_smallest &&
primme->target != primme_largest &&
(primme->projectionParams.projection == primme_proj_RR ||
primme->target == primme_closest_geq ||
primme->target == primme_closest_leq)) {
break;
}
}
/* Compute W = A*V for the orthogonalized corrections */
CHKERR(matrixMatvec_Sprimme(V, primme->nLocal, ldV, W, ldW,
basisSize, blockSize, ctx));
if (numQR) {
int nQ = basisSize;
CHKERR(update_Q_Sprimme(BV ? BV : V, primme->nLocal, ldBV, W,
ldW, Q, ldQ, R, primme->maxBasisSize, QtQ, ldQtQ, fQtQ,
ldfQtQ, primme->targetShifts[targetShiftIndex], basisSize,
blockSize, &nQ, ctx));
if (basisSize + blockSize != nQ) {
blockSize = 0;
reset = 1;
break;
}
}
/* Extend H by blockSize columns and rows and solve the */
/* eigenproblem for the new H. */
if (H)
CHKERR(update_projection_Sprimme(V, ldV, W, ldW, H,
primme->maxBasisSize, primme->nLocal, basisSize, blockSize,
KIND(1 /*symmetric*/, 0 /* unsymmetric */), ctx));
if (QtV)
CHKERR(update_projection_Sprimme(
Q, ldQ, V, ldV, QtV, primme->maxBasisSize, primme->nLocal,
basisSize, blockSize, 0 /*unsymmetric*/, ctx));
/* Copy hVecs into prevhVecs */
CHKERR(Num_copy_matrix_SHprimme(hVecs, basisSize, basisSize,
basisSize, prevhVecs, primme->maxBasisSize, ctx));
CHKERR(Num_zero_matrix_SHprimme(&prevhVecs[basisSize],
primme->maxBasisSize - basisSize, basisSize,
primme->maxBasisSize, ctx));
nprevhVecs = basisSize;
basisSize += blockSize;
blockSize = 0;
solve_timer = primme_wTimer();
CHKERR(solve_H_SHprimme(H, basisSize, primme->maxBasisSize,
VtBV
? &VtBV[(primme->numOrthoConst + numLocked) * ldVtBV +
primme->numOrthoConst + numLocked]
: NULL,
ldVtBV, R, primme->maxBasisSize, QtV, primme->maxBasisSize,
QtQ, ldQtQ, hU, basisSize, hVecs, basisSize, hVals, hSVals,
numConverged, ctx));
primme->stats.timeRR += primme_wTimer() - solve_timer;
numArbitraryVecs = 0;
candidates_prepared = 0;
/* If QR decomposition accumulates so much error, force it to */
/* reset by setting targetShiftIndex to -1. We use the next */
/* heuristic. Note that if (s_0, u_0, y_0) is the smallest triplet*/
/* of R, (A-tau*I)*V = Q*R, and l_0 is the Rayleigh quotient of */
/* V*y_0, then */
/* |l_0-tau| = |y_0'*V'*(A-tau*I)*V*y_0| = |y_0'*V'*Q*R*y_0| = */
/* = |y_0'*V'*Q*u_0*s_0| <= s_0. */
/* So when |l_0-tau|-machEps*|A| > s_0, we consider to reset the */
/* QR factorization. machEps*|A| is the error computing l_0. */
/* NOTE: rarely observed |l_0-tau|-machEps*|A| slightly greater */
/* than s_0 after resetting. The condition restartsSinceReset > 0 */
/* avoids infinite loop in those cases. */
double eps_orth;
CHKERR(machineEpsOrth_Sprimme(&eps_orth, ctx));
if (primme->projectionParams.projection == primme_proj_refined &&
basisSize > 0 && restartsSinceReset > 1 &&
targetShiftIndex >= 0 &&
ABS((HSCALAR)primme->targetShifts[targetShiftIndex] -
hVals[0]) -
max(primme->aNorm,
primme->stats.estimateLargestSVal) *
eps_orth >
hSVals[0]) {
reset = 2;
PRINTF(5, "Detected some errors in QR");
break;
}
/* --------------------------------------------------------------- */
} /* while (basisSize<maxBasisSize && basisSize<n-orthoConst-numLocked)
* --------------------------------------------------------------- */
/* If wholeSpace, reset if the accumulated errors on the residual */
/* vectors are too large */
if (basisSize >= primme->n - primme->numOrthoConst - numLocked) {
if (primme->stats.maxConvTol < primme->stats.estimateResidualError)
reset = 1;
}
if (reset > 0) break;
/* ----------------------------------------------------------------- */
/* Restart basis will need the final coefficient vectors in hVecs */
/* to lock out converged vectors and to compute X and R for the next */
/* iteration. The function prepare_vecs will make sure that hVecs */
/* has proper coefficient vectors. */
/* ----------------------------------------------------------------- */
if (!candidates_prepared) {
/* -------------------------------------------------------------- */
/* TODO: merge all this logic in the regular main loop. After all */
/* restarting is a regular step that also shrinks the basis. */
/* -------------------------------------------------------------- */
int j,k,l,m;
double *dummySmallestResNorm, dummyZero = 0.0;
if (blockSize > 0) {
availableBlockSize = blockSize;
maxRecentlyConverged = 0;
}
/* When there are more than one target shift, */
/* limit blockSize and the converged values to one. */
else if (primme->numTargetShifts > numConverged+1) {
if (primme->locking) {
maxRecentlyConverged =
max(min(primme->numEvals, numLocked+1) - numConverged, 0);
}
else {
maxRecentlyConverged =
max(min(primme->numEvals, numConverged+1) - numConverged, 0);
}
availableBlockSize = maxRecentlyConverged;
}
else {
maxRecentlyConverged = max(0, primme->numEvals-numConverged);
/* Limit blockSize to vacant vectors in the basis */
availableBlockSize = max(0, min(primme->maxBlockSize, primme->maxBasisSize-(numConverged-numLocked)));
/* Limit blockSize to remaining values to converge plus one */
availableBlockSize = min(availableBlockSize, maxRecentlyConverged+1);
}
/* Limit basisSize to the matrix dimension */
availableBlockSize =
max(0, min(availableBlockSize,
primme->n - numLocked - primme->numOrthoConst));
// Calling prepare_candidates before restarting is compulsory when:
// - the order of the Ritz pairs depends on the residual vector
// norms (e.g. when using refine extraction, or target being