static PetscErrorCode VecDotNorm2_Nest(Vec x,Vec y,PetscScalar *dp, PetscScalar *nm) { Vec_Nest *bx = (Vec_Nest*)x->data; Vec_Nest *by = (Vec_Nest*)y->data; PetscInt i,nr; PetscScalar x_dot_y,_dp,_nm; PetscReal norm2_y; PetscErrorCode ierr; PetscFunctionBegin; nr = bx->nb; _dp = 0.0; _nm = 0.0; for (i=0; i<nr; i++) { ierr = VecDotNorm2(bx->v[i],by->v[i],&x_dot_y,&norm2_y);CHKERRQ(ierr); _dp += x_dot_y; _nm += norm2_y; } *dp = _dp; *nm = _nm; PetscFunctionReturn(0); }
PetscErrorCode KSPSolve_GCR_cycle( KSP ksp ) { KSP_GCR *ctx = (KSP_GCR*)ksp->data; PetscErrorCode ierr; PetscScalar r_dot_v; Mat A, B; PC pc; Vec s,v,r; PetscReal norm_r,nrm; PetscInt k, i, restart; Vec x; PetscReal res; PetscFunctionBegin; restart = ctx->restart; ierr = KSPGetPC( ksp, &pc );CHKERRQ(ierr); ierr = KSPGetOperators( ksp, &A, &B, 0 );CHKERRQ(ierr); x = ksp->vec_sol; r = ctx->R; for ( k=0; k<restart; k++ ) { v = ctx->VV[k]; s = ctx->SS[k]; if (ctx->modifypc) { ierr = (*ctx->modifypc)(ksp,ksp->its,ksp->rnorm,ctx->modifypc_ctx);CHKERRQ(ierr); } ierr = PCApply( pc, r, s );CHKERRQ(ierr); /* s = B^{-1} r */ ierr = MatMult( A, s, v );CHKERRQ(ierr); /* v = A s */ ierr = VecMDot( v,k, ctx->VV, ctx->val );CHKERRQ(ierr); for (i=0; i<k; i++) ctx->val[i] = -ctx->val[i]; ierr = VecMAXPY(v,k,ctx->val,ctx->VV);CHKERRQ(ierr); /* v = v - sum_{i=0}^{k-1} alpha_i v_i */ ierr = VecMAXPY(s,k,ctx->val,ctx->SS);CHKERRQ(ierr); /* s = s - sum_{i=0}^{k-1} alpha_i s_i */ ierr = VecDotNorm2(r,v,&r_dot_v,&nrm);CHKERRQ(ierr); nrm = PetscSqrtReal(nrm); r_dot_v = r_dot_v/nrm; ierr = VecScale( v, 1.0/nrm );CHKERRQ(ierr); ierr = VecScale( s, 1.0/nrm );CHKERRQ(ierr); ierr = VecAXPY( x, r_dot_v, s );CHKERRQ(ierr); ierr = VecAXPY( r, -r_dot_v, v );CHKERRQ(ierr); if (ksp->its > ksp->chknorm ) { ierr = VecNorm( r, NORM_2, &norm_r );CHKERRQ(ierr); } /* update the local counter and the global counter */ ksp->its++; res = norm_r; ksp->rnorm = res; KSPLogResidualHistory(ksp,res); ierr = KSPMonitor(ksp,ksp->its,res);CHKERRQ(ierr); if ( ksp->its > ksp->chknorm ) { ierr = (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; } if ( ksp->its >= ksp->max_it ) { ksp->reason = KSP_CONVERGED_ITS; break; } } ctx->n_restarts++; PetscFunctionReturn(0); }
PetscErrorCode test_axpy_dot_max( void ) { Vec x1,y1, x2,y2; Vec tmp_buf[2]; Vec X, Y; PetscReal real,real2; PetscScalar scalar; PetscInt index; PetscErrorCode ierr; PetscFunctionBegin; PetscPrintf( PETSC_COMM_WORLD, "\n\n============== %s ==============\n", PETSC_FUNCTION_NAME ); gen_test_vector( PETSC_COMM_WORLD, 4, 0, 1, &x1 ); gen_test_vector( PETSC_COMM_WORLD, 5, 10, 2, &x2 ); gen_test_vector( PETSC_COMM_WORLD, 4, 4, 3, &y1 ); gen_test_vector( PETSC_COMM_WORLD, 5, 5, 1, &y2 ); tmp_buf[0] = x1; tmp_buf[1] = x2; ierr = VecCreateNest(PETSC_COMM_WORLD,2,PETSC_NULL,tmp_buf,&X);CHKERRQ(ierr); ierr = VecAssemblyBegin(X);CHKERRQ(ierr); ierr = VecAssemblyEnd(X);CHKERRQ(ierr); ierr = VecDestroy(&x1);CHKERRQ(ierr); ierr = VecDestroy(&x2);CHKERRQ(ierr); tmp_buf[0] = y1; tmp_buf[1] = y2; ierr = VecCreateNest(PETSC_COMM_WORLD,2,PETSC_NULL,tmp_buf,&Y);CHKERRQ(ierr); ierr = VecAssemblyBegin(Y);CHKERRQ(ierr); ierr = VecAssemblyEnd(Y);CHKERRQ(ierr); ierr = VecDestroy(&y1);CHKERRQ(ierr); ierr = VecDestroy(&y2);CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "VecAXPY \n"); ierr = VecAXPY( Y, 1.0, X ); /* Y <- a X + Y */ ierr = VecNestGetSubVec( Y, 0, &y1 );CHKERRQ(ierr); ierr = VecNestGetSubVec( Y, 1, &y2 );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "(1) y1 = \n" ); ierr = VecView( y1, PETSC_VIEWER_STDOUT_WORLD );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "(1) y2 = \n" ); ierr = VecView( y2, PETSC_VIEWER_STDOUT_WORLD );CHKERRQ(ierr); ierr = VecDot( X,Y, &scalar );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "X.Y = %lf + %lfi \n", PetscRealPart(scalar), PetscImaginaryPart(scalar) ); ierr = VecDotNorm2( X,Y, &scalar, &real2 );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "X.Y = %lf + %lfi norm2(Y) = %lf\n", PetscRealPart(scalar), PetscImaginaryPart(scalar), real2); ierr = VecAXPY( Y, 1.0, X ); /* Y <- a X + Y */ ierr = VecNestGetSubVec( Y, 0, &y1 );CHKERRQ(ierr); ierr = VecNestGetSubVec( Y, 1, &y2 );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "(2) y1 = \n" ); ierr = VecView( y1, PETSC_VIEWER_STDOUT_WORLD );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "(2) y2 = \n" ); ierr = VecView( y2, PETSC_VIEWER_STDOUT_WORLD );CHKERRQ(ierr); ierr = VecDot( X,Y, &scalar );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "X.Y = %lf + %lfi \n", PetscRealPart(scalar), PetscImaginaryPart(scalar) ); ierr = VecDotNorm2( X,Y, &scalar, &real2 );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "X.Y = %lf + %lfi norm2(Y) = %lf\n", PetscRealPart(scalar), PetscImaginaryPart(scalar), real2); ierr = VecMax( X, &index, &real );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "(max-X) = %f : index = %d \n", real, index ); ierr = VecMin( X, &index, &real );CHKERRQ(ierr); PetscPrintf( PETSC_COMM_WORLD, "(min-X) = %f : index = %d \n", real, index ); ierr = VecDestroy(&X);CHKERRQ(ierr); ierr = VecDestroy(&Y);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode KSPSolve_Richardson(KSP ksp) { PetscErrorCode ierr; PetscInt i,maxit; PetscReal rnorm = 0.0,abr; PetscScalar scale,rdot; Vec x,b,r,z,w = NULL,y = NULL; PetscInt xs, ws; Mat Amat,Pmat; KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data; PetscBool exists,diagonalscale; PetscFunctionBegin; ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr); if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name); ksp->its = 0; ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); x = ksp->vec_sol; b = ksp->vec_rhs; ierr = VecGetSize(x,&xs);CHKERRQ(ierr); ierr = VecGetSize(ksp->work[0],&ws);CHKERRQ(ierr); if (xs != ws) { if (richardsonP->selfscale) { ierr = KSPSetWorkVecs(ksp,4);CHKERRQ(ierr); } else { ierr = KSPSetWorkVecs(ksp,2);CHKERRQ(ierr); } } r = ksp->work[0]; z = ksp->work[1]; if (richardsonP->selfscale) { w = ksp->work[2]; y = ksp->work[3]; } maxit = ksp->max_it; /* if user has provided fast Richardson code use that */ ierr = PCApplyRichardsonExists(ksp->pc,&exists);CHKERRQ(ierr); if (exists && !ksp->numbermonitors && !ksp->transpose_solve & !ksp->nullsp) { PCRichardsonConvergedReason reason; ierr = PCApplyRichardson(ksp->pc,b,x,r,ksp->rtol,ksp->abstol,ksp->divtol,maxit,ksp->guess_zero,&ksp->its,&reason);CHKERRQ(ierr); ksp->reason = (KSPConvergedReason)reason; PetscFunctionReturn(0); } scale = richardsonP->scale; if (!ksp->guess_zero) { /* r <- b - A x */ ierr = KSP_MatMult(ksp,Amat,x,r);CHKERRQ(ierr); ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); } else { ierr = VecCopy(b,r);CHKERRQ(ierr); } ksp->its = 0; if (richardsonP->selfscale) { ierr = KSP_PCApply(ksp,r,z);CHKERRQ(ierr); /* z <- B r */ for (i=0; i<maxit; i++) { if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) { ierr = VecNorm(r,NORM_2,&rnorm);CHKERRQ(ierr); /* rnorm <- r'*r */ ierr = KSPMonitor(ksp,i,rnorm);CHKERRQ(ierr); ksp->rnorm = rnorm; ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) { ierr = VecNorm(z,NORM_2,&rnorm);CHKERRQ(ierr); /* rnorm <- z'*z */ ierr = KSPMonitor(ksp,i,rnorm);CHKERRQ(ierr); ksp->rnorm = rnorm; ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; } ierr = KSP_PCApplyBAorAB(ksp,z,y,w);CHKERRQ(ierr); /* y = BAz = BABr */ ierr = VecDotNorm2(z,y,&rdot,&abr);CHKERRQ(ierr); /* rdot = (Br)^T(BABR); abr = (BABr)^T (BABr) */ scale = rdot/abr; ierr = PetscInfo1(ksp,"Self-scale factor %g\n",(double)PetscRealPart(scale));CHKERRQ(ierr); ierr = VecAXPY(x,scale,z);CHKERRQ(ierr); /* x <- x + scale z */ ierr = VecAXPY(r,-scale,w);CHKERRQ(ierr); /* r <- r - scale*Az */ ierr = VecAXPY(z,-scale,y);CHKERRQ(ierr); /* z <- z - scale*y */ ksp->its++; } } else { for (i=0; i<maxit; i++) { if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) { ierr = VecNorm(r,NORM_2,&rnorm);CHKERRQ(ierr); /* rnorm <- r'*r */ ierr = KSPMonitor(ksp,i,rnorm);CHKERRQ(ierr); ksp->rnorm = rnorm; ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; } ierr = KSP_PCApply(ksp,r,z);CHKERRQ(ierr); /* z <- B r */ if (ksp->normtype == KSP_NORM_PRECONDITIONED) { ierr = VecNorm(z,NORM_2,&rnorm);CHKERRQ(ierr); /* rnorm <- z'*z */ ierr = KSPMonitor(ksp,i,rnorm);CHKERRQ(ierr); ksp->rnorm = rnorm; ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; } ierr = VecAXPY(x,scale,z);CHKERRQ(ierr); /* x <- x + scale z */ ksp->its++; if (i+1 < maxit || ksp->normtype != KSP_NORM_NONE) { ierr = KSP_MatMult(ksp,Amat,x,r);CHKERRQ(ierr); /* r <- b - Ax */ ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); } } } if (!ksp->reason) { if (ksp->normtype != KSP_NORM_NONE) { if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) { ierr = VecNorm(r,NORM_2,&rnorm);CHKERRQ(ierr); /* rnorm <- r'*r */ } else { ierr = KSP_PCApply(ksp,r,z);CHKERRQ(ierr); /* z <- B r */ ierr = VecNorm(z,NORM_2,&rnorm);CHKERRQ(ierr); /* rnorm <- z'*z */ } ksp->rnorm = rnorm; ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = KSPMonitor(ksp,i,rnorm);CHKERRQ(ierr); } if (ksp->its >= ksp->max_it) { if (ksp->normtype != KSP_NORM_NONE) { ierr = (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS; } else { ksp->reason = KSP_CONVERGED_ITS; } } } PetscFunctionReturn(0); }