PetscErrorCode EPSSetUp_Lanczos(EPS eps) { EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data; BVOrthogRefineType refine; PetscReal eta; PetscErrorCode ierr; PetscFunctionBegin; ierr = EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);CHKERRQ(ierr); if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd"); if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv); if (!eps->which) { ierr = EPSSetWhichEigenpairs_Default(eps);CHKERRQ(ierr); } switch (eps->which) { case EPS_LARGEST_IMAGINARY: case EPS_SMALLEST_IMAGINARY: case EPS_TARGET_IMAGINARY: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which"); default: ; /* default case to remove warning */ } if (!eps->extraction) { ierr = EPSSetExtraction(eps,EPS_RITZ);CHKERRQ(ierr); } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type"); if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver"); ierr = EPSAllocateSolution(eps,1);CHKERRQ(ierr); ierr = EPS_SetInnerProduct(eps);CHKERRQ(ierr); if (lanczos->reorthog != EPS_LANCZOS_REORTHOG_FULL) { ierr = BVGetOrthogonalization(eps->V,NULL,&refine,&eta);CHKERRQ(ierr); ierr = BVSetOrthogonalization(eps->V,BV_ORTHOG_MGS,refine,eta);CHKERRQ(ierr); ierr = PetscInfo(eps,"Switching to MGS orthogonalization\n");CHKERRQ(ierr); } if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) { ierr = BVDuplicate(eps->V,&lanczos->AV);CHKERRQ(ierr); } ierr = DSSetType(eps->ds,DSHEP);CHKERRQ(ierr); ierr = DSSetCompact(eps->ds,PETSC_TRUE);CHKERRQ(ierr); ierr = DSAllocate(eps->ds,eps->ncv+1);CHKERRQ(ierr); if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) { ierr = EPSSetWorkVecs(eps,1);CHKERRQ(ierr); } /* dispatch solve method */ if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems"); if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems"); eps->ops->solve = EPSSolve_Lanczos; PetscFunctionReturn(0); }
PetscErrorCode SVDSolve_TRLanczos(SVD svd) { PetscErrorCode ierr; SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data; PetscReal *alpha,*beta,lastbeta,norm; PetscScalar *Q,*swork=NULL,*w; PetscInt i,k,l,nv,ld; Mat U,VT; PetscBool conv; BVOrthogType orthog; PetscFunctionBegin; /* allocate working space */ ierr = DSGetLeadingDimension(svd->ds,&ld);CHKERRQ(ierr); ierr = BVGetOrthogonalization(svd->V,&orthog,NULL,NULL);CHKERRQ(ierr); ierr = PetscMalloc1(ld,&w);CHKERRQ(ierr); if (lanczos->oneside && orthog == BV_ORTHOG_CGS) { ierr = PetscMalloc1(svd->ncv+1,&swork);CHKERRQ(ierr); } /* normalize start vector */ if (!svd->nini) { ierr = BVSetRandomColumn(svd->V,0,svd->rand);CHKERRQ(ierr); ierr = BVNormColumn(svd->V,0,NORM_2,&norm);CHKERRQ(ierr); ierr = BVScaleColumn(svd->V,0,1.0/norm);CHKERRQ(ierr); } l = 0; while (svd->reason == SVD_CONVERGED_ITERATING) { svd->its++; /* inner loop */ nv = PetscMin(svd->nconv+svd->mpd,svd->ncv); ierr = BVSetActiveColumns(svd->V,svd->nconv,nv);CHKERRQ(ierr); ierr = BVSetActiveColumns(svd->U,svd->nconv,nv);CHKERRQ(ierr); ierr = DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);CHKERRQ(ierr); beta = alpha + ld; if (lanczos->oneside) { if (orthog == BV_ORTHOG_MGS) { ierr = SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv);CHKERRQ(ierr); } else { ierr = SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);CHKERRQ(ierr); } } else { ierr = SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv+l,nv);CHKERRQ(ierr); } lastbeta = beta[nv-1]; ierr = DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);CHKERRQ(ierr); ierr = BVScaleColumn(svd->V,nv,1.0/lastbeta);CHKERRQ(ierr); /* compute SVD of general matrix */ ierr = DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);CHKERRQ(ierr); if (l==0) { ierr = DSSetState(svd->ds,DS_STATE_INTERMEDIATE);CHKERRQ(ierr); } else { ierr = DSSetState(svd->ds,DS_STATE_RAW);CHKERRQ(ierr); } ierr = DSSolve(svd->ds,w,NULL);CHKERRQ(ierr); ierr = DSSort(svd->ds,w,NULL,NULL,NULL,NULL);CHKERRQ(ierr); /* compute error estimates */ k = 0; conv = PETSC_TRUE; ierr = DSGetArray(svd->ds,DS_MAT_U,&Q);CHKERRQ(ierr); ierr = DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);CHKERRQ(ierr); beta = alpha + ld; for (i=svd->nconv;i<nv;i++) { svd->sigma[i] = PetscRealPart(w[i]); beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta; svd->errest[i] = PetscAbsScalar(beta[i]); if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i]; if (conv) { if (svd->errest[i] < svd->tol) k++; else conv = PETSC_FALSE; } } ierr = DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);CHKERRQ(ierr); ierr = DSRestoreArray(svd->ds,DS_MAT_U,&Q);CHKERRQ(ierr); /* check convergence and update l */ if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS; if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL; if (svd->reason != SVD_CONVERGED_ITERATING) l = 0; else l = PetscMax((nv-svd->nconv-k)/2,0); /* compute converged singular vectors and restart vectors */ ierr = DSGetMat(svd->ds,DS_MAT_VT,&VT);CHKERRQ(ierr); ierr = BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k+l);CHKERRQ(ierr); ierr = MatDestroy(&VT);CHKERRQ(ierr); ierr = DSGetMat(svd->ds,DS_MAT_U,&U);CHKERRQ(ierr); ierr = BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k+l);CHKERRQ(ierr); ierr = MatDestroy(&U);CHKERRQ(ierr); /* copy the last vector to be the next initial vector */ if (svd->reason == SVD_CONVERGED_ITERATING) { ierr = BVCopyColumn(svd->V,nv,svd->nconv+k+l);CHKERRQ(ierr); } svd->nconv += k; ierr = SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);CHKERRQ(ierr); } /* orthonormalize U columns in one side method */ if (lanczos->oneside) { for (i=0;i<svd->nconv;i++) { ierr = BVOrthogonalizeColumn(svd->U,i,NULL,&norm,NULL);CHKERRQ(ierr); ierr = BVScaleColumn(svd->U,i,1.0/norm);CHKERRQ(ierr); } } /* free working space */ ierr = PetscFree(w);CHKERRQ(ierr); if (swork) { ierr = PetscFree(swork);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode EPSSolve_Arnoldi(EPS eps) { PetscErrorCode ierr; PetscInt k,nv,ld; Mat U; PetscScalar *H,*X; PetscReal beta,gamma=1.0; PetscBool breakdown,harmonic,refined; BVOrthogRefineType orthog_ref; EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data; PetscFunctionBegin; ierr = DSGetLeadingDimension(eps->ds,&ld);CHKERRQ(ierr); ierr = DSGetRefined(eps->ds,&refined);CHKERRQ(ierr); harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE; ierr = BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL);CHKERRQ(ierr); /* Get the starting Arnoldi vector */ ierr = EPSGetStartVector(eps,0,NULL);CHKERRQ(ierr); /* Restart loop */ while (eps->reason == EPS_CONVERGED_ITERATING) { eps->its++; /* Compute an nv-step Arnoldi factorization */ nv = PetscMin(eps->nconv+eps->mpd,eps->ncv); ierr = DSSetDimensions(eps->ds,nv,0,eps->nconv,0);CHKERRQ(ierr); ierr = DSGetArray(eps->ds,DS_MAT_A,&H);CHKERRQ(ierr); if (!arnoldi->delayed) { ierr = EPSBasicArnoldi(eps,PETSC_FALSE,H,ld,eps->nconv,&nv,&beta,&breakdown);CHKERRQ(ierr); } else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Not implemented"); /*if (orthog_ref == BV_ORTHOG_REFINE_NEVER) { ierr = EPSDelayedArnoldi1(eps,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);CHKERRQ(ierr); } else { ierr = EPSDelayedArnoldi(eps,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);CHKERRQ(ierr); }*/ ierr = DSRestoreArray(eps->ds,DS_MAT_A,&H);CHKERRQ(ierr); ierr = DSSetState(eps->ds,DS_STATE_INTERMEDIATE);CHKERRQ(ierr); ierr = BVSetActiveColumns(eps->V,eps->nconv,nv);CHKERRQ(ierr); /* Compute translation of Krylov decomposition if harmonic extraction used */ if (harmonic) { ierr = DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,NULL,&gamma);CHKERRQ(ierr); } /* Solve projected problem */ ierr = DSSolve(eps->ds,eps->eigr,eps->eigi);CHKERRQ(ierr); ierr = DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);CHKERRQ(ierr); ierr = DSUpdateExtraRow(eps->ds);CHKERRQ(ierr); /* Check convergence */ ierr = EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);CHKERRQ(ierr); if (refined) { ierr = DSGetArray(eps->ds,DS_MAT_X,&X);CHKERRQ(ierr); ierr = BVMultColumn(eps->V,1.0,0.0,k,X+k*ld);CHKERRQ(ierr); ierr = DSRestoreArray(eps->ds,DS_MAT_X,&X);CHKERRQ(ierr); ierr = DSGetMat(eps->ds,DS_MAT_Q,&U);CHKERRQ(ierr); ierr = BVMultInPlace(eps->V,U,eps->nconv,nv);CHKERRQ(ierr); ierr = MatDestroy(&U);CHKERRQ(ierr); ierr = BVOrthogonalizeColumn(eps->V,k,NULL,NULL,NULL);CHKERRQ(ierr); } else { ierr = DSGetMat(eps->ds,DS_MAT_Q,&U);CHKERRQ(ierr); ierr = BVMultInPlace(eps->V,U,eps->nconv,nv);CHKERRQ(ierr); ierr = MatDestroy(&U);CHKERRQ(ierr); } eps->nconv = k; ierr = EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);CHKERRQ(ierr); if (breakdown && k<eps->nev) { ierr = PetscInfo2(eps,"Breakdown in Arnoldi method (it=%D norm=%g)\n",eps->its,(double)beta);CHKERRQ(ierr); ierr = EPSGetStartVector(eps,k,&breakdown);CHKERRQ(ierr); if (breakdown) { eps->reason = EPS_DIVERGED_BREAKDOWN; ierr = PetscInfo(eps,"Unable to generate more start vectors\n");CHKERRQ(ierr); } } if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS; if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL; } /* truncate Schur decomposition and change the state to raw so that PSVectors() computes eigenvectors from scratch */ ierr = DSSetDimensions(eps->ds,eps->nconv,0,0,0);CHKERRQ(ierr); ierr = DSSetState(eps->ds,DS_STATE_RAW);CHKERRQ(ierr); PetscFunctionReturn(0); }
static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work) { PetscErrorCode ierr; PetscReal a,b,eta; PetscScalar gamma; PetscInt i,j,k=nconv+l; Vec ui,ui1,vi; BVOrthogRefineType refine; PetscFunctionBegin; ierr = BVGetColumn(V,k,&vi);CHKERRQ(ierr); ierr = BVGetColumn(U,k,&ui);CHKERRQ(ierr); ierr = SVDMatMult(svd,PETSC_FALSE,vi,ui);CHKERRQ(ierr); ierr = BVRestoreColumn(V,k,&vi);CHKERRQ(ierr); ierr = BVRestoreColumn(U,k,&ui);CHKERRQ(ierr); if (l>0) { ierr = BVMultColumn(U,-1.0,1.0,k,&gamma);CHKERRQ(ierr); beta[nconv] = PetscRealPart(gamma); } ierr = BVGetOrthogonalization(V,NULL,&refine,&eta);CHKERRQ(ierr); for (i=k+1;i<n;i++) { ierr = BVGetColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVGetColumn(U,i-1,&ui1);CHKERRQ(ierr); ierr = SVDMatMult(svd,PETSC_TRUE,ui1,vi);CHKERRQ(ierr); ierr = BVRestoreColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVRestoreColumn(U,i-1,&ui1);CHKERRQ(ierr); ierr = BVNormColumn(U,i-1,NORM_2,&a);CHKERRQ(ierr); if (refine == BV_ORTHOG_REFINE_IFNEEDED) { ierr = BVSetActiveColumns(V,0,i+1);CHKERRQ(ierr); ierr = BVGetColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVDotVec(V,vi,work);CHKERRQ(ierr); ierr = BVRestoreColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVSetActiveColumns(V,0,i);CHKERRQ(ierr); } else { ierr = BVSetActiveColumns(V,0,i);CHKERRQ(ierr); ierr = BVDotColumn(V,i,work);CHKERRQ(ierr); } ierr = BVScaleColumn(U,i-1,1.0/a);CHKERRQ(ierr); for (j=0;j<i;j++) work[j] = work[j] / a; ierr = BVMultColumn(V,-1.0,1.0/a,i,work);CHKERRQ(ierr); ierr = SVDOrthogonalizeCGS(V,i,work,a,refine,eta,&b);CHKERRQ(ierr); ierr = BVScaleColumn(V,i,1.0/b);CHKERRQ(ierr); ierr = BVGetColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVGetColumn(U,i,&ui);CHKERRQ(ierr); ierr = BVGetColumn(U,i-1,&ui1);CHKERRQ(ierr); ierr = SVDMatMult(svd,PETSC_FALSE,vi,ui);CHKERRQ(ierr); ierr = VecAXPY(ui,-b,ui1);CHKERRQ(ierr); ierr = BVRestoreColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVRestoreColumn(U,i,&ui);CHKERRQ(ierr); ierr = BVRestoreColumn(U,i-1,&ui1);CHKERRQ(ierr); alpha[i-1] = a; beta[i-1] = b; } ierr = BVGetColumn(V,n,&vi);CHKERRQ(ierr); ierr = BVGetColumn(U,n-1,&ui1);CHKERRQ(ierr); ierr = SVDMatMult(svd,PETSC_TRUE,ui1,vi);CHKERRQ(ierr); ierr = BVRestoreColumn(V,n,&vi);CHKERRQ(ierr); ierr = BVRestoreColumn(U,n-1,&ui1);CHKERRQ(ierr); ierr = BVNormColumn(svd->U,n-1,NORM_2,&a);CHKERRQ(ierr); if (refine == BV_ORTHOG_REFINE_IFNEEDED) { ierr = BVSetActiveColumns(V,0,n+1);CHKERRQ(ierr); ierr = BVGetColumn(V,n,&vi);CHKERRQ(ierr); ierr = BVDotVec(V,vi,work);CHKERRQ(ierr); ierr = BVRestoreColumn(V,n,&vi);CHKERRQ(ierr); } else { ierr = BVSetActiveColumns(V,0,n);CHKERRQ(ierr); ierr = BVDotColumn(V,n,work);CHKERRQ(ierr); } ierr = BVScaleColumn(U,n-1,1.0/a);CHKERRQ(ierr); for (j=0;j<n;j++) work[j] = work[j] / a; ierr = BVMultColumn(V,-1.0,1.0/a,n,work);CHKERRQ(ierr); ierr = SVDOrthogonalizeCGS(V,n,work,a,refine,eta,&b);CHKERRQ(ierr); ierr = BVSetActiveColumns(V,nconv,n);CHKERRQ(ierr); alpha[n-1] = a; beta[n-1] = b; PetscFunctionReturn(0); }