Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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"
12 :
13 : Method: Arnoldi with Eiermann-Ernst restart
14 :
15 : Algorithm:
16 :
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.
20 :
21 : References:
22 :
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.
27 :
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 : */
32 :
33 : #include <slepc/private/lmeimpl.h>
34 : #include <slepcblaslapack.h>
35 :
36 8 : static PetscErrorCode LMESetUp_Krylov(LME lme)
37 : {
38 8 : PetscInt N;
39 :
40 8 : PetscFunctionBegin;
41 8 : PetscCall(MatGetSize(lme->A,&N,NULL));
42 8 : if (lme->ncv==PETSC_DEFAULT) lme->ncv = PetscMin(30,N);
43 8 : if (lme->max_it==PETSC_DEFAULT) lme->max_it = 100;
44 8 : PetscCall(LMEAllocateSolution(lme,1));
45 8 : PetscFunctionReturn(PETSC_SUCCESS);
46 : }
47 :
48 57 : static PetscErrorCode LMESolve_Krylov_Lyapunov_Vec(LME lme,Vec b,PetscBool fixed,PetscInt rrank,BV C1,BV *X1,PetscInt *col,PetscBool *fail,PetscInt *totalits)
49 : {
50 57 : PetscInt n=0,m,ldh,ldg=0,i,j,rank=0,lrank,pass,nouter=0,its;
51 57 : PetscReal bnorm,beta,errest;
52 57 : PetscBool breakdown;
53 57 : PetscScalar *Harray,*G=NULL,*Gnew=NULL,*L=NULL,*U=NULL,*r,*Qarray,sone=1.0,zero=0.0;
54 57 : PetscBLASInt n_,m_,rk_;
55 57 : Mat Q,H;
56 :
57 57 : PetscFunctionBegin;
58 57 : *fail = PETSC_FALSE;
59 57 : its = 0;
60 57 : m = lme->ncv;
61 57 : ldh = m+1;
62 57 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ldh,m,NULL,&H));
63 57 : PetscCall(MatDenseGetArray(H,&Harray));
64 :
65 57 : PetscCall(VecNorm(b,NORM_2,&bnorm));
66 57 : PetscCheck(bnorm,PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_WRONG,"Cannot process a zero vector in the right-hand side");
67 :
68 171 : for (pass=0;pass<2;pass++) {
69 :
70 : /* set initial vector to b/||b|| */
71 114 : PetscCall(BVInsertVec(lme->V,0,b));
72 114 : PetscCall(BVScaleColumn(lme->V,0,1.0/bnorm));
73 :
74 : /* Restart loop */
75 183 : while ((pass==0 && !*fail) || (pass==1 && its+1<nouter)) {
76 126 : its++;
77 :
78 : /* compute Arnoldi factorization */
79 126 : PetscCall(BVMatArnoldi(lme->V,lme->A,H,0,&m,&beta,&breakdown));
80 126 : PetscCall(BVSetActiveColumns(lme->V,0,m));
81 :
82 126 : if (pass==0) {
83 : /* glue together the previous H and the new H obtained with Arnoldi */
84 63 : ldg = n+m+1;
85 63 : PetscCall(PetscCalloc1(ldg*(n+m),&Gnew));
86 1845 : for (j=0;j<m;j++) PetscCall(PetscArraycpy(Gnew+n+(j+n)*ldg,Harray+j*ldh,m));
87 63 : Gnew[n+m+(n+m-1)*ldg] = beta;
88 63 : if (G) {
89 132 : for (j=0;j<n;j++) PetscCall(PetscArraycpy(Gnew+j*ldg,G+j*(n+1),n+1));
90 6 : PetscCall(PetscFree(G));
91 : }
92 63 : G = Gnew;
93 63 : n += m;
94 : } else {
95 : /* update Z = Z + V(:,1:m)*Q with Q=U(blk,:)*P(1:nrk,:)' */
96 63 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,m,*col+rank,NULL,&Q));
97 63 : PetscCall(MatDenseGetArray(Q,&Qarray));
98 63 : PetscCall(PetscBLASIntCast(m,&m_));
99 63 : PetscCall(PetscBLASIntCast(n,&n_));
100 63 : PetscCall(PetscBLASIntCast(rank,&rk_));
101 63 : PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&m_,&rk_,&rk_,&sone,U+its*m,&n_,L,&n_,&zero,Qarray+(*col)*m,&m_));
102 63 : PetscCall(MatDenseRestoreArray(Q,&Qarray));
103 63 : PetscCall(BVSetActiveColumns(*X1,*col,*col+rank));
104 63 : PetscCall(BVMult(*X1,1.0,1.0,lme->V,Q));
105 63 : PetscCall(MatDestroy(&Q));
106 : }
107 :
108 126 : if (pass==0) {
109 : /* solve compressed Lyapunov equation */
110 63 : PetscCall(PetscCalloc1(n,&r));
111 63 : PetscCall(PetscCalloc1(n*n,&L));
112 63 : r[0] = bnorm;
113 63 : errest = PetscAbsScalar(G[n+(n-1)*ldg]);
114 63 : PetscCall(LMEDenseHessLyapunovChol(lme,n,G,ldg,1,r,n,L,n,&errest));
115 63 : PetscCall(LMEMonitor(lme,*totalits+its,errest));
116 63 : PetscCall(PetscFree(r));
117 :
118 : /* check convergence */
119 63 : if (errest<lme->tol) {
120 57 : lme->errest += errest;
121 57 : PetscCall(PetscMalloc1(n*n,&U));
122 : /* transpose L */
123 1839 : for (j=0;j<n;j++) {
124 29133 : for (i=j+1;i<n;i++) {
125 27351 : L[i+j*n] = PetscConj(L[j+i*n]);
126 27351 : L[j+i*n] = 0.0;
127 : }
128 : }
129 57 : PetscCall(LMEDenseRankSVD(lme,n,L,n,U,n,&lrank));
130 57 : PetscCall(PetscInfo(lme,"Rank of the Cholesky factor = %" PetscInt_FMT "\n",lrank));
131 57 : nouter = its;
132 57 : its = -1;
133 57 : if (!fixed) { /* X1 was not set by user, allocate it with rank columns */
134 8 : rank = lrank;
135 8 : if (*col) PetscCall(BVResize(*X1,*col+rank,PETSC_TRUE));
136 4 : else PetscCall(BVDuplicateResize(C1,rank,X1));
137 49 : } else rank = PetscMin(lrank,rrank);
138 57 : PetscCall(PetscFree(G));
139 : break;
140 : } else {
141 6 : PetscCall(PetscFree(L));
142 6 : if (*totalits+its>=lme->max_it) *fail = PETSC_TRUE;
143 : }
144 : }
145 :
146 : /* restart with vector v_{m+1} */
147 252 : if (!*fail) PetscCall(BVCopyColumn(lme->V,m,0));
148 : }
149 : }
150 :
151 57 : *col += rank;
152 57 : *totalits += its+1;
153 57 : PetscCall(MatDenseRestoreArray(H,&Harray));
154 57 : PetscCall(MatDestroy(&H));
155 57 : if (L) PetscCall(PetscFree(L));
156 57 : if (U) PetscCall(PetscFree(U));
157 57 : PetscFunctionReturn(PETSC_SUCCESS);
158 : }
159 :
160 33 : static PetscErrorCode LMESolve_Krylov_Lyapunov(LME lme)
161 : {
162 33 : PetscBool fail,fixed = lme->X? PETSC_TRUE: PETSC_FALSE;
163 33 : PetscInt i,k,rank=0,col=0;
164 33 : Vec b;
165 33 : BV X1=NULL,C1;
166 33 : Mat X1m,X1t,C1m;
167 :
168 33 : PetscFunctionBegin;
169 33 : PetscCall(MatLRCGetMats(lme->C,NULL,&C1m,NULL,NULL));
170 33 : PetscCall(BVCreateFromMat(C1m,&C1));
171 33 : PetscCall(BVSetFromOptions(C1));
172 33 : PetscCall(BVGetActiveColumns(C1,NULL,&k));
173 33 : if (fixed) {
174 29 : PetscCall(MatLRCGetMats(lme->X,NULL,&X1m,NULL,NULL));
175 29 : PetscCall(BVCreateFromMat(X1m,&X1));
176 29 : PetscCall(BVSetFromOptions(X1));
177 29 : PetscCall(BVGetActiveColumns(X1,NULL,&rank));
178 29 : rank = rank/k;
179 : }
180 90 : for (i=0;i<k;i++) {
181 57 : PetscCall(BVGetColumn(C1,i,&b));
182 57 : PetscCall(LMESolve_Krylov_Lyapunov_Vec(lme,b,fixed,rank,C1,&X1,&col,&fail,&lme->its));
183 57 : PetscCall(BVRestoreColumn(C1,i,&b));
184 57 : if (fail) {
185 0 : lme->reason = LME_DIVERGED_ITS;
186 0 : break;
187 : }
188 : }
189 33 : if (lme->reason==LME_CONVERGED_ITERATING) lme->reason = LME_CONVERGED_TOL;
190 33 : PetscCall(BVCreateMat(X1,&X1t));
191 33 : if (fixed) PetscCall(MatCopy(X1t,X1m,SAME_NONZERO_PATTERN));
192 4 : else PetscCall(MatCreateLRC(NULL,X1t,NULL,NULL,&lme->X));
193 33 : PetscCall(MatDestroy(&X1t));
194 33 : PetscCall(BVDestroy(&C1));
195 33 : PetscCall(BVDestroy(&X1));
196 33 : PetscFunctionReturn(PETSC_SUCCESS);
197 : }
198 :
199 8 : SLEPC_EXTERN PetscErrorCode LMECreate_Krylov(LME lme)
200 : {
201 8 : PetscFunctionBegin;
202 8 : lme->ops->solve[LME_LYAPUNOV] = LMESolve_Krylov_Lyapunov;
203 8 : lme->ops->setup = LMESetUp_Krylov;
204 8 : PetscFunctionReturn(PETSC_SUCCESS);
205 : }
|