/* BVOrthogonalizeCGS1 - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization */ PetscErrorCode BVOrthogonalizeCGS1(BV bv,PetscInt j,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm) { PetscErrorCode ierr; PetscInt i; PetscReal sum,nrm,beta; Vec w=v; PetscFunctionBegin; /* h = W^* v ; alpha = (v, v) */ bv->k = j; if (onorm || norm) { if (!v) { bv->k++; ierr = BVGetColumn(bv,j,&w);CHKERRQ(ierr); } ierr = BVDotVec(bv,w,H);CHKERRQ(ierr); if (!v) { ierr = BVRestoreColumn(bv,j,&w);CHKERRQ(ierr); bv->k--; beta = PetscSqrtReal(PetscRealPart(H[bv->nc+j])); } else { ierr = BVNormVec(bv,w,NORM_2,&beta);CHKERRQ(ierr); } } else { if (!v) { ierr = BVDotColumn(bv,j,H);CHKERRQ(ierr); } else { ierr = BVDotVec(bv,w,H);CHKERRQ(ierr); } } /* q = v - V h */ if (bv->indef) { for (i=0;i<bv->nc+j;i++) H[i] /= bv->omega[i]; /* apply inverse of signature */ } if (!v) { ierr = BVMultColumn(bv,-1.0,1.0,j,H);CHKERRQ(ierr); } else { ierr = BVMultVec(bv,-1.0,1.0,w,H);CHKERRQ(ierr); } if (bv->indef) { for (i=0;i<bv->nc+j;i++) H[i] *= bv->omega[i]; /* revert signature */ } /* compute |v| */ if (onorm) *onorm = beta; if (bv->indef) { if (!v) { ierr = BVNormColumn(bv,j,NORM_2,&nrm);CHKERRQ(ierr); } else { ierr = BVNormVec(bv,w,NORM_2,&nrm);CHKERRQ(ierr); } if (norm) *norm = nrm; bv->omega[bv->nc+j] = (nrm<0.0)? -1.0: 1.0; } else if (norm) { /* estimate |v'| from |v| */ sum = 0.0; for (i=0;i<bv->nc+j;i++) sum += PetscRealPart(H[i]*PetscConj(H[i])); *norm = beta*beta-sum; if (*norm <= 0.0) { if (!v) { ierr = BVNormColumn(bv,j,NORM_2,norm);CHKERRQ(ierr); } else { ierr = BVNormVec(bv,w,NORM_2,norm);CHKERRQ(ierr); } } else *norm = PetscSqrtReal(*norm); } PetscFunctionReturn(0); }
static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n) { PetscErrorCode ierr; PetscReal a,b; PetscScalar gamma; PetscInt i,k=nconv+l; Vec ui,ui1,vi; 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 = BVNormColumn(U,k,NORM_2,&a);CHKERRQ(ierr); ierr = BVScaleColumn(U,k,1.0/a);CHKERRQ(ierr); alpha[k] = a; 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 = BVOrthogonalizeColumn(V,i,NULL,&b,NULL);CHKERRQ(ierr); ierr = BVScaleColumn(V,i,1.0/b);CHKERRQ(ierr); beta[i-1] = b; ierr = BVGetColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVGetColumn(U,i,&ui);CHKERRQ(ierr); ierr = SVDMatMult(svd,PETSC_FALSE,vi,ui);CHKERRQ(ierr); ierr = BVRestoreColumn(V,i,&vi);CHKERRQ(ierr); ierr = BVGetColumn(U,i-1,&ui1);CHKERRQ(ierr); ierr = VecAXPY(ui,-b,ui1);CHKERRQ(ierr); ierr = BVRestoreColumn(U,i-1,&ui1);CHKERRQ(ierr); ierr = BVRestoreColumn(U,i,&ui);CHKERRQ(ierr); ierr = BVNormColumn(U,i,NORM_2,&a);CHKERRQ(ierr); ierr = BVScaleColumn(U,i,1.0/a);CHKERRQ(ierr); alpha[i] = a; } 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 = BVOrthogonalizeColumn(V,n,NULL,&b,NULL);CHKERRQ(ierr); beta[n-1] = b; PetscFunctionReturn(0); }
/* BVOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt */ static PetscErrorCode BVOrthogonalizeCGS(BV bv,PetscInt j,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep) { PetscErrorCode ierr; PetscReal onrm,nrm; PetscInt i,k,l; PetscFunctionBegin; if (v) k = bv->k; else k = j; switch (bv->orthog_ref) { case BV_ORTHOG_REFINE_IFNEEDED: ierr = BVOrthogonalizeCGS1(bv,k,v,bv->h,&onrm,&nrm);CHKERRQ(ierr); /* ||q|| < eta ||h|| */ l = 1; while (l<3 && nrm && nrm < bv->orthog_eta*onrm) { l++; ierr = BVOrthogonalizeCGS1(bv,k,v,bv->c,&onrm,&nrm);CHKERRQ(ierr); for (i=0;i<bv->nc+k;i++) bv->h[i] += bv->c[i]; } if (norm) *norm = nrm; if (lindep) { if (nrm < bv->orthog_eta*onrm) *lindep = PETSC_TRUE; else *lindep = PETSC_FALSE; } break; case BV_ORTHOG_REFINE_NEVER: ierr = BVOrthogonalizeCGS1(bv,k,v,bv->h,NULL,NULL);CHKERRQ(ierr); /* compute |v| */ if (norm || lindep) { if (v) { ierr = BVNormVec(bv,v,NORM_2,&nrm);CHKERRQ(ierr); } else { ierr = BVNormColumn(bv,k,NORM_2,&nrm);CHKERRQ(ierr); } } if (norm) *norm = nrm; /* linear dependence check: just test for exactly zero norm */ if (lindep) *lindep = nrm? PETSC_FALSE: PETSC_TRUE; break; case BV_ORTHOG_REFINE_ALWAYS: ierr = BVOrthogonalizeCGS1(bv,k,v,bv->h,NULL,NULL);CHKERRQ(ierr); if (lindep) { ierr = BVOrthogonalizeCGS1(bv,k,v,bv->c,&onrm,&nrm);CHKERRQ(ierr); if (norm) *norm = nrm; if (nrm==0.0 || nrm < bv->orthog_eta*onrm) *lindep = PETSC_TRUE; else *lindep = PETSC_FALSE; } else { ierr = BVOrthogonalizeCGS1(bv,k,v,bv->c,NULL,norm);CHKERRQ(ierr); } for (i=0;i<bv->nc+k;i++) bv->h[i] += bv->c[i]; break; } PetscFunctionReturn(0); }
/* Custom CGS orthogonalization, preprocess after first orthogonalization */ static PetscErrorCode SVDOrthogonalizeCGS(BV V,PetscInt i,PetscScalar* h,PetscReal a,BVOrthogRefineType refine,PetscReal eta,PetscReal *norm) { PetscErrorCode ierr; PetscReal sum,onorm; PetscScalar dot; PetscInt j; PetscFunctionBegin; switch (refine) { case BV_ORTHOG_REFINE_NEVER: ierr = BVNormColumn(V,i,NORM_2,norm);CHKERRQ(ierr); break; case BV_ORTHOG_REFINE_ALWAYS: ierr = BVSetActiveColumns(V,0,i);CHKERRQ(ierr); ierr = BVDotColumn(V,i,h);CHKERRQ(ierr); ierr = BVMultColumn(V,-1.0,1.0,i,h);CHKERRQ(ierr); ierr = BVNormColumn(V,i,NORM_2,norm);CHKERRQ(ierr); break; case BV_ORTHOG_REFINE_IFNEEDED: dot = h[i]; onorm = PetscSqrtReal(PetscRealPart(dot)) / a; sum = 0.0; for (j=0;j<i;j++) { sum += PetscRealPart(h[j] * PetscConj(h[j])); } *norm = PetscRealPart(dot)/(a*a) - sum; if (*norm>0.0) *norm = PetscSqrtReal(*norm); else { ierr = BVNormColumn(V,i,NORM_2,norm);CHKERRQ(ierr); } if (*norm < eta*onorm) { ierr = BVSetActiveColumns(V,0,i);CHKERRQ(ierr); ierr = BVDotColumn(V,i,h);CHKERRQ(ierr); ierr = BVMultColumn(V,-1.0,1.0,i,h);CHKERRQ(ierr); ierr = BVNormColumn(V,i,NORM_2,norm);CHKERRQ(ierr); } break; } PetscFunctionReturn(0); }
int main(int argc,char **argv) { PetscErrorCode ierr; Vec t,v; Mat Q,M; BV X,Y; PetscInt i,j,n=10,kx=6,lx=3,ky=5,ly=2; PetscScalar *q,*z; PetscReal nrm; PetscViewer view; PetscBool verbose,trans; SlepcInitialize(&argc,&argv,(char*)0,help); ierr = PetscOptionsGetInt(NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,"-kx",&kx,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,"-lx",&lx,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,"-ky",&ky,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,"-ly",&ly,NULL);CHKERRQ(ierr); ierr = PetscOptionsHasName(NULL,"-verbose",&verbose);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"First BV with %D active columns (%D leading columns) of dimension %D.\n",kx,lx,n);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Second BV with %D active columns (%D leading columns) of dimension %D.\n",ky,ly,n);CHKERRQ(ierr); /* Create template vector */ ierr = VecCreate(PETSC_COMM_WORLD,&t);CHKERRQ(ierr); ierr = VecSetSizes(t,PETSC_DECIDE,n);CHKERRQ(ierr); ierr = VecSetFromOptions(t);CHKERRQ(ierr); /* Create BV object X */ ierr = BVCreate(PETSC_COMM_WORLD,&X);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)X,"X");CHKERRQ(ierr); ierr = BVSetSizesFromVec(X,t,kx+2);CHKERRQ(ierr); /* two extra columns to test active columns */ ierr = BVSetFromOptions(X);CHKERRQ(ierr); ierr = BVSetActiveColumns(X,lx,kx);CHKERRQ(ierr); /* Set up viewer */ ierr = PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);CHKERRQ(ierr); if (verbose) { ierr = PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); } /* Fill X entries */ for (j=0;j<kx+2;j++) { ierr = BVGetColumn(X,j,&v);CHKERRQ(ierr); ierr = VecZeroEntries(v);CHKERRQ(ierr); for (i=0;i<4;i++) { if (i+j<n) { ierr = VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);CHKERRQ(ierr); } } ierr = VecAssemblyBegin(v);CHKERRQ(ierr); ierr = VecAssemblyEnd(v);CHKERRQ(ierr); ierr = BVRestoreColumn(X,j,&v);CHKERRQ(ierr); } if (verbose) { ierr = BVView(X,view);CHKERRQ(ierr); } /* Create BV object Y */ ierr = BVCreate(PETSC_COMM_WORLD,&Y);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Y,"Y");CHKERRQ(ierr); ierr = BVSetSizesFromVec(Y,t,ky+1);CHKERRQ(ierr); ierr = BVSetFromOptions(Y);CHKERRQ(ierr); ierr = BVSetActiveColumns(Y,ly,ky);CHKERRQ(ierr); /* Fill Y entries */ for (j=0;j<ky+1;j++) { ierr = BVGetColumn(Y,j,&v);CHKERRQ(ierr); ierr = VecSet(v,(PetscScalar)(j+1)/4.0);CHKERRQ(ierr); ierr = BVRestoreColumn(Y,j,&v);CHKERRQ(ierr); } if (verbose) { ierr = BVView(Y,view);CHKERRQ(ierr); } /* Create Mat */ ierr = MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Q,"Q");CHKERRQ(ierr); ierr = MatDenseGetArray(Q,&q);CHKERRQ(ierr); for (i=0;i<kx;i++) for (j=0;j<ky;j++) q[i+j*kx] = (i<j)? 2.0: -0.5; ierr = MatDenseRestoreArray(Q,&q);CHKERRQ(ierr); if (verbose) { ierr = MatView(Q,NULL);CHKERRQ(ierr); } /* Test BVMult */ ierr = BVMult(Y,2.0,1.0,X,Q);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");CHKERRQ(ierr); ierr = BVView(Y,view);CHKERRQ(ierr); } /* Test BVMultVec */ ierr = BVGetColumn(Y,0,&v);CHKERRQ(ierr); ierr = PetscMalloc1(kx-lx,&z);CHKERRQ(ierr); z[0] = 2.0; for (i=1;i<kx-lx;i++) z[i] = -0.5*z[i-1]; ierr = BVMultVec(X,-1.0,1.0,v,z);CHKERRQ(ierr); ierr = PetscFree(z);CHKERRQ(ierr); ierr = BVRestoreColumn(Y,0,&v);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");CHKERRQ(ierr); ierr = BVView(Y,view);CHKERRQ(ierr); } /* Test BVDot */ ierr = MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)M,"M");CHKERRQ(ierr); ierr = BVDot(X,Y,M);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");CHKERRQ(ierr); ierr = MatView(M,NULL);CHKERRQ(ierr); } /* Test BVDotVec */ ierr = BVGetColumn(Y,0,&v);CHKERRQ(ierr); ierr = PetscMalloc1(kx-lx,&z);CHKERRQ(ierr); ierr = BVDotVec(X,v,z);CHKERRQ(ierr); ierr = BVRestoreColumn(Y,0,&v);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");CHKERRQ(ierr); ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,kx-lx,z,&v);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)v,"z");CHKERRQ(ierr); ierr = VecView(v,view);CHKERRQ(ierr); ierr = VecDestroy(&v);CHKERRQ(ierr); } ierr = PetscFree(z);CHKERRQ(ierr); /* Test BVMultInPlace and BVScale */ ierr = PetscOptionsHasName(NULL,"-trans",&trans);CHKERRQ(ierr); if (trans) { Mat Qt; ierr = MatTranspose(Q,MAT_INITIAL_MATRIX,&Qt);CHKERRQ(ierr); ierr = BVMultInPlaceTranspose(X,Qt,lx+1,ky);CHKERRQ(ierr); ierr = MatDestroy(&Qt);CHKERRQ(ierr); } else { ierr = BVMultInPlace(X,Q,lx+1,ky);CHKERRQ(ierr); } ierr = BVScale(X,2.0);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");CHKERRQ(ierr); ierr = BVView(X,view);CHKERRQ(ierr); } /* Test BVNorm */ ierr = BVNormColumn(X,lx,NORM_2,&nrm);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"2-Norm or X[%D] = %g\n",lx,(double)nrm);CHKERRQ(ierr); ierr = BVNorm(X,NORM_FROBENIUS,&nrm);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm or X = %g\n",(double)nrm);CHKERRQ(ierr); ierr = BVDestroy(&X);CHKERRQ(ierr); ierr = BVDestroy(&Y);CHKERRQ(ierr); ierr = MatDestroy(&Q);CHKERRQ(ierr); ierr = MatDestroy(&M);CHKERRQ(ierr); ierr = VecDestroy(&t);CHKERRQ(ierr); ierr = SlepcFinalize(); return 0; }
int main(int argc,char **argv) { PetscErrorCode ierr; Vec t,v; Mat B,M; BV X; PetscInt i,j,n=10,k=5,Istart,Iend,col[3]; PetscScalar value[3],alpha; PetscReal nrm; PetscViewer view; PetscBool verbose,FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; SlepcInitialize(&argc,&argv,(char*)0,help); ierr = PetscOptionsGetInt(NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,"-k",&k,NULL);CHKERRQ(ierr); ierr = PetscOptionsHasName(NULL,"-verbose",&verbose);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Test BV with non-standard inner product (n=%D, k=%D).\n",n,k);CHKERRQ(ierr); /* Create inner product matrix */ ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); ierr = MatSetFromOptions(B);CHKERRQ(ierr); ierr = MatSetUp(B);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)B,"B");CHKERRQ(ierr); ierr = MatGetOwnershipRange(B,&Istart,&Iend);CHKERRQ(ierr); if (Istart==0) FirstBlock=PETSC_TRUE; if (Iend==n) LastBlock=PETSC_TRUE; value[0]=-1.0; value[1]=2.0; value[2]=-1.0; for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) { col[0]=i-1; col[1]=i; col[2]=i+1; ierr = MatSetValues(B,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); } if (LastBlock) { i=n-1; col[0]=n-2; col[1]=n-1; ierr = MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); } if (FirstBlock) { i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0; ierr = MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); } ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatGetVecs(B,&t,NULL);CHKERRQ(ierr); /* Create BV object X */ ierr = BVCreate(PETSC_COMM_WORLD,&X);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)X,"X");CHKERRQ(ierr); ierr = BVSetSizesFromVec(X,t,k);CHKERRQ(ierr); ierr = BVSetFromOptions(X);CHKERRQ(ierr); ierr = BVSetMatrix(X,B,PETSC_FALSE);CHKERRQ(ierr); /* Set up viewer */ ierr = PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);CHKERRQ(ierr); if (verbose) { ierr = PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); } /* Fill X entries */ for (j=0;j<k;j++) { ierr = BVGetColumn(X,j,&v);CHKERRQ(ierr); ierr = VecZeroEntries(v);CHKERRQ(ierr); for (i=0;i<4;i++) { if (i+j<n) { ierr = VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);CHKERRQ(ierr); } } ierr = VecAssemblyBegin(v);CHKERRQ(ierr); ierr = VecAssemblyEnd(v);CHKERRQ(ierr); ierr = BVRestoreColumn(X,j,&v);CHKERRQ(ierr); } if (verbose) { ierr = MatView(B,view);CHKERRQ(ierr); ierr = BVView(X,view);CHKERRQ(ierr); } /* Test BVNormColumn */ ierr = BVNormColumn(X,0,NORM_2,&nrm);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"B-Norm or X[0] = %g\n",(double)nrm);CHKERRQ(ierr); /* Test BVOrthogonalizeColumn */ for (j=0;j<k;j++) { ierr = BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL);CHKERRQ(ierr); alpha = 1.0/nrm; ierr = BVScaleColumn(X,j,alpha);CHKERRQ(ierr); } if (verbose) { ierr = BVView(X,view);CHKERRQ(ierr); } /* Check orthogonality */ ierr = MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);CHKERRQ(ierr); ierr = BVDot(X,X,M);CHKERRQ(ierr); ierr = MatShift(M,-1.0);CHKERRQ(ierr); ierr = MatNorm(M,NORM_1,&nrm);CHKERRQ(ierr); if (nrm<100*PETSC_MACHINE_EPSILON) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");CHKERRQ(ierr); } else { ierr = PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm);CHKERRQ(ierr); } ierr = BVDestroy(&X);CHKERRQ(ierr); ierr = MatDestroy(&M);CHKERRQ(ierr); ierr = MatDestroy(&B);CHKERRQ(ierr); ierr = VecDestroy(&t);CHKERRQ(ierr); ierr = SlepcFinalize(); return 0; }
PetscErrorCode NEPSolve_NArnoldi(NEP nep) { PetscErrorCode ierr; Mat T=nep->function,Tsigma; Vec f,r=nep->work[0],x=nep->work[1],w=nep->work[2]; PetscScalar *X,lambda; PetscReal beta,resnorm=0.0,nrm; PetscInt n; PetscBool breakdown; KSPConvergedReason kspreason; PetscFunctionBegin; /* get initial space and shift */ ierr = NEPGetDefaultShift(nep,&lambda);CHKERRQ(ierr); if (!nep->nini) { ierr = BVSetRandomColumn(nep->V,0,nep->rand);CHKERRQ(ierr); ierr = BVNormColumn(nep->V,0,NORM_2,&nrm);CHKERRQ(ierr); ierr = BVScaleColumn(nep->V,0,1.0/nrm);CHKERRQ(ierr); n = 1; } else n = nep->nini; /* build projected matrices for initial space */ ierr = DSSetDimensions(nep->ds,n,0,0,0);CHKERRQ(ierr); ierr = NEPProjectOperator(nep,0,n);CHKERRQ(ierr); /* prepare linear solver */ ierr = NEPComputeFunction(nep,lambda,T,T);CHKERRQ(ierr); ierr = MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);CHKERRQ(ierr); ierr = KSPSetOperators(nep->ksp,Tsigma,Tsigma);CHKERRQ(ierr); /* Restart loop */ while (nep->reason == NEP_CONVERGED_ITERATING) { nep->its++; /* solve projected problem */ ierr = DSSetDimensions(nep->ds,n,0,0,0);CHKERRQ(ierr); ierr = DSSetState(nep->ds,DS_STATE_RAW);CHKERRQ(ierr); ierr = DSSolve(nep->ds,nep->eigr,NULL);CHKERRQ(ierr); lambda = nep->eigr[0]; /* compute Ritz vector, x = V*s */ ierr = DSGetArray(nep->ds,DS_MAT_X,&X);CHKERRQ(ierr); ierr = BVSetActiveColumns(nep->V,0,n);CHKERRQ(ierr); ierr = BVMultVec(nep->V,1.0,0.0,x,X);CHKERRQ(ierr); ierr = DSRestoreArray(nep->ds,DS_MAT_X,&X);CHKERRQ(ierr); /* compute the residual, r = T(lambda)*x */ ierr = NEPApplyFunction(nep,lambda,x,w,r,NULL,NULL);CHKERRQ(ierr); /* convergence test */ ierr = VecNorm(r,NORM_2,&resnorm);CHKERRQ(ierr); nep->errest[nep->nconv] = resnorm; if (resnorm<=nep->rtol) { ierr = BVInsertVec(nep->V,nep->nconv,x);CHKERRQ(ierr); nep->nconv = nep->nconv + 1; nep->reason = NEP_CONVERGED_FNORM_RELATIVE; } ierr = NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->errest,1);CHKERRQ(ierr); if (nep->reason == NEP_CONVERGED_ITERATING) { /* continuation vector: f = T(sigma)\r */ ierr = BVGetColumn(nep->V,n,&f);CHKERRQ(ierr); ierr = NEP_KSPSolve(nep,r,f);CHKERRQ(ierr); ierr = BVRestoreColumn(nep->V,n,&f);CHKERRQ(ierr); ierr = KSPGetConvergedReason(nep->ksp,&kspreason);CHKERRQ(ierr); if (kspreason<0) { ierr = PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);CHKERRQ(ierr); nep->reason = NEP_DIVERGED_LINEAR_SOLVE; break; } /* orthonormalize */ ierr = BVOrthogonalizeColumn(nep->V,n,NULL,&beta,&breakdown);CHKERRQ(ierr); if (breakdown || beta==0.0) { ierr = PetscInfo1(nep,"iter=%D, orthogonalization failed, stopping solve\n",nep->its);CHKERRQ(ierr); nep->reason = NEP_DIVERGED_BREAKDOWN; break; } ierr = BVScaleColumn(nep->V,n,1.0/beta);CHKERRQ(ierr); /* update projected matrices */ ierr = DSSetDimensions(nep->ds,n+1,0,0,0);CHKERRQ(ierr); ierr = NEPProjectOperator(nep,n,n+1);CHKERRQ(ierr); n++; } if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT; } ierr = MatDestroy(&Tsigma);CHKERRQ(ierr); 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); }
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); }