Actual source code: mfnkrylov.c
slepc-3.17.0 2022-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: {
32: PetscInt N;
34: MatGetSize(mfn->A,&N,NULL);
35: if (mfn->ncv==PETSC_DEFAULT) mfn->ncv = PetscMin(30,N);
36: if (mfn->max_it==PETSC_DEFAULT) mfn->max_it = 100;
37: MFNAllocateSolution(mfn,1);
38: PetscFunctionReturn(0);
39: }
41: PetscErrorCode MFNSolve_Krylov(MFN mfn,Vec b,Vec x)
42: {
43: PetscInt n=0,m,ld,ldh,j;
44: PetscBLASInt m_,inc=1;
45: Mat M,G=NULL,H=NULL;
46: Vec F=NULL;
47: PetscScalar *marray,*farray,*harray;
48: const PetscScalar *garray;
49: PetscReal beta,betaold=0.0,nrm=1.0;
50: PetscBool breakdown;
52: m = mfn->ncv;
53: ld = m+1;
54: MatCreateSeqDense(PETSC_COMM_SELF,ld,m,NULL,&M);
55: MatDenseGetArray(M,&marray);
57: /* set initial vector to b/||b|| */
58: BVInsertVec(mfn->V,0,b);
59: BVScaleColumn(mfn->V,0,1.0/mfn->bnorm);
60: VecSet(x,0.0);
62: /* Restart loop */
63: while (mfn->reason == MFN_CONVERGED_ITERATING) {
64: mfn->its++;
66: /* compute Arnoldi factorization */
67: BVMatArnoldi(mfn->V,mfn->transpose_solve?mfn->AT:mfn->A,M,0,&m,&beta,&breakdown);
69: /* save previous Hessenberg matrix in G; allocate new storage for H and f(H) */
70: if (mfn->its>1) { G = H; H = NULL; }
71: ldh = n+m;
72: MFN_CreateVec(ldh,&F);
73: MFN_CreateDenseMat(ldh,&H);
75: /* glue together the previous H and the new H obtained with Arnoldi */
76: MatDenseGetArray(H,&harray);
77: for (j=0;j<m;j++) PetscArraycpy(harray+n+(j+n)*ldh,marray+j*ld,m);
78: if (mfn->its>1) {
79: MatDenseGetArrayRead(G,&garray);
80: for (j=0;j<n;j++) PetscArraycpy(harray+j*ldh,garray+j*n,n);
81: MatDenseRestoreArrayRead(G,&garray);
82: MatDestroy(&G);
83: harray[n+(n-1)*ldh] = betaold;
84: }
85: MatDenseRestoreArray(H,&harray);
87: if (mfn->its==1) {
88: /* set symmetry flag of H from A */
89: MatPropagateSymmetryOptions(mfn->A,H);
90: }
92: /* evaluate f(H) */
93: FNEvaluateFunctionMatVec(mfn->fn,H,F);
95: /* x += ||b||*V*f(H)*e_1 */
96: VecGetArray(F,&farray);
97: PetscBLASIntCast(m,&m_);
98: nrm = BLASnrm2_(&m_,farray+n,&inc); /* relative norm of the update ||u||/||b|| */
99: MFNMonitor(mfn,mfn->its,nrm);
100: for (j=0;j<m;j++) farray[j+n] *= mfn->bnorm;
101: BVSetActiveColumns(mfn->V,0,m);
102: BVMultVec(mfn->V,1.0,1.0,x,farray+n);
103: VecRestoreArray(F,&farray);
105: /* check convergence */
106: if (mfn->its >= mfn->max_it) mfn->reason = MFN_DIVERGED_ITS;
107: if (mfn->its>1) {
108: if (m<mfn->ncv || breakdown || beta==0.0 || nrm<mfn->tol) mfn->reason = MFN_CONVERGED_TOL;
109: }
111: /* restart with vector v_{m+1} */
112: if (mfn->reason == MFN_CONVERGED_ITERATING) {
113: BVCopyColumn(mfn->V,m,0);
114: n += m;
115: betaold = beta;
116: }
117: }
119: MatDestroy(&H);
120: MatDestroy(&G);
121: VecDestroy(&F);
122: MatDenseRestoreArray(M,&marray);
123: MatDestroy(&M);
124: PetscFunctionReturn(0);
125: }
127: SLEPC_EXTERN PetscErrorCode MFNCreate_Krylov(MFN mfn)
128: {
129: mfn->ops->solve = MFNSolve_Krylov;
130: mfn->ops->setup = MFNSetUp_Krylov;
131: PetscFunctionReturn(0);
132: }