Actual source code: mfnkrylov.c
slepc-3.14.0 2020-09-30
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: SLEPc matrix function solver: "krylov"
13: Method: Arnoldi with Eiermann-Ernst restart
15: Algorithm:
17: Build Arnoldi approximations using f(H) for the Hessenberg matrix H,
18: restart by discarding the Krylov basis but keeping H.
20: References:
22: [1] M. Eiermann and O. Ernst, "A restarted Krylov subspace method
23: for the evaluation of matrix functions", SIAM J. Numer. Anal.
24: 44(6):2481-2504, 2006.
25: */
27: #include <slepc/private/mfnimpl.h>
28: #include <slepcblaslapack.h>
30: PetscErrorCode MFNSetUp_Krylov(MFN mfn)
31: {
33: PetscInt N;
36: MatGetSize(mfn->A,&N,NULL);
37: if (mfn->ncv==PETSC_DEFAULT) mfn->ncv = PetscMin(30,N);
38: if (mfn->max_it==PETSC_DEFAULT) mfn->max_it = 100;
39: MFNAllocateSolution(mfn,1);
40: return(0);
41: }
43: PetscErrorCode MFNSolve_Krylov(MFN mfn,Vec b,Vec x)
44: {
46: PetscInt n=0,m,ld,ldh,j;
47: PetscBLASInt m_,inc=1;
48: Mat G=NULL,H=NULL;
49: Vec F=NULL;
50: PetscScalar *array,*farray,*garray,*harray;
51: PetscReal beta,betaold=0.0,nrm=1.0;
52: PetscBool breakdown,set,flg,symm=PETSC_FALSE;
55: m = mfn->ncv;
56: ld = m+1;
57: PetscCalloc1(ld*ld,&array);
59: /* set initial vector to b/||b|| */
60: BVInsertVec(mfn->V,0,b);
61: BVScaleColumn(mfn->V,0,1.0/mfn->bnorm);
62: VecSet(x,0.0);
64: /* Restart loop */
65: while (mfn->reason == MFN_CONVERGED_ITERATING) {
66: mfn->its++;
68: /* compute Arnoldi factorization */
69: BVMatArnoldi(mfn->V,mfn->transpose_solve?mfn->AT:mfn->A,array,ld,0,&m,&beta,&breakdown);
71: /* save previous Hessenberg matrix in G; allocate new storage for H and f(H) */
72: if (mfn->its>1) { G = H; H = NULL; }
73: ldh = n+m;
74: MFN_CreateVec(ldh,&F);
75: MFN_CreateDenseMat(ldh,&H);
77: /* glue together the previous H and the new H obtained with Arnoldi */
78: MatDenseGetArray(H,&harray);
79: for (j=0;j<m;j++) {
80: PetscArraycpy(harray+n+(j+n)*ldh,array+j*ld,m);
81: }
82: if (mfn->its>1) {
83: MatDenseGetArray(G,&garray);
84: for (j=0;j<n;j++) {
85: PetscArraycpy(harray+j*ldh,garray+j*n,n);
86: }
87: MatDenseRestoreArray(G,&garray);
88: MatDestroy(&G);
89: harray[n+(n-1)*ldh] = betaold;
90: }
91: MatDenseRestoreArray(H,&harray);
93: if (mfn->its==1) {
94: /* set symmetry flag of H from A */
95: MatIsHermitianKnown(mfn->A,&set,&flg);
96: symm = set? flg: PETSC_FALSE;
97: if (symm) {
98: MatSetOption(H,MAT_HERMITIAN,PETSC_TRUE);
99: }
100: }
102: /* evaluate f(H) */
103: FNEvaluateFunctionMatVec(mfn->fn,H,F);
105: /* x += ||b||*V*f(H)*e_1 */
106: VecGetArray(F,&farray);
107: PetscBLASIntCast(m,&m_);
108: nrm = BLASnrm2_(&m_,farray+n,&inc); /* relative norm of the update ||u||/||b|| */
109: MFNMonitor(mfn,mfn->its,nrm);
110: for (j=0;j<m;j++) farray[j+n] *= mfn->bnorm;
111: BVSetActiveColumns(mfn->V,0,m);
112: BVMultVec(mfn->V,1.0,1.0,x,farray+n);
113: VecRestoreArray(F,&farray);
115: /* check convergence */
116: if (mfn->its >= mfn->max_it) mfn->reason = MFN_DIVERGED_ITS;
117: if (mfn->its>1) {
118: if (m<mfn->ncv || breakdown || beta==0.0 || nrm<mfn->tol) mfn->reason = MFN_CONVERGED_TOL;
119: }
121: /* restart with vector v_{m+1} */
122: if (mfn->reason == MFN_CONVERGED_ITERATING) {
123: BVCopyColumn(mfn->V,m,0);
124: n += m;
125: betaold = beta;
126: }
127: }
129: MatDestroy(&H);
130: MatDestroy(&G);
131: VecDestroy(&F);
132: PetscFree(array);
133: return(0);
134: }
136: SLEPC_EXTERN PetscErrorCode MFNCreate_Krylov(MFN mfn)
137: {
139: mfn->ops->solve = MFNSolve_Krylov;
140: mfn->ops->setup = MFNSetUp_Krylov;
141: return(0);
142: }