Actual source code: mfnkrylov.c
slepc-3.16.1 2021-11-17
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: 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: {
45: PetscErrorCode ierr;
46: PetscInt n=0,m,ld,ldh,j;
47: PetscBLASInt m_,inc=1;
48: Mat M,G=NULL,H=NULL;
49: Vec F=NULL;
50: PetscScalar *marray,*farray,*harray;
51: const PetscScalar *garray;
52: PetscReal beta,betaold=0.0,nrm=1.0;
53: PetscBool breakdown;
56: m = mfn->ncv;
57: ld = m+1;
58: MatCreateSeqDense(PETSC_COMM_SELF,ld,m,NULL,&M);
59: MatDenseGetArray(M,&marray);
61: /* set initial vector to b/||b|| */
62: BVInsertVec(mfn->V,0,b);
63: BVScaleColumn(mfn->V,0,1.0/mfn->bnorm);
64: VecSet(x,0.0);
66: /* Restart loop */
67: while (mfn->reason == MFN_CONVERGED_ITERATING) {
68: mfn->its++;
70: /* compute Arnoldi factorization */
71: BVMatArnoldi(mfn->V,mfn->transpose_solve?mfn->AT:mfn->A,M,0,&m,&beta,&breakdown);
73: /* save previous Hessenberg matrix in G; allocate new storage for H and f(H) */
74: if (mfn->its>1) { G = H; H = NULL; }
75: ldh = n+m;
76: MFN_CreateVec(ldh,&F);
77: MFN_CreateDenseMat(ldh,&H);
79: /* glue together the previous H and the new H obtained with Arnoldi */
80: MatDenseGetArray(H,&harray);
81: for (j=0;j<m;j++) {
82: PetscArraycpy(harray+n+(j+n)*ldh,marray+j*ld,m);
83: }
84: if (mfn->its>1) {
85: MatDenseGetArrayRead(G,&garray);
86: for (j=0;j<n;j++) {
87: PetscArraycpy(harray+j*ldh,garray+j*n,n);
88: }
89: MatDenseRestoreArrayRead(G,&garray);
90: MatDestroy(&G);
91: harray[n+(n-1)*ldh] = betaold;
92: }
93: MatDenseRestoreArray(H,&harray);
95: if (mfn->its==1) {
96: /* set symmetry flag of H from A */
97: MatPropagateSymmetryOptions(mfn->A,H);
98: }
100: /* evaluate f(H) */
101: FNEvaluateFunctionMatVec(mfn->fn,H,F);
103: /* x += ||b||*V*f(H)*e_1 */
104: VecGetArray(F,&farray);
105: PetscBLASIntCast(m,&m_);
106: nrm = BLASnrm2_(&m_,farray+n,&inc); /* relative norm of the update ||u||/||b|| */
107: MFNMonitor(mfn,mfn->its,nrm);
108: for (j=0;j<m;j++) farray[j+n] *= mfn->bnorm;
109: BVSetActiveColumns(mfn->V,0,m);
110: BVMultVec(mfn->V,1.0,1.0,x,farray+n);
111: VecRestoreArray(F,&farray);
113: /* check convergence */
114: if (mfn->its >= mfn->max_it) mfn->reason = MFN_DIVERGED_ITS;
115: if (mfn->its>1) {
116: if (m<mfn->ncv || breakdown || beta==0.0 || nrm<mfn->tol) mfn->reason = MFN_CONVERGED_TOL;
117: }
119: /* restart with vector v_{m+1} */
120: if (mfn->reason == MFN_CONVERGED_ITERATING) {
121: BVCopyColumn(mfn->V,m,0);
122: n += m;
123: betaold = beta;
124: }
125: }
127: MatDestroy(&H);
128: MatDestroy(&G);
129: VecDestroy(&F);
130: MatDenseRestoreArray(M,&marray);
131: MatDestroy(&M);
132: return(0);
133: }
135: SLEPC_EXTERN PetscErrorCode MFNCreate_Krylov(MFN mfn)
136: {
138: mfn->ops->solve = MFNSolve_Krylov;
139: mfn->ops->setup = MFNSetUp_Krylov;
140: return(0);
141: }