/*
I should not modify setup called!!
This is handled via petsc.
*/
PetscErrorCode BSSCR_PCSetUp_ScGtKG( PC pc )
{
	PC_SC_GtKG    ctx = (PC_SC_GtKG)pc->data;
	
	if( ctx->F == PETSC_NULL ) {   Stg_SETERRQ( PETSC_ERR_SUP, "gtkg: F not set" );   }
	if( ctx->Bt == PETSC_NULL ) {  Stg_SETERRQ( PETSC_ERR_SUP, "gtkg: Bt not set" );  }
	
	if(ctx->ksp_BBt==PETSC_NULL) {
		BSSCR_PCScGtKGUseStandardBBtOperator( pc ) ;
	}
	
	BSSCR_PCScGtKGBBtContainsConstantNullSpace( pc, &ctx->BBt_has_cnst_nullspace );
	if( ctx->BBt_has_cnst_nullspace == PETSC_TRUE ) {
		PetscPrintf( PETSC_COMM_WORLD, "\t* Detected prescence of constant nullspace in BBt-C\n" );
	}
	
	PetscFunctionReturn(0);
}
PetscErrorCode _Stokes_SLE_UzawaSolver_FormResidual( void *stokesSLE, void *solver, Vec r )
{
        Stokes_SLE_UzawaSolver *self = (Stokes_SLE_UzawaSolver*)solver;
        Stokes_SLE *sle = (Stokes_SLE*)stokesSLE;
 	/* stg linear algebra objects */
        Mat         A12 = sle->gStiffMat->matrix;
        Mat         A22 = NULL;
        Vec         x2 = sle->pSolnVec->vector;
        Vec         f_star = self->fTempVec;
        Vec         u_star = self->vStarVec;
	Vec         q_star = self->pTempVec;
        KSP     A11_solver = self->velSolver;

	/* petsc objects */
	KSP ksp_A11;
	PetscInt r_N, x2_N;

        /* check operations will be valid */
	if (A11_solver==NULL){   Stg_SETERRQ( PETSC_ERR_ARG_NULL, "vel_solver is NULL" ); }
        if (sle->dStiffMat!=NULL) {  Stg_SETERRQ( PETSC_ERR_SUP, "A21 must be NULL" ); }
        if (A12==NULL){          Stg_SETERRQ( PETSC_ERR_ARG_NULL, "A12 is NULL" ); }
        if (x2==NULL){         Stg_SETERRQ( PETSC_ERR_ARG_NULL, "x2 is NULL" ); }
        if (u_star==NULL){     Stg_SETERRQ( PETSC_ERR_ARG_NULL, "u* is NULL" ); }
        if (f_star==NULL){     Stg_SETERRQ( PETSC_ERR_ARG_NULL, "f* is NULL" ); }
	if (q_star==NULL) {    Stg_SETERRQ( PETSC_ERR_ARG_NULL, "q* is NULL" ); }

	A22 = PETSC_NULL;
	if (sle->cStiffMat!=NULL) {
		A22 = sle->cStiffMat->matrix;
	}


        /* Extract petsc objects */
        ksp_A11 =     A11_solver;

	/* Check sizes match */
	VecGetSize( r, &r_N );
	VecGetSize( x2, &x2_N );
	if (r_N!=x2_N) {
		Stg_SETERRQ2( PETSC_ERR_ARG_SIZ, "Solution vector for pressure (N=%D) is not compatible with residual vector (N=%D)", x2_N, r_N );
	}	
	
	/* r = f_hat - (G^T K^{-1} G - M) p */
	_Stokes_SLE_UzawaSolver_GetRhs( stokesSLE, solver, r );  /* r <- f_hat */

	/* correct for non zero A22 */
	if (A22!=PETSC_NULL) {
		MatMultAdd( A22, x2, r, r );  /* r <- r + A22 p */
	}

	/* make correction r <- r - G^T K^{-1} G */
	MatMult( A12, x2, f_star );                 /* f* <- A12 x2 */
	KSPSolve( ksp_A11, f_star, u_star );        /* u* <- A11^{-1} f* */
	MatMultTranspose( A12, u_star, q_star );    /* q* <- A12 u* */
	
	VecAXPY( r, -1.0, q_star );  /* r <- r - q* */
}
PetscErrorCode _Stokes_SLE_UzawaSolver_GetRhs( void *stokesSLE, void *solver, Vec rhs )
{
        Stokes_SLE_UzawaSolver *self = (Stokes_SLE_UzawaSolver*)solver;
	Stokes_SLE *sle = (Stokes_SLE*)stokesSLE;
	/* stg linear algebra */
	KSP     A11_solver = self->velSolver;
        Mat         A12 = sle->gStiffMat->matrix;
        Vec         b1  = sle->fForceVec->vector;
        Vec         b2  = sle->hForceVec->vector;
        Vec         u_star = self->vStarVec;
	
	/* petsc variables */
	KSP ksp_A11;

	/* check operations will be valid */
	if (sle->dStiffMat!=NULL) {   Stg_SETERRQ( PETSC_ERR_SUP, "A21 must be NULL" ); }
	if (A11_solver==NULL){    Stg_SETERRQ( PETSC_ERR_ARG_NULL, "vel_solver is NULL" ); }
        if (A12==NULL){           Stg_SETERRQ( PETSC_ERR_ARG_NULL, "A12 is NULL" ); }
        if (b1==NULL){       Stg_SETERRQ( PETSC_ERR_ARG_NULL, "b1 is NULL" ); }
        if (b2==NULL){       Stg_SETERRQ( PETSC_ERR_ARG_NULL, "b2 is NULL" ); }
        if (u_star==NULL){   Stg_SETERRQ( PETSC_ERR_ARG_NULL, "u* is NULL" ); }

	/* Extract petsc objects */
	ksp_A11 = A11_solver;

	/* compute rhs = A12^T A11^{-1} b1 - b2 */
	KSPSolve( ksp_A11, b1, u_star );         /* u* = A11^{-1} b1 */
	MatMultTranspose( A12, u_star, rhs );    /* b2 = A12^T u* */
	VecAXPY( rhs, -1.0, b2 );  /* rhs <- rhs - b2 */
}
/*
F & Bt must different to PETSC_NULL
B may be PETSC_NULL
C can be PETSC_NULL
*/
PetscErrorCode BSSCR_PCScGtKGSetOperators( PC pc, Mat F, Mat Bt, Mat B, Mat C )
{
	PC_SC_GtKG ctx = (PC_SC_GtKG)pc->data;
	BSSCR_BSSCR_pc_error_ScGtKG( pc, __func__ );
	
	
	ctx->F  = F;
	ctx->Bt = Bt;
	ctx->B  = B;
	ctx->C  = C;
	
	if( C != PETSC_NULL ) {
		pc->ops->apply          = BSSCR_BSSCR_PCApply_ScGtKG_C;
		pc->ops->applytranspose = PETSC_NULL;
	}
	
	/* Create vectors */
	if( ctx->s == PETSC_NULL ) {  MatGetVecs( ctx->F, &ctx->s, PETSC_NULL );  }
	if( ctx->X == PETSC_NULL ) {  MatGetVecs( ctx->F, PETSC_NULL, &ctx->X );  }
	if( ctx->t == PETSC_NULL ) {  MatGetVecs( ctx->Bt, &ctx->t, PETSC_NULL ); }
	
	if( ctx->F == PETSC_NULL ) {   Stg_SETERRQ( PETSC_ERR_SUP, "gtkg: F not set" );   }
	if( ctx->Bt == PETSC_NULL ) {  Stg_SETERRQ( PETSC_ERR_SUP, "gtkg: Bt not set" );  }
	
	if( ctx->X1 == PETSC_NULL ) {  MatGetVecs( ctx->F, &ctx->X1, PETSC_NULL );  }
	if( ctx->Y1 == PETSC_NULL ) {  MatGetVecs( ctx->F, &ctx->Y1, PETSC_NULL );  }
	if( ctx->X2 == PETSC_NULL ) {  MatGetVecs( ctx->Bt, &ctx->X2, PETSC_NULL ); }
	if( ctx->Y2 == PETSC_NULL ) {  MatGetVecs( ctx->Bt, &ctx->Y2, PETSC_NULL ); }
	
	VecSet( ctx->X1, 1.0 );
	VecSet( ctx->Y1, 1.0 );
	VecSet( ctx->X2, 1.0 );
	VecSet( ctx->Y2, 1.0 );
	
	
	PetscFunctionReturn(0);
}
// updated
PetscErrorCode BSSCR_MatMVBlock_ConstructScaling( MatStokesBlockScaling BA, Mat A, Vec b, Vec x, PetscTruth sym )
{
	
	
	if( BA->scaling_exists == PETSC_FALSE ) {
		PetscInt M,N;
		PetscTruth is_block;
		
		/* check A is 2x2 block matrix */
		Stg_PetscObjectTypeCompare( (PetscObject)A, "block", &is_block );
		if (is_block==PETSC_FALSE) {
			Stg_SETERRQ( PETSC_ERR_SUP, "Only valid for MatType = block" );
		}
		MatGetSize( A, &M, &N );
		if ( (M!=2) || (N!=2) ) {
			Stg_SETERRQ2( PETSC_ERR_SUP, "Only valid for 2x2 block. Yours has dimension %Dx%D", M,N );
		}
		
		
		
		VecDuplicate( x, &BA->Lz ); 
		VecDuplicate( x, &BA->Rz );
		
		BA->scaling_exists = PETSC_TRUE;
	}
	
//	if( BA->user_build_scaling != PETSC_NULL ) {
//		BA->user_build_scaling( A,b,x,BA->scaling_ctx);
//	}
//	else {
	BSSCR_MatStokesMVBlockDefaultBuildScaling(BA,A,b,x, sym);
//	}
	
	BA->scalings_have_been_inverted = PETSC_FALSE;
	
	PetscFunctionReturn(0);
}
PetscErrorCode BSSCR_KSPNormInfConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
{
    PetscErrorCode         ierr;
    KSPPWConvergedCtx     *cctx = (KSPPWConvergedCtx*)ctx;
    KSPNormType            normtype;
    PetscReal              min, max, R_max, R_min, R_Ninf;
    Vec                    R, work, w1,w2;
        
    PetscFunctionBegin;
    PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
    PetscValidPointer(reason,4);
    *reason = KSP_CONVERGED_ITERATING;
                        
    ierr = VecDuplicate(ksp->vec_rhs,&work);CHKERRQ(ierr);
    ierr = VecDuplicate(ksp->vec_rhs,&w1);CHKERRQ(ierr);
    ierr = VecDuplicate(ksp->vec_rhs,&w2);CHKERRQ(ierr);      
                
    KSPBuildResidual( ksp, w1,w2, &R );
    VecNorm( R, NORM_INFINITY, &R_Ninf );
                
    //PetscPrintf( PETSC_COMM_WORLD, "Norm inf convergence %s\n ", ksp->prefix);
        
    cctx->pointwise_max = R_Ninf;
        
    ierr = KSPGetNormType(ksp,&normtype); 
        
    CHKERRQ(ierr);
    if (normtype == KSP_NORM_NO) 
        Stg_SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Use BSSCR_KSPSkipConverged() with KSPNormType of KSP_NORM_NO");
        
    if (!cctx) 
        Stg_SETERRQ(PETSC_ERR_ARG_NULL,"Convergence context must have been created with BSSCR_KSPDefaultConvergedCreate()");
        
    if (!n) {
        /* if user gives initial guess need to compute norm of b */
        if (!ksp->guess_zero && !cctx->initialrtol) {
            PetscReal      snorm;
            if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
                ierr = PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of RHS\n");
                CHKERRQ(ierr);
                ierr = VecNorm(ksp->vec_rhs,NORM_INFINITY,&snorm);CHKERRQ(ierr);        /*     <- b'*b */       
                                
                PetscPrintf( PETSC_COMM_WORLD,  "Non Zero Guess; RHS - %g\n", snorm);
                                                
            } 
            else {
                Vec z;
                if (!cctx->work) {
                    ierr = VecDuplicate(ksp->vec_rhs,&cctx->work);CHKERRQ(ierr);
                }
                z = cctx->work;
                ierr = KSP_PCApply(ksp,ksp->vec_rhs,z);CHKERRQ(ierr);
                if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
                    ierr = PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n");CHKERRQ(ierr);
                    ierr = VecNorm(z,NORM_INFINITY,&snorm);CHKERRQ(ierr);                 /*    dp <- b'*B'*B*b */

                } 
                else if (ksp->normtype == KSP_NORM_NATURAL) {
                    PetscScalar norm;
                    Vec bz;
                    ierr = PetscInfo(ksp,"user has provided nonzero initial guess, computing natural norm of RHS\n");CHKERRQ(ierr);
                    //        ierr  = VecDot(ksp->vec_rhs,z,&norm);
                    //        snorm = sqrt(PetscAbsScalar(norm));                            /*    dp <- b'*B*b */
                    VecDuplicate( z, &bz );
                    VecPointwiseMult( bz, ksp->vec_rhs, z );
                    ierr = VecNorm(bz,NORM_INFINITY,&snorm);CHKERRQ(ierr);                 
                    Stg_VecDestroy(&bz);
                }
            }
            /* handle special case of zero RHS and nonzero guess */
            if (!snorm) {
                ierr = PetscInfo(ksp,"Special case, user has provided nonzero initial guess and zero RHS\n");CHKERRQ(ierr);
                snorm = rnorm;
            }
            if (cctx->mininitialrtol) {
                ksp->rnorm0 = PetscMin(snorm,rnorm);
            } else {
                ksp->rnorm0 = snorm;
            }
        } else {
            ksp->rnorm0 = rnorm;
        }
        ksp->ttol   = PetscMax(ksp->rtol*ksp->rnorm0,ksp->abstol);
    }
        
    // if (n <= ksp->chknorm) PetscFunctionReturn(0);
        
    if ( R_Ninf != R_Ninf ) {
        ierr = PetscInfo(ksp,"Linear solver has created a not a number (NaN) as the pointwise residual norm, declaring divergence \n");CHKERRQ(ierr);
        *reason = KSP_DIVERGED_NAN;
    } 
    else 
        if (R_Ninf <= ksp->ttol) {
            if (R_Ninf < ksp->abstol) {
                ierr = PetscInfo3(ksp,"Linear solver has converged. Pointwise residual %G is less than absolute tolerance %G at iteration %D\n",R_Ninf,ksp->abstol,n);
                CHKERRQ(ierr);
                *reason = KSP_CONVERGED_ATOL;
            } 
            else {
                if (cctx->initialrtol) {
                    ierr = PetscInfo4(ksp,"Linear solver has converged. Norm_infinity %G is less than relative tolerance %G times initial Norm_infinity %G at iteration %D\n",R_Ninf,ksp->rtol,ksp->rnorm0,n);
                    CHKERRQ(ierr);
                } 
                else {
                    ierr = PetscInfo4(ksp,"Linear solver has converged. Norm_infinity %G is less than relative tolerance %G times initial norm_infinity right hand side %G at iteration %D\n",R_Ninf,ksp->rtol,ksp->rnorm0,n);CHKERRQ(ierr);
                }
                *reason = KSP_CONVERGED_RTOL;
            }
        } 
        else 
            if (R_Ninf >= ksp->divtol*ksp->rnorm0) {
                ierr = PetscInfo3(ksp,"Linear solver is diverging. Initial right hand size Norm_infinity value %G, current residual norm %G at iteration %D\n",ksp->rnorm0,R_Ninf,n);CHKERRQ(ierr);
                *reason = KSP_DIVERGED_DTOL;
            }
                        
    /* trash all work vectors here */

    Stg_VecDestroy(&work);
    Stg_VecDestroy(&w1);
    Stg_VecDestroy(&w2);
        
    PetscFunctionReturn(0);
}
PetscErrorCode BSSCR_MatStokesMVBlockReportOperatorScales( Mat A, PetscTruth sym )
{
	Vec rA, rG;
	PetscInt loc,M,N;
	PetscReal min, max;
	Mat K,G,D,C;
	PetscTruth is_block;
	
	
	/* check A is 2x2 block matrix */
	Stg_PetscObjectTypeCompare( (PetscObject)A, "block", &is_block );
	if (is_block==PETSC_FALSE) {
		Stg_SETERRQ( PETSC_ERR_SUP, "Only valid for MatType = block" );
	}
	MatGetSize( A, &M, &N );
	if ( (M!=2) || (N!=2) ) {
		Stg_SETERRQ2( PETSC_ERR_SUP, "Only valid for 2x2 block. Yours has dimension %Dx%D", M,N );
	}
	
	
	MatNestGetSubMat( A, 0,0, &K );
	MatNestGetSubMat( A, 0,1, &G );
	MatNestGetSubMat( A, 1,0, &D );
	MatNestGetSubMat( A, 1,1, &C );
	
	
	MatGetVecs( K, PETSC_NULL, &rA );
	VecDuplicate( rA, &rG );
	
	/* Report the row max and mins */
	if (K!=PETSC_NULL) {
		MatGetRowMax( K, rA, PETSC_NULL );
		VecMax( rA, &loc, &max );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_max(K) = %g \n", max );
		
		MatGetRowMinAbs( K, rA, PETSC_NULL );
		VecMin( rA, &loc, &min );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_min(K) = %g \n\n", min );
	}
	
	if( G != PETSC_NULL ) {       
		MatGetRowMax( G, rG, PETSC_NULL );
		VecMax( rG, &loc, &max );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_max(G) = %g \n", max );
		
		MatGetRowMinAbs( G, rG, PETSC_NULL );
		VecMin( rG, &loc, &min );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_min(G) = %g \n", min );
	}
	
	if( D != PETSC_NULL && !sym ) {
                Vec rD;

                MatGetVecs( D, PETSC_NULL, &rD );
		MatGetRowMax( D, rD, PETSC_NULL );
		VecMax( rD, &loc, &max );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_max(D) = %g \n", max );
		
		MatGetRowMinAbs( D, rD, PETSC_NULL );
		VecMin( rD, &loc, &min );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_min(D) = %g \n", min );

                Stg_VecDestroy(&rD );
	}
	
	if( C != PETSC_NULL ) {
		Vec cG;

		MatGetVecs( G, &cG, PETSC_NULL );
		MatGetRowMax( C, cG, PETSC_NULL );
		VecMax( cG, &loc, &max );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_max(C) = %g \n", max );
		
		MatGetRowMin( C, cG, PETSC_NULL );
		VecMin( cG, &loc, &min );
		PetscPrintf( PETSC_COMM_WORLD, "Sup_min(C) = %g \n\n", min );
	
		Stg_VecDestroy(&cG);
	}
	
	
	Stg_VecDestroy(&rA );
	Stg_VecDestroy(&rG );
	
	
	PetscFunctionReturn(0);
}
예제 #8
0
/*
I should not modify setup called!!
This is handled via petsc.
*/
PetscErrorCode BSSCR_PCSetUp_GtKG( PC pc )
{
	PC_GtKG    ctx = (PC_GtKG)pc->data;
	PetscReal  fill;
	Mat        Ident;
	Vec        diag;
	PetscInt   M,N, m,n;
	MPI_Comm   comm;
	PetscInt   nnz_I, nnz_G;
	const MatType mtype;
	const char *prefix;
	PetscTruth wasSetup;
	
	
	
	if( ctx->K == PETSC_NULL ) {
		Stg_SETERRQ( PETSC_ERR_SUP, "gtkg: K not set" );
	}
	if( ctx->G == PETSC_NULL ) {
		Stg_SETERRQ( PETSC_ERR_SUP, "gtkg: G not set" );
	}
	
	PetscObjectGetComm( (PetscObject)ctx->K, &comm ); 
	
	
	/* Check for existence of objects and trash any which exist */
	if( ctx->form_GtG == PETSC_TRUE && ctx->GtG != PETSC_NULL ) {
		Stg_MatDestroy(&ctx->GtG );
		ctx->GtG = PETSC_NULL;
	}
	
	if( ctx->s != PETSC_NULL ) {
		Stg_VecDestroy(&ctx->s );
		ctx->s = PETSC_NULL;
	}
	if( ctx->X != PETSC_NULL ) {
		Stg_VecDestroy(&ctx->X );
		ctx->X = PETSC_NULL;
	}
	if( ctx->t != PETSC_NULL ) {
		Stg_VecDestroy(&ctx->t );
		ctx->t = PETSC_NULL;
	}
	if( ctx->inv_diag_M != PETSC_NULL ) {
		Stg_VecDestroy(&ctx->inv_diag_M );
		ctx->inv_diag_M = PETSC_NULL;
	}
	
	
	
	/* Create vectors */
	MatGetVecs( ctx->K, &ctx->s, &ctx->X );
	MatGetVecs( ctx->G, &ctx->t, PETSC_NULL );
	
	if( ctx->M != PETSC_NULL ) {
		MatGetVecs( ctx->K, &ctx->inv_diag_M, PETSC_NULL );
		MatGetDiagonal( ctx->M, ctx->inv_diag_M );
		VecReciprocal( ctx->inv_diag_M );
		
		/* change the pc_apply routines */
		pc->ops->apply          = BSSCR_BSSCR_PCApply_GtKG_diagonal_scaling;
		pc->ops->applytranspose = BSSCR_BSSCR_PCApplyTranspose_GtKG_diagonal_scaling;
	}
	
	
	/* Assemble GtG */
	MatGetSize( ctx->G, &M, &N );
	MatGetLocalSize( ctx->G, &m, &n );
	
	MatGetVecs( ctx->G, PETSC_NULL, &diag );
	VecSet( diag, 1.0 );
	
	MatCreate( comm, &Ident );
	MatSetSizes( Ident, m,m , M, M );
#if (((PETSC_VERSION_MAJOR==3) && (PETSC_VERSION_MINOR>=3)) || (PETSC_VERSION_MAJOR>3) )
        MatSetUp(Ident);
#endif

	MatGetType( ctx->G, &mtype );
	MatSetType( Ident, mtype );
	
	if( ctx->M == PETSC_NULL ) {
		MatDiagonalSet( Ident, diag, INSERT_VALUES );
	}
	else {
		MatDiagonalSet( Ident, ctx->inv_diag_M, INSERT_VALUES );
	}
	
	BSSCR_get_number_nonzeros_AIJ( Ident, &nnz_I );
	BSSCR_get_number_nonzeros_AIJ( ctx->G, &nnz_G );
	//fill = 1.0;
	/* 
	Not sure the best way to estimate the fill factor.
	GtG is a laplacian on the pressure space. 
	This might tell us something useful...
	*/
	fill = (PetscReal)(nnz_G)/(PetscReal)( nnz_I );
	MatPtAP( Ident, ctx->G, MAT_INITIAL_MATRIX, fill, &ctx->GtG );
	
	Stg_MatDestroy(&Ident);
	Stg_VecDestroy(&diag );
	
	
	Stg_KSPSetOperators( ctx->ksp, ctx->GtG, ctx->GtG, SAME_NONZERO_PATTERN );
	
	if (!pc->setupcalled) {	
		wasSetup = PETSC_FALSE;
		
		PCGetOptionsPrefix( pc,&prefix );
		KSPSetOptionsPrefix( ctx->ksp, prefix );
		KSPAppendOptionsPrefix( ctx->ksp, "pc_gtkg_" ); /* -pc_GtKG_ksp_type <type>, -ksp_GtKG_pc_type <type> */
	}
	else {
		wasSetup = PETSC_TRUE;
	}	
	
	
//	if (!wasSetup && pc->setfromoptionscalled) {
	if (!wasSetup) {
		KSPSetFromOptions(ctx->ksp);
	}
	
	
	PetscFunctionReturn(0);
}