/* PEPConvergedNormRelative - Checks convergence relative to the matrix norms. */ PetscErrorCode PEPConvergedNormRelative(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx) { PetscErrorCode ierr; PetscInt i; PetscReal w=0.0; PetscScalar t[20],*vals=t,*ivals=NULL; #if !defined(PETSC_USE_COMPLEX) PetscScalar it[20]; #endif PetscFunctionBegin; #if !defined(PETSC_USE_COMPLEX) ivals = it; #endif if (pep->nmat>20) { #if !defined(PETSC_USE_COMPLEX) ierr = PetscMalloc2(pep->nmat,&vals,pep->nmat,&ivals);CHKERRQ(ierr); #else ierr = PetscMalloc1(pep->nmat,&vals);CHKERRQ(ierr); #endif } ierr = PEPEvaluateBasis(pep,eigr,eigi,vals,ivals);CHKERRQ(ierr); for (i=0;i<pep->nmat;i++) w += SlepcAbsEigenvalue(vals[i],ivals[i])*pep->nrma[i]; *errest = res/w; if (pep->nmat>20) { #if !defined(PETSC_USE_COMPLEX) ierr = PetscFree2(vals,ivals);CHKERRQ(ierr); #else ierr = PetscFree(vals);CHKERRQ(ierr); #endif } PetscFunctionReturn(0); }
/* PEPComputeResidualNorm_Private - Computes the norm of the residual vector associated with an eigenpair. */ PetscErrorCode PEPComputeResidualNorm_Private(PEP pep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm) { PetscErrorCode ierr; Vec u,w; Mat *A=pep->A; PetscInt i,nmat=pep->nmat; PetscScalar t[20],*vals=t,*ivals=NULL; #if !defined(PETSC_USE_COMPLEX) Vec ui,wi; PetscReal ni; PetscBool imag; PetscScalar it[20]; #endif PetscFunctionBegin; ierr = BVGetVec(pep->V,&u);CHKERRQ(ierr); ierr = BVGetVec(pep->V,&w);CHKERRQ(ierr); ierr = VecZeroEntries(u);CHKERRQ(ierr); #if !defined(PETSC_USE_COMPLEX) ivals = it; #endif if (nmat>20) { ierr = PetscMalloc(nmat*sizeof(PetscScalar),&vals);CHKERRQ(ierr); #if !defined(PETSC_USE_COMPLEX) ierr = PetscMalloc(nmat*sizeof(PetscScalar),&ivals);CHKERRQ(ierr); #endif } ierr = PEPEvaluateBasis(pep,kr,ki,vals,ivals);CHKERRQ(ierr); #if !defined(PETSC_USE_COMPLEX) if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) imag = PETSC_FALSE; else { imag = PETSC_TRUE; ierr = VecDuplicate(u,&ui);CHKERRQ(ierr); ierr = VecDuplicate(u,&wi);CHKERRQ(ierr); ierr = VecZeroEntries(ui);CHKERRQ(ierr); } #endif for (i=0;i<nmat;i++) { if (vals[i]!=0.0) { ierr = MatMult(A[i],xr,w);CHKERRQ(ierr); ierr = VecAXPY(u,vals[i],w);CHKERRQ(ierr); } #if !defined(PETSC_USE_COMPLEX) if (imag) { if (ivals[i]!=0 || vals[i]!=0) { ierr = MatMult(A[i],xi,wi);CHKERRQ(ierr); if (vals[i]==0) { ierr = MatMult(A[i],xr,w);CHKERRQ(ierr); } } if (ivals[i]!=0){ ierr = VecAXPY(u,-ivals[i],wi);CHKERRQ(ierr); ierr = VecAXPY(ui,ivals[i],w);CHKERRQ(ierr); } if (vals[i]!=0) { ierr = VecAXPY(ui,vals[i],wi);CHKERRQ(ierr); } } #endif } ierr = VecNorm(u,NORM_2,norm);CHKERRQ(ierr); #if !defined(PETSC_USE_COMPLEX) if (imag) { ierr = VecNorm(ui,NORM_2,&ni);CHKERRQ(ierr); *norm = SlepcAbsEigenvalue(*norm,ni); } #endif ierr = VecDestroy(&w);CHKERRQ(ierr); ierr = VecDestroy(&u);CHKERRQ(ierr); if (nmat>20) { ierr = PetscFree(vals);CHKERRQ(ierr); #if !defined(PETSC_USE_COMPLEX) ierr = PetscFree(ivals);CHKERRQ(ierr); #endif } #if !defined(PETSC_USE_COMPLEX) if (imag) { ierr = VecDestroy(&wi);CHKERRQ(ierr); ierr = VecDestroy(&ui);CHKERRQ(ierr); } #endif PetscFunctionReturn(0); }