Actual source code: lmekrylov.c
slepc-3.17.1 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc matrix equation solver: "krylov"
13: Method: Arnoldi with Eiermann-Ernst restart
15: Algorithm:
17: Project the equation onto the Arnoldi basis and solve the compressed
18: equation the Hessenberg matrix H, restart by discarding the Krylov
19: basis but keeping H.
21: References:
23: [1] Y. Saad, "Numerical solution of large Lyapunov equations", in
24: Signal processing, scattering and operator theory, and numerical
25: methods, vol. 5 of Progr. Systems Control Theory, pages 503-511,
26: 1990.
28: [2] D. Kressner, "Memory-efficient Krylov subspace techniques for
29: solving large-scale Lyapunov equations", in 2008 IEEE Int. Conf.
30: Computer-Aided Control Systems, pages 613-618, 2008.
31: */
33: #include <slepc/private/lmeimpl.h>
34: #include <slepcblaslapack.h>
36: PetscErrorCode LMESetUp_Krylov(LME lme)
37: {
38: PetscInt N;
40: MatGetSize(lme->A,&N,NULL);
41: if (lme->ncv==PETSC_DEFAULT) lme->ncv = PetscMin(30,N);
42: if (lme->max_it==PETSC_DEFAULT) lme->max_it = 100;
43: LMEAllocateSolution(lme,1);
44: PetscFunctionReturn(0);
45: }
47: PetscErrorCode LMESolve_Krylov_Lyapunov_Vec(LME lme,Vec b,PetscBool fixed,PetscInt rrank,BV C1,BV *X1,PetscInt *col,PetscBool *fail,PetscInt *totalits)
48: {
49: PetscInt n=0,m,ldh,ldg=0,i,j,rank=0,lrank,pass,nouter=0,its;
50: PetscReal bnorm,beta,errest;
51: PetscBool breakdown;
52: PetscScalar *Harray,*G=NULL,*Gnew=NULL,*L,*U,*r,*Qarray,sone=1.0,zero=0.0;
53: PetscBLASInt n_,m_,rk_;
54: Mat Q,H;
56: *fail = PETSC_FALSE;
57: its = 0;
58: m = lme->ncv;
59: ldh = m+1;
60: MatCreateSeqDense(PETSC_COMM_SELF,ldh,m,NULL,&H);
61: MatDenseGetArray(H,&Harray);
63: VecNorm(b,NORM_2,&bnorm);
66: for (pass=0;pass<2;pass++) {
68: /* set initial vector to b/||b|| */
69: BVInsertVec(lme->V,0,b);
70: BVScaleColumn(lme->V,0,1.0/bnorm);
72: /* Restart loop */
73: while ((pass==0 && !*fail) || (pass==1 && its+1<nouter)) {
74: its++;
76: /* compute Arnoldi factorization */
77: BVMatArnoldi(lme->V,lme->A,H,0,&m,&beta,&breakdown);
78: BVSetActiveColumns(lme->V,0,m);
80: if (pass==0) {
81: /* glue together the previous H and the new H obtained with Arnoldi */
82: ldg = n+m+1;
83: PetscCalloc1(ldg*(n+m),&Gnew);
84: for (j=0;j<m;j++) PetscArraycpy(Gnew+n+(j+n)*ldg,Harray+j*ldh,m);
85: Gnew[n+m+(n+m-1)*ldg] = beta;
86: if (G) {
87: for (j=0;j<n;j++) PetscArraycpy(Gnew+j*ldg,G+j*(n+1),n+1);
88: PetscFree(G);
89: }
90: G = Gnew;
91: n += m;
92: } else {
93: /* update Z = Z + V(:,1:m)*Q with Q=U(blk,:)*P(1:nrk,:)' */
94: MatCreateSeqDense(PETSC_COMM_SELF,m,*col+rank,NULL,&Q);
95: MatDenseGetArray(Q,&Qarray);
96: PetscBLASIntCast(m,&m_);
97: PetscBLASIntCast(n,&n_);
98: PetscBLASIntCast(rank,&rk_);
99: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m_,&rk_,&rk_,&sone,U+its*m,&n_,L,&n_,&zero,Qarray+(*col)*m,&m_));
100: MatDenseRestoreArray(Q,&Qarray);
101: BVSetActiveColumns(*X1,*col,*col+rank);
102: BVMult(*X1,1.0,1.0,lme->V,Q);
103: MatDestroy(&Q);
104: }
106: if (pass==0) {
107: /* solve compressed Lyapunov equation */
108: PetscCalloc1(n,&r);
109: PetscCalloc1(n*n,&L);
110: r[0] = bnorm;
111: errest = PetscAbsScalar(G[n+(n-1)*ldg]);
112: LMEDenseHessLyapunovChol(lme,n,G,ldg,1,r,n,L,n,&errest);
113: LMEMonitor(lme,*totalits+its,errest);
114: PetscFree(r);
116: /* check convergence */
117: if (errest<lme->tol) {
118: lme->errest += errest;
119: PetscMalloc1(n*n,&U);
120: /* transpose L */
121: for (j=0;j<n;j++) {
122: for (i=j+1;i<n;i++) {
123: L[i+j*n] = PetscConj(L[j+i*n]);
124: L[j+i*n] = 0.0;
125: }
126: }
127: LMEDenseRankSVD(lme,n,L,n,U,n,&lrank);
128: PetscInfo(lme,"Rank of the Cholesky factor = %" PetscInt_FMT "\n",lrank);
129: nouter = its;
130: its = -1;
131: if (!fixed) { /* X1 was not set by user, allocate it with rank columns */
132: rank = lrank;
133: if (*col) BVResize(*X1,*col+rank,PETSC_TRUE);
134: else BVDuplicateResize(C1,rank,X1);
135: } else rank = PetscMin(lrank,rrank);
136: PetscFree(G);
137: break;
138: } else {
139: PetscFree(L);
140: if (*totalits+its>=lme->max_it) *fail = PETSC_TRUE;
141: }
142: }
144: /* restart with vector v_{m+1} */
145: if (!*fail) BVCopyColumn(lme->V,m,0);
146: }
147: }
149: *col += rank;
150: *totalits += its+1;
151: MatDenseRestoreArray(H,&Harray);
152: MatDestroy(&H);
153: if (L) PetscFree(L);
154: if (U) PetscFree(U);
155: PetscFunctionReturn(0);
156: }
158: PetscErrorCode LMESolve_Krylov_Lyapunov(LME lme)
159: {
160: PetscBool fail,fixed = lme->X? PETSC_TRUE: PETSC_FALSE;
161: PetscInt i,k,rank=0,col=0;
162: Vec b;
163: BV X1=NULL,C1;
164: Mat X1m,X1t,C1m;
166: MatLRCGetMats(lme->C,NULL,&C1m,NULL,NULL);
167: BVCreateFromMat(C1m,&C1);
168: BVSetFromOptions(C1);
169: BVGetActiveColumns(C1,NULL,&k);
170: if (fixed) {
171: MatLRCGetMats(lme->X,NULL,&X1m,NULL,NULL);
172: BVCreateFromMat(X1m,&X1);
173: BVSetFromOptions(X1);
174: BVGetActiveColumns(X1,NULL,&rank);
175: rank = rank/k;
176: }
177: for (i=0;i<k;i++) {
178: BVGetColumn(C1,i,&b);
179: LMESolve_Krylov_Lyapunov_Vec(lme,b,fixed,rank,C1,&X1,&col,&fail,&lme->its);
180: BVRestoreColumn(C1,i,&b);
181: if (fail) {
182: lme->reason = LME_DIVERGED_ITS;
183: break;
184: }
185: }
186: if (lme->reason==LME_CONVERGED_ITERATING) lme->reason = LME_CONVERGED_TOL;
187: BVCreateMat(X1,&X1t);
188: if (fixed) MatCopy(X1t,X1m,SAME_NONZERO_PATTERN);
189: else MatCreateLRC(NULL,X1t,NULL,NULL,&lme->X);
190: MatDestroy(&X1t);
191: BVDestroy(&C1);
192: BVDestroy(&X1);
193: PetscFunctionReturn(0);
194: }
196: SLEPC_EXTERN PetscErrorCode LMECreate_Krylov(LME lme)
197: {
198: lme->ops->solve[LME_LYAPUNOV] = LMESolve_Krylov_Lyapunov;
199: lme->ops->setup = LMESetUp_Krylov;
200: PetscFunctionReturn(0);
201: }