Ejemplo n.º 1
0
/*
   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);
}
Ejemplo n.º 2
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);
}
Ejemplo n.º 3
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);
}
Ejemplo n.º 4
0
PetscErrorCode EPSSolve_Lanczos(EPS eps)
{
  EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
  PetscErrorCode ierr;
  PetscInt       nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
  Vec            vi,vj,w;
  Mat            U;
  PetscScalar    *Y,*ritz,stmp;
  PetscReal      *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
  PetscBool      breakdown;
  char           *conv,ctmp;

  PetscFunctionBegin;
  ierr = DSGetLeadingDimension(eps->ds,&ld);CHKERRQ(ierr);
  ierr = PetscMalloc4(ncv,&ritz,ncv,&bnd,ncv,&perm,ncv,&conv);CHKERRQ(ierr);

  /* The first Lanczos vector is the normalized initial vector */
  ierr = EPSGetStartVector(eps,0,NULL);CHKERRQ(ierr);

  anorm = -1.0;
  nconv = 0;

  /* Restart loop */
  while (eps->reason == EPS_CONVERGED_ITERATING) {
    eps->its++;

    /* Compute an ncv-step Lanczos factorization */
    n = PetscMin(nconv+eps->mpd,ncv);
    ierr = DSGetArrayReal(eps->ds,DS_MAT_T,&d);CHKERRQ(ierr);
    e = d + ld;
    ierr = EPSBasicLanczos(eps,d,e,nconv,&n,&breakdown,anorm);CHKERRQ(ierr);
    beta = e[n-1];
    ierr = DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);CHKERRQ(ierr);
    ierr = DSSetDimensions(eps->ds,n,0,nconv,0);CHKERRQ(ierr);
    ierr = DSSetState(eps->ds,DS_STATE_INTERMEDIATE);CHKERRQ(ierr);
    ierr = BVSetActiveColumns(eps->V,nconv,n);CHKERRQ(ierr);

    /* Solve projected problem */
    ierr = DSSolve(eps->ds,ritz,NULL);CHKERRQ(ierr);
    ierr = DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL);CHKERRQ(ierr);

    /* Estimate ||A|| */
    for (i=nconv;i<n;i++)
      anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));

    /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
    ierr = DSGetArray(eps->ds,DS_MAT_Q,&Y);CHKERRQ(ierr);
    for (i=nconv;i<n;i++) {
      resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
      ierr = (*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx);CHKERRQ(ierr);
      if (bnd[i]<eps->tol) conv[i] = 'C';
      else conv[i] = 'N';
    }
    ierr = DSRestoreArray(eps->ds,DS_MAT_Q,&Y);CHKERRQ(ierr);

    /* purge repeated ritz values */
    if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
      for (i=nconv+1;i<n;i++) {
        if (conv[i] == 'C' && PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol) conv[i] = 'R';
      }
    }

    /* Compute restart vector */
    if (breakdown) {
      ierr = PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%g)\n",eps->its,(double)beta);CHKERRQ(ierr);
    } else {
      restart = nconv;
      while (restart<n && conv[restart] != 'N') restart++;
      if (restart >= n) {
        breakdown = PETSC_TRUE;
      } else {
        for (i=restart+1;i<n;i++) {
          if (conv[i] == 'N') {
            ierr = SlepcSCCompare(eps->sc,ritz[restart],0.0,ritz[i],0.0,&r);CHKERRQ(ierr);
            if (r>0) restart = i;
          }
        }
        ierr = DSGetArray(eps->ds,DS_MAT_Q,&Y);CHKERRQ(ierr);
        ierr = BVMultColumn(eps->V,1.0,0.0,n,Y+restart*ld+nconv);CHKERRQ(ierr);
        ierr = DSRestoreArray(eps->ds,DS_MAT_Q,&Y);CHKERRQ(ierr);
      }
    }

    /* Count and put converged eigenvalues first */
    for (i=nconv;i<n;i++) perm[i] = i;
    for (k=nconv;k<n;k++) {
      if (conv[perm[k]] != 'C') {
        j = k + 1;
        while (j<n && conv[perm[j]] != 'C') j++;
        if (j>=n) break;
        l = perm[k]; perm[k] = perm[j]; perm[j] = l;
      }
    }

    /* Sort eigenvectors according to permutation */
    ierr = DSGetArray(eps->ds,DS_MAT_Q,&Y);CHKERRQ(ierr);
    for (i=nconv;i<k;i++) {
      x = perm[i];
      if (x != i) {
        j = i + 1;
        while (perm[j] != i) j++;
        /* swap eigenvalues i and j */
        stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
        rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
        ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
        perm[j] = x; perm[i] = i;
        /* swap eigenvectors i and j */
        for (l=0;l<n;l++) {
          stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
        }
      }
    }
    ierr = DSRestoreArray(eps->ds,DS_MAT_Q,&Y);CHKERRQ(ierr);

    /* compute converged eigenvectors */
    ierr = DSGetMat(eps->ds,DS_MAT_Q,&U);CHKERRQ(ierr);
    ierr = BVMultInPlace(eps->V,U,nconv,k);CHKERRQ(ierr);
    ierr = MatDestroy(&U);CHKERRQ(ierr);

    /* purge spurious ritz values */
    if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
      for (i=nconv;i<k;i++) {
        ierr = BVGetColumn(eps->V,i,&vi);CHKERRQ(ierr);
        ierr = VecNorm(vi,NORM_2,&norm);CHKERRQ(ierr);
        ierr = VecScale(vi,1.0/norm);CHKERRQ(ierr);
        w = eps->work[0];
        ierr = STApply(eps->st,vi,w);CHKERRQ(ierr);
        ierr = VecAXPY(w,-ritz[i],vi);CHKERRQ(ierr);
        ierr = BVRestoreColumn(eps->V,i,&vi);CHKERRQ(ierr);
        ierr = VecNorm(w,NORM_2,&norm);CHKERRQ(ierr);
        ierr = (*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx);CHKERRQ(ierr);
        if (bnd[i]>=eps->tol) conv[i] = 'S';
      }
      for (i=nconv;i<k;i++) {
        if (conv[i] != 'C') {
          j = i + 1;
          while (j<k && conv[j] != 'C') j++;
          if (j>=k) break;
          /* swap eigenvalues i and j */
          stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
          rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
          ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
          /* swap eigenvectors i and j */
          ierr = BVGetColumn(eps->V,i,&vi);CHKERRQ(ierr);
          ierr = BVGetColumn(eps->V,j,&vj);CHKERRQ(ierr);
          ierr = VecSwap(vi,vj);CHKERRQ(ierr);
          ierr = BVRestoreColumn(eps->V,i,&vi);CHKERRQ(ierr);
          ierr = BVRestoreColumn(eps->V,j,&vj);CHKERRQ(ierr);
        }
      }
      k = i;
    }

    /* store ritz values and estimated errors */
    for (i=nconv;i<n;i++) {
      eps->eigr[i] = ritz[i];
      eps->errest[i] = bnd[i];
    }
    ierr = EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);CHKERRQ(ierr);
    nconv = k;
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
    if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;

    if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
      ierr = BVCopyColumn(eps->V,n,nconv);CHKERRQ(ierr);
      if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
        /* Reorthonormalize restart vector */
        ierr = BVOrthogonalizeColumn(eps->V,nconv,NULL,&norm,&breakdown);CHKERRQ(ierr);
        ierr = BVScaleColumn(eps->V,nconv,1.0/norm);CHKERRQ(ierr);
      }
      if (breakdown) {
        /* Use random vector for restarting */
        ierr = PetscInfo(eps,"Using random vector for restart\n");CHKERRQ(ierr);
        ierr = EPSGetStartVector(eps,nconv,&breakdown);CHKERRQ(ierr);
      }
      if (breakdown) { /* give up */
        eps->reason = EPS_DIVERGED_BREAKDOWN;
        ierr = PetscInfo(eps,"Unable to generate more start vectors\n");CHKERRQ(ierr);
      }
    }
  }
  eps->nconv = nconv;

  ierr = PetscFree4(ritz,bnd,perm,conv);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}
Ejemplo n.º 5
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);
}
Ejemplo n.º 6
0
PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
{
  PetscErrorCode  ierr;
  EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
  PetscInt        i,j,*pj,k,l,nv,ld;
  Mat             U;
  PetscScalar     *S,*Q,*g;
  PetscReal       beta,gamma=1.0;
  PetscBool       breakdown,harmonic;

  PetscFunctionBegin;
  ierr = DSGetLeadingDimension(eps->ds,&ld);CHKERRQ(ierr);
  harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
  if (harmonic) { ierr = PetscMalloc1(ld,&g);CHKERRQ(ierr); }
  if (eps->arbitrary) pj = &j;
  else pj = NULL;

  /* Get the starting Arnoldi vector */
  ierr = EPSGetStartVector(eps,0,NULL);CHKERRQ(ierr);
  l = 0;

  /* 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 = DSGetArray(eps->ds,DS_MAT_A,&S);CHKERRQ(ierr);
    ierr = EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);CHKERRQ(ierr);
    ierr = DSRestoreArray(eps->ds,DS_MAT_A,&S);CHKERRQ(ierr);
    ierr = DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);CHKERRQ(ierr);
    if (l==0) {
      ierr = DSSetState(eps->ds,DS_STATE_INTERMEDIATE);CHKERRQ(ierr);
    } else {
      ierr = DSSetState(eps->ds,DS_STATE_RAW);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,g,&gamma);CHKERRQ(ierr);
    }

    /* Solve projected problem */
    ierr = DSSolve(eps->ds,eps->eigr,eps->eigi);CHKERRQ(ierr);
    if (eps->arbitrary) {
      ierr = EPSGetArbitraryValues(eps,eps->rr,eps->ri);CHKERRQ(ierr);
      j=1;
    }
    ierr = DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);CHKERRQ(ierr);

    /* Check convergence */
    ierr = EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);CHKERRQ(ierr);
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;

    /* Update l */
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
    else {
      l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
#if !defined(PETSC_USE_COMPLEX)
      ierr = DSGetArray(eps->ds,DS_MAT_A,&S);CHKERRQ(ierr);
      if (S[k+l+(k+l-1)*ld] != 0.0) {
        if (k+l<nv-1) l = l+1;
        else l = l-1;
      }
      ierr = DSRestoreArray(eps->ds,DS_MAT_A,&S);CHKERRQ(ierr);
#endif
    }

    if (eps->reason == EPS_CONVERGED_ITERATING) {
      if (breakdown) {
        /* Start a new Arnoldi factorization */
        ierr = PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);CHKERRQ(ierr);
        if (k<eps->nev) {
          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);
          }
        }
      } else {
        /* Undo translation of Krylov decomposition */
        if (harmonic) {
          ierr = DSSetDimensions(eps->ds,nv,0,k,l);CHKERRQ(ierr);
          ierr = DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);CHKERRQ(ierr);
          /* gamma u^ = u - U*g~ */
          ierr = BVMultColumn(eps->V,-1.0,1.0,nv,g);CHKERRQ(ierr);
          ierr = BVScaleColumn(eps->V,nv,1.0/gamma);CHKERRQ(ierr);
        }
        /* Prepare the Rayleigh quotient for restart */
        ierr = DSGetArray(eps->ds,DS_MAT_A,&S);CHKERRQ(ierr);
        ierr = DSGetArray(eps->ds,DS_MAT_Q,&Q);CHKERRQ(ierr);
        for (i=k;i<k+l;i++) {
          S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
        }
        ierr = DSRestoreArray(eps->ds,DS_MAT_A,&S);CHKERRQ(ierr);
        ierr = DSRestoreArray(eps->ds,DS_MAT_Q,&Q);CHKERRQ(ierr);
      }
    }
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
    ierr = DSGetMat(eps->ds,DS_MAT_Q,&U);CHKERRQ(ierr);
    ierr = BVMultInPlace(eps->V,U,eps->nconv,k+l);CHKERRQ(ierr);
    ierr = MatDestroy(&U);CHKERRQ(ierr);

    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
      ierr = BVCopyColumn(eps->V,nv,k+l);CHKERRQ(ierr);
    }
    eps->nconv = k;
    ierr = EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);CHKERRQ(ierr);
  }

  if (harmonic) { ierr = PetscFree(g);CHKERRQ(ierr); }
  /* 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);
}
Ejemplo n.º 7
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);
}