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> /*I "slepcbv.h" I*/
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: BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta,&lindep);
85: if (lindep) {
86: *m = j+1;
87: break;
88: }
89: }
90: if (breakdown) *breakdown = lindep;
91: /* extract Hessenberg matrix from the BV object */
92: BVGetBufferVec(V,&buf);
93: VecGetArray(buf,&a);
94: for (j=k;j<*m;j++) {
95: PetscArraycpy(H+j*ldh,a+V->nc+(j+1)*(V->nc+V->m),j+2);
96: }
97: VecRestoreArray(buf,&a);
99: PetscObjectStateIncrease((PetscObject)V);
100: return(0);
101: }
103: /*@C
104: BVMatLanczos - Computes a Lanczos factorization associated with a matrix.
106: Collective on V
108: Input Parameters:
109: + V - basis vectors context
110: . A - the matrix
111: . alpha - diagonal entries of tridiagonal matrix
112: . beta - subdiagonal entries of tridiagonal matrix
113: . k - number of locked columns
114: - m - dimension of the Lanczos basis
116: Output Parameters:
117: + m - the modified dimension
118: - breakdown - (optional) flag indicating that breakdown occurred
120: Notes:
121: Computes an m-step Lanczos factorization for matrix A, with full
122: reorthogonalization. At each Lanczos step, the corresponding Lanczos
123: vector is orthogonalized with respect to all previous Lanczos vectors.
124: This is equivalent to computing an m-step Arnoldi factorization and
125: exploting symmetry of the operator.
127: The first k columns are assumed to be locked and therefore they are
128: not modified. On exit, the following relation is satisfied:
130: A * V - V * T = beta_m*v_m * e_m^T
132: where the columns of V are the Lanczos vectors (which are B-orthonormal),
133: T is a real symmetric tridiagonal matrix, and e_m is the m-th vector of
134: the canonical basis. The tridiagonal is stored as two arrays: alpha
135: contains the diagonal elements, beta the off-diagonal. On exit, the last
136: element of beta contains the B-norm of V[m] before normalization.
137: The basis V must have at least m+1 columns, while the arrays alpha and
138: beta must have space for at least m elements.
140: The breakdown flag indicates that orthogonalization failed, see
141: BVOrthonormalizeColumn(). In that case, on exit m contains the index of
142: the column that failed.
144: The values of k and m are not restricted to the active columns of V.
146: To create a Lanczos factorization from scratch, set k=0 and make sure the
147: first column contains the normalized initial vector.
149: Level: advanced
151: .seealso: BVMatArnoldi(), BVSetActiveColumns(), BVOrthonormalizeColumn()
152: @*/
153: PetscErrorCode BVMatLanczos(BV V,Mat A,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown)154: {
156: PetscScalar *a;
157: PetscInt j;
158: PetscBool lindep=PETSC_FALSE;
159: Vec buf;
170: BVCheckSizes(V,1);
174: 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);
175: 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);
176: if (*m<=k) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");
178: BVSetActiveColumns(V,0,*m);
179: for (j=k;j<*m;j++) {
180: BVMatMultColumn(V,A,j);
181: BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta+j,&lindep);
182: if (lindep) {
183: *m = j+1;
184: break;
185: }
186: }
187: if (breakdown) *breakdown = lindep;
189: /* extract Hessenberg matrix from the BV buffer */
190: BVGetBufferVec(V,&buf);
191: VecGetArray(buf,&a);
192: for (j=k;j<*m;j++) alpha[j] = PetscRealPart(a[V->nc+j+(j+1)*(V->nc+V->m)]);
193: VecRestoreArray(buf,&a);
195: PetscObjectStateIncrease((PetscObject)V);
196: return(0);
197: }