Actual source code: bvkrylov.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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:    BV routines related to Krylov decompositions
 12: */

 14: #include <slepc/private/bvimpl.h>

 16: /*@
 17:    BVMatArnoldi - Computes an Arnoldi factorization associated with a matrix.

 19:    Collective on V

 21:    Input Parameters:
 22: +  V - basis vectors context
 23: .  A - the matrix
 24: .  H - (optional) the upper Hessenberg matrix
 25: .  k - number of locked columns
 26: -  m - dimension of the Arnoldi basis, may be modified

 28:    Output Parameters:
 29: +  beta - (optional) norm of last vector before normalization
 30: -  breakdown - (optional) flag indicating that breakdown occurred

 32:    Notes:
 33:    Computes an m-step Arnoldi factorization for matrix A. The first k columns
 34:    are assumed to be locked and therefore they are not modified. On exit, the
 35:    following relation is satisfied:

 37:                     A * V - V * H = beta*v_m * e_m^T

 39:    where the columns of V are the Arnoldi vectors (which are orthonormal), H is
 40:    an upper Hessenberg matrix, e_m is the m-th vector of the canonical basis.
 41:    On exit, beta contains the norm of V[m] before normalization.

 43:    The breakdown flag indicates that orthogonalization failed, see
 44:    BVOrthonormalizeColumn(). In that case, on exit m contains the index of
 45:    the column that failed.

 47:    The values of k and m are not restricted to the active columns of V.

 49:    To create an Arnoldi factorization from scratch, set k=0 and make sure the
 50:    first column contains the normalized initial vector.

 52:    Level: advanced

 54: .seealso: BVMatLanczos(), BVSetActiveColumns(), BVOrthonormalizeColumn()
 55: @*/
 56: PetscErrorCode BVMatArnoldi(BV V,Mat A,Mat H,PetscInt k,PetscInt *m,PetscReal *beta,PetscBool *breakdown)
 57: {
 58:   PetscErrorCode    ierr;
 59:   PetscScalar       *h;
 60:   const PetscScalar *a;
 61:   PetscInt          j,ldh,rows,cols;
 62:   PetscBool         lindep=PETSC_FALSE;
 63:   Vec               buf;

 72:   BVCheckSizes(V,1);

 76:   if (k<0 || k>V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument k has wrong value %D, should be between 0 and %D",k,V->m);
 77:   if (*m<1 || *m>V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m has wrong value %D, should be between 1 and %D",*m,V->m);
 78:   if (*m<=k) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");
 79:   if (H) {
 83:     MatGetSize(H,&rows,&cols);
 84:     MatDenseGetLDA(H,&ldh);
 85:     if (rows<*m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Matrix H has %D rows, should have at least %D",rows,*m);
 86:     if (cols<*m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Matrix H has %D columns, should have at least %D",cols,*m);
 87:   }

 89:   for (j=k;j<*m;j++) {
 90:     BVMatMultColumn(V,A,j);
 91:     if (PetscUnlikely(j==V->N-1)) {   /* safeguard in case the full basis is requested */
 92:       BV_OrthogonalizeColumn_Safe(V,j+1,NULL,beta,&lindep);
 93:     } else {
 94:       BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta,&lindep);
 95:     }
 96:     if (lindep) {
 97:       *m = j+1;
 98:       break;
 99:     }
100:   }
101:   if (breakdown) *breakdown = lindep;
102:   if (lindep) { PetscInfo1(V,"Arnoldi finished early at m=%D\n",*m); }

104:   if (H) {
105:     MatDenseGetArray(H,&h);
106:     BVGetBufferVec(V,&buf);
107:     VecGetArrayRead(buf,&a);
108:     for (j=k;j<*m-1;j++) {
109:       PetscArraycpy(h+j*ldh,a+V->nc+(j+1)*(V->nc+V->m),j+2);
110:     }
111:     PetscArraycpy(h+(*m-1)*ldh,a+V->nc+(*m)*(V->nc+V->m),*m);
112:     if (ldh>*m) h[(*m)+(*m-1)*ldh] = a[V->nc+(*m)+(*m)*(V->nc+V->m)];
113:     VecRestoreArrayRead(buf,&a);
114:     MatDenseRestoreArray(H,&h);
115:   }

117:   PetscObjectStateIncrease((PetscObject)V);
118:   return(0);
119: }

121: /*@C
122:    BVMatLanczos - Computes a Lanczos factorization associated with a matrix.

124:    Collective on V

126:    Input Parameters:
127: +  V - basis vectors context
128: .  A - the matrix
129: .  alpha - diagonal entries of tridiagonal matrix
130: .  beta - subdiagonal entries of tridiagonal matrix
131: -  k - number of locked columns

133:    Input/Output Parameter:
134: .  m - dimension of the Lanczos basis, may be modified

136:    Output Parameter:
137: .  breakdown - (optional) flag indicating that breakdown occurred

139:    Notes:
140:    Computes an m-step Lanczos factorization for matrix A, with full
141:    reorthogonalization. At each Lanczos step, the corresponding Lanczos
142:    vector is orthogonalized with respect to all previous Lanczos vectors.
143:    This is equivalent to computing an m-step Arnoldi factorization and
144:    exploting symmetry of the operator.

146:    The first k columns are assumed to be locked and therefore they are
147:    not modified. On exit, the following relation is satisfied:

149:                     A * V - V * T = beta_m*v_m * e_m^T

151:    where the columns of V are the Lanczos vectors (which are B-orthonormal),
152:    T is a real symmetric tridiagonal matrix, and e_m is the m-th vector of
153:    the canonical basis. The tridiagonal is stored as two arrays: alpha
154:    contains the diagonal elements, beta the off-diagonal. On exit, the last
155:    element of beta contains the B-norm of V[m] before normalization.
156:    The basis V must have at least m+1 columns, while the arrays alpha and
157:    beta must have space for at least m elements.

159:    The breakdown flag indicates that orthogonalization failed, see
160:    BVOrthonormalizeColumn(). In that case, on exit m contains the index of
161:    the column that failed.

163:    The values of k and m are not restricted to the active columns of V.

165:    To create a Lanczos factorization from scratch, set k=0 and make sure the
166:    first column contains the normalized initial vector.

168:    Level: advanced

170: .seealso: BVMatArnoldi(), BVSetActiveColumns(), BVOrthonormalizeColumn()
171: @*/
172: PetscErrorCode BVMatLanczos(BV V,Mat A,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown)
173: {
175:   PetscScalar    *a;
176:   PetscInt       j;
177:   PetscBool      lindep=PETSC_FALSE;
178:   Vec            buf;

189:   BVCheckSizes(V,1);

193:   if (k<0 || k>V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument k has wrong value %D, should be between 0 and %D",k,V->m);
194:   if (*m<1 || *m>V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m has wrong value %D, should be between 1 and %D",*m,V->m);
195:   if (*m<=k) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");

197:   for (j=k;j<*m;j++) {
198:     BVMatMultColumn(V,A,j);
199:     if (PetscUnlikely(j==V->N-1)) {   /* safeguard in case the full basis is requested */
200:       BV_OrthogonalizeColumn_Safe(V,j+1,NULL,beta+j,&lindep);
201:     } else {
202:       BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta+j,&lindep);
203:     }
204:     if (lindep) {
205:       *m = j+1;
206:       break;
207:     }
208:   }
209:   if (breakdown) *breakdown = lindep;
210:   if (lindep) { PetscInfo1(V,"Lanczos finished early at m=%D\n",*m); }

212:   /* extract Hessenberg matrix from the BV buffer */
213:   BVGetBufferVec(V,&buf);
214:   VecGetArray(buf,&a);
215:   for (j=k;j<*m;j++) alpha[j] = PetscRealPart(a[V->nc+j+(j+1)*(V->nc+V->m)]);
216:   VecRestoreArray(buf,&a);

218:   PetscObjectStateIncrease((PetscObject)V);
219:   return(0);
220: }