Actual source code: mfnkrylov.c
slepc-3.7.0 2016-05-16
1: /*
3: SLEPc matrix function solver: "krylov"
5: Method: Arnoldi with Eiermann-Ernst restart
7: Algorithm:
9: Build Arnoldi approximations using f(H) for the Hessenberg matrix H,
10: restart by discarding the Krylov basis but keeping H.
12: References:
14: [1] M. Eiermann and O. Ernst, "A restarted Krylov subspace method
15: for the evaluation of matrix functions", SIAM J. Numer. Anal.
16: 44(6):2481-2504, 2006.
18: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: SLEPc - Scalable Library for Eigenvalue Problem Computations
20: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
22: This file is part of SLEPc.
24: SLEPc is free software: you can redistribute it and/or modify it under the
25: terms of version 3 of the GNU Lesser General Public License as published by
26: the Free Software Foundation.
28: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
29: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
30: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
31: more details.
33: You should have received a copy of the GNU Lesser General Public License
34: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: */
38: #include <slepc/private/mfnimpl.h>
39: #include <slepcblaslapack.h>
43: PetscErrorCode MFNSetUp_Krylov(MFN mfn)
44: {
46: PetscInt N;
49: MatGetSize(mfn->A,&N,NULL);
50: if (!mfn->ncv) mfn->ncv = PetscMin(30,N);
51: if (!mfn->max_it) mfn->max_it = 100;
52: MFNAllocateSolution(mfn,1);
53: return(0);
54: }
58: PetscErrorCode MFNBasicArnoldi(MFN mfn,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscReal *beta,PetscBool *breakdown)
59: {
61: PetscInt j,m = *M;
62: Vec vj,vj1;
65: BVSetActiveColumns(mfn->V,0,m);
66: for (j=k;j<m;j++) {
67: BVGetColumn(mfn->V,j,&vj);
68: BVGetColumn(mfn->V,j+1,&vj1);
69: MatMult(mfn->A,vj,vj1);
70: BVRestoreColumn(mfn->V,j,&vj);
71: BVRestoreColumn(mfn->V,j+1,&vj1);
72: BVOrthogonalizeColumn(mfn->V,j+1,H+ldh*j,beta,breakdown);
73: H[j+1+ldh*j] = *beta;
74: if (*breakdown) {
75: *M = j+1;
76: break;
77: } else {
78: BVScaleColumn(mfn->V,j+1,1.0/(*beta));
79: }
80: }
81: return(0);
82: }
86: PetscErrorCode MFNSolve_Krylov(MFN mfn,Vec b,Vec x)
87: {
89: PetscInt n=0,m,ld,ldh,j;
90: PetscBLASInt m_,inc=1;
91: Mat G=NULL,H=NULL;
92: Vec F=NULL;
93: PetscScalar *array,*farray,*garray,*harray;
94: PetscReal beta,nrm=1.0;
95: PetscBool breakdown,set,flg,symm=PETSC_FALSE;
98: m = mfn->ncv;
99: ld = m+1;
100: PetscCalloc1(ld*ld,&array);
101: VecSet(x,0.0);
103: /* set initial vector to b/||b|| */
104: BVInsertVec(mfn->V,0,b);
105: BVScaleColumn(mfn->V,0,1.0/mfn->bnorm);
107: /* Restart loop */
108: while (mfn->reason == MFN_CONVERGED_ITERATING) {
109: mfn->its++;
111: /* compute Arnoldi factorization */
112: MFNBasicArnoldi(mfn,array,ld,0,&m,&beta,&breakdown);
114: /* save previous Hessenberg matrix in G; allocate new storage for H and f(H) */
115: if (mfn->its>1) { G = H; H = NULL; }
116: ldh = n+m;
117: MFN_CreateVec(ldh,&F);
118: MFN_CreateDenseMat(ldh,&H);
120: /* glue together the previous H and the new H obtained with Arnoldi */
121: MatDenseGetArray(H,&harray);
122: for (j=0;j<m;j++) {
123: PetscMemcpy(harray+n+(j+n)*ldh,array+j*ld,m*sizeof(PetscScalar));
124: }
125: if (mfn->its>1) {
126: MatDenseGetArray(G,&garray);
127: for (j=0;j<n;j++) {
128: PetscMemcpy(harray+j*ldh,garray+j*n,n*sizeof(PetscScalar));
129: }
130: MatDenseRestoreArray(G,&garray);
131: MatDestroy(&G);
132: harray[n+(n-1)*ldh] = beta;
133: }
134: MatDenseRestoreArray(H,&harray);
136: if (mfn->its==1) {
137: /* set symmetry flag of H from A */
138: MatIsHermitianKnown(mfn->A,&set,&flg);
139: symm = set? flg: PETSC_FALSE;
140: if (symm) {
141: MatSetOption(H,MAT_HERMITIAN,PETSC_TRUE);
142: }
143: }
145: /* evaluate f(H) */
146: FNEvaluateFunctionMatVec(mfn->fn,H,F);
148: /* x += ||b||*V*f(H)*e_1 */
149: VecGetArray(F,&farray);
150: if (mfn->its>1) {
151: PetscBLASIntCast(m,&m_);
152: nrm = BLASnrm2_(&m_,farray+n,&inc); /* relative norm of the update ||u||/||b|| */
153: MFNMonitor(mfn,mfn->its,nrm);
154: } else {
155: MFNMonitor(mfn,1,1.0); /* no error estimate available */
156: }
157: for (j=0;j<m;j++) farray[j+n] *= mfn->bnorm;
158: BVSetActiveColumns(mfn->V,0,m);
159: BVMultVec(mfn->V,1.0,1.0,x,farray+n);
160: VecRestoreArray(F,&farray);
162: /* check convergence */
163: if (mfn->its>1) {
164: if (mfn->its >= mfn->max_it) mfn->reason = MFN_DIVERGED_ITS;
165: if (m<mfn->ncv || breakdown || beta==0.0 || nrm<mfn->tol) mfn->reason = MFN_CONVERGED_TOL;
166: }
168: /* restart with vector v_{m+1} */
169: if (mfn->reason == MFN_CONVERGED_ITERATING) {
170: BVCopyColumn(mfn->V,m,0);
171: n += m;
172: }
173: }
175: MatDestroy(&H);
176: MatDestroy(&G);
177: VecDestroy(&F);
178: PetscFree(array);
179: return(0);
180: }
184: PETSC_EXTERN PetscErrorCode MFNCreate_Krylov(MFN mfn)
185: {
187: mfn->ops->solve = MFNSolve_Krylov;
188: mfn->ops->setup = MFNSetUp_Krylov;
189: return(0);
190: }