Actual source code: bvkrylov.c

slepc-3.13.3 2020-06-14
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 - the upper Hessenberg matrix
 25: .  ldh - leading dimension of H
 26: .  k - number of locked columns
 27: -  m - dimension of the Arnoldi basis

 29:    Output Parameters:
 30: +  m - the modified dimension
 31: .  beta - (optional) norm of last vector before normalization
 32: -  breakdown - (optional) flag indicating that breakdown occurred

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

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

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

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

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

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

 54:    Level: advanced

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

 73:   BVCheckSizes(V,1);

 77:   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);
 78:   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);
 79:   if (*m<=k) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");

 81:   BVSetActiveColumns(V,0,*m);
 82:   for (j=k;j<*m;j++) {
 83:     BVMatMultColumn(V,A,j);
 84:     if (PetscUnlikely(j==V->N-1)) {   /* safeguard in case the full basis is requested */
 85:       BV_OrthogonalizeColumn_Safe(V,j+1,NULL,beta,&lindep);
 86:     } else {
 87:       BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta,&lindep);
 88:     }
 89:     if (lindep) {
 90:       *m = j+1;
 91:       break;
 92:     }
 93:   }
 94:   if (breakdown) *breakdown = lindep;
 95:   /* extract Hessenberg matrix from the BV object */
 96:   BVGetBufferVec(V,&buf);
 97:   VecGetArray(buf,&a);
 98:   for (j=k;j<*m;j++) {
 99:     PetscArraycpy(H+j*ldh,a+V->nc+(j+1)*(V->nc+V->m),j+2);
100:   }
101:   VecRestoreArray(buf,&a);

103:   PetscObjectStateIncrease((PetscObject)V);
104:   return(0);
105: }

107: /*@C
108:    BVMatLanczos - Computes a Lanczos factorization associated with a matrix.

110:    Collective on V

112:    Input Parameters:
113: +  V - basis vectors context
114: .  A - the matrix
115: .  alpha - diagonal entries of tridiagonal matrix
116: .  beta - subdiagonal entries of tridiagonal matrix
117: .  k - number of locked columns
118: -  m - dimension of the Lanczos basis

120:    Output Parameters:
121: +  m - the modified dimension
122: -  breakdown - (optional) flag indicating that breakdown occurred

124:    Notes:
125:    Computes an m-step Lanczos factorization for matrix A, with full
126:    reorthogonalization. At each Lanczos step, the corresponding Lanczos
127:    vector is orthogonalized with respect to all previous Lanczos vectors.
128:    This is equivalent to computing an m-step Arnoldi factorization and
129:    exploting symmetry of the operator.

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

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

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

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

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

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

153:    Level: advanced

155: .seealso: BVMatArnoldi(), BVSetActiveColumns(), BVOrthonormalizeColumn()
156: @*/
157: PetscErrorCode BVMatLanczos(BV V,Mat A,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown)
158: {
160:   PetscScalar    *a;
161:   PetscInt       j;
162:   PetscBool      lindep=PETSC_FALSE;
163:   Vec            buf;

174:   BVCheckSizes(V,1);

178:   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);
179:   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);
180:   if (*m<=k) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");

182:   BVSetActiveColumns(V,0,*m);
183:   for (j=k;j<*m;j++) {
184:     BVMatMultColumn(V,A,j);
185:     if (PetscUnlikely(j==V->N-1)) {   /* safeguard in case the full basis is requested */
186:       BV_OrthogonalizeColumn_Safe(V,j+1,NULL,beta+j,&lindep);
187:     } else {
188:       BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta+j,&lindep);
189:     }
190:     if (lindep) {
191:       *m = j+1;
192:       break;
193:     }
194:   }
195:   if (breakdown) *breakdown = lindep;

197:   /* extract Hessenberg matrix from the BV buffer */
198:   BVGetBufferVec(V,&buf);
199:   VecGetArray(buf,&a);
200:   for (j=k;j<*m;j++) alpha[j] = PetscRealPart(a[V->nc+j+(j+1)*(V->nc+V->m)]);
201:   VecRestoreArray(buf,&a);

203:   PetscObjectStateIncrease((PetscObject)V);
204:   return(0);
205: }