/* 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); }
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; }
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); }