PetscErrorCode bsscr_DGMiGtD( Mat *_K2, Mat K, Mat G, Mat M){ Mat K2; Vec diag; Vec Mdiag; Mat MinvGt; Mat Gtrans; PetscErrorCode ierr; PetscFunctionBegin; MatGetVecs( K, &diag, PETSC_NULL ); MatGetDiagonal( K, diag ); VecSqrt(diag); //VecReciprocal(diag);/* trying something different here */ MatGetVecs( M, &Mdiag, PETSC_NULL ); MatGetDiagonal( M, Mdiag ); VecReciprocal(Mdiag); #if( PETSC_VERSION_MAJOR <= 2 ) ierr=MatTranspose(G, &Gtrans);CHKERRQ(ierr); #else ierr=MatTranspose(G, MAT_INITIAL_MATRIX,&Gtrans);CHKERRQ(ierr); #endif ierr=MatConvert(Gtrans, MATSAME, MAT_INITIAL_MATRIX, &MinvGt);CHKERRQ(ierr);/* copy Gtrans -> MinvGt */ MatDiagonalScale(MinvGt, Mdiag, PETSC_NULL);/* Minv*Gtrans */ /* MAT_INITIAL_MATRIX -> creates K2 matrix : PETSC_DEFAULT for fill ratio: run with -info to find what it should be*/ ierr=MatMatMult( G, MinvGt, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &K2);CHKERRQ(ierr);/* K2 = G*Minv*Gtrans */ MatDiagonalScale(K2, diag, diag );/* K2 = D*K2*D = D*G*Minv*Gtrans*D */ Stg_MatDestroy(&Gtrans); Stg_MatDestroy(&MinvGt); Stg_VecDestroy(&diag); Stg_VecDestroy(&Mdiag); *_K2=K2; PetscFunctionReturn(0); }
// updated PetscErrorCode BSSCR_MatStokesMVBlockDefaultBuildScaling( MatStokesBlockScaling BA, Mat A, Vec b, Vec x, PetscTruth sym ) { Mat K,G,D,C; Vec rG; PetscScalar rg2, rg, ra; PetscInt N; Vec rA, rC; Vec L1,L2, R1,R2; Mat S; VecNestGetSubVec( BA->Lz, 0, &L1 ); VecNestGetSubVec( BA->Lz, 1, &L2 ); VecNestGetSubVec( BA->Rz, 0, &R1 ); VecNestGetSubVec( BA->Rz, 1, &R2 ); rA = L1; rC = L2; MatNestGetSubMat( A, 0,0, &K ); MatNestGetSubMat( A, 0,1, &G ); MatNestGetSubMat( A, 1,0, &D ); MatNestGetSubMat( A, 1,1, &C ); VecDuplicate( rA, &rG ); /* Get magnitude of K */ //px_MatGetAbsRowSum( K, rA ); //MatGetRowMax( K, rA, PETSC_NULL ); MatGetDiagonal( K, rA); VecSqrt( rA ); VecReciprocal( rA ); /* VecDot( rA,rA, &ra ); */ /* VecGetSize( rA, &N ); */ /* ra = PetscSqrtScalar( ra/N ); */ /* Get magnitude of G */ //px_MatGetAbsRowSum( G, rG ); //MatGetRowMax( G, rG, PETSC_NULL ); Mat A21_cpy; Mat Shat; Vec diag; /* same as rA*rA */ MatGetVecs( K, &diag, PETSC_NULL ); MatGetDiagonal( K, diag ); VecReciprocal( diag ); if( sym ) #if( PETSC_VERSION_MAJOR <= 2 ) MatTranspose( G, &A21_cpy ); #else MatTranspose( G, MAT_INITIAL_MATRIX, &A21_cpy ); #endif else {
PetscErrorCode ComputeShearStress(UserContext *uc) { PetscErrorCode ierr; Vec u2, v2; PetscReal scale = VISCOSITY / DELTA_Z; // mu * (u^2 + v^2)^0.5 / dz PetscFunctionBegin; PetscLogEventBegin(EVENT_ComputeShearStress,0,0,0,0); VecDuplicate(uc->v, &uc->ss); VecDuplicate(uc->v, &u2); VecDuplicate(uc->v, &v2); VecPointwiseMult( u2, uc->u, uc->u); VecPointwiseMult( v2, uc->v, uc->v); VecWAXPY(uc->ss, one, u2, v2); VecSqrt(uc->ss); VecScale(uc->ss, scale); VecDestroy(u2); VecDestroy(v2); PetscLogEventEnd(EVENT_ComputeShearStress,0,0,0,0); PetscFunctionReturn(0); }
PetscErrorCode BSSCR_PCScGtKGUseStandardScaling( PC pc ) { PC_SC_GtKG ctx = (PC_SC_GtKG)pc->data; Mat K,G,D,C; Vec rG; PetscScalar rg2, rg, ra; PetscInt N; Vec rA, rC; Vec L1,L2, R1,R2; BSSCR_BSSCR_pc_error_ScGtKG( pc, __func__ ); L1 = ctx->X1; L2 = ctx->X2; R1 = ctx->Y1; R2 = ctx->Y2; rA = L1; rC = L2; K = ctx->F; G = ctx->Bt; D = ctx->B; C = ctx->C; VecDuplicate( rA, &rG ); /* Get magnitude of K */ MatGetRowMax( K, rA, PETSC_NULL ); VecSqrt( rA ); VecReciprocal( rA ); VecDot( rA,rA, &ra ); VecGetSize( rA, &N ); ra = PetscSqrtScalar( ra/N ); /* Get magnitude of G */ MatGetRowMax( G, rG, PETSC_NULL ); VecDot( rG, rG, &rg2 ); VecGetSize( rG, &N ); rg = PetscSqrtScalar(rg2/N); // printf("rg = %f \n", rg ); VecSet( rC, 1.0/(rg*ra) ); Stg_VecDestroy(&rG ); VecCopy( L1, R1 ); VecCopy( L2, R2 ); PetscFunctionReturn(0); }
PetscErrorCode PressurePoissonUpdate( ) { PetscErrorCode ierr; PetscFunctionBegin; PetscLogEventBegin(EVENT_PressurePoissonUpdate,0,0,0,0); // PetscLogEventRegister(&EVENT_PressurePoissonUpdate,"PressurePoissonUpdate", 0); VecSet(rhsC->v, 0.); // TODO: this is important when reusing vecs, need to reset back to zero after time step ierr = IIMComputeCorrection2D(iim,JumpConditionPressure,ls,rhsC); CHKERRQ(ierr); IIMDiscreteCompatabilityCondition(ls, rhsC); ierr = KSPSolve(ksp, rhsC->v, p->v); CHKERRQ(ierr); IrregularNodeListWrite(ls,t); ierr = IIMPressureGradient(ls,p,px,py); CHKERRQ(ierr); ierr = IIMComputeCorrection2D(iim,JumpConditionXVelocity,ls,px); CHKERRQ(ierr); px->v2[d1/2][d2/2] = 0.; ierr = KSPSolve(ksp, px->v, u->v); CHKERRQ(ierr); ierr = IIMComputeCorrection2D(iim,JumpConditionYVelocity,ls,py); CHKERRQ(ierr); py->v2[d1/2][d2/2] = 0.; ierr = KSPSolve(ksp, py->v, v->v); CHKERRQ(ierr); VecPointwiseMult(u2, u->v, u->v); VecPointwiseMult(v2, u->v, u->v); VecWAXPY(magVel, 1, u2, v2); VecSqrt( magVel ); VecNorm(magVel, NORM_INFINITY, &maxVel); dt = PetscMin(1/(2*maxVel),.1); LevelSetAdvect2D(dt, u, v, ls, lstemp); printf("\tmaxVel: %f",maxVel); printf("\tdt: %f\n",dt); PetscLogEventEnd(EVENT_PressurePoissonUpdate,0,0,0,0); PetscFunctionReturn(0); }
END_TEST START_TEST( test_IIMContext ) { IIM iim; int d1=63, d2=d1; Grid2D rhsC, p, px, py, u, v; LevelSet2D ls, lstemp; PetscErrorCode ierr; Mat m; KSP ksp; PC pc; Generate2DLaplacianPeriodicBC( d1, d2, &m); ierr = KSPCreate(PETSC_COMM_SELF, &ksp); CHKERRQ(ierr); KSPGetPC(ksp, &pc); // KSPSetType(ksp, KSPCG); // PCSetType(pc, PCICC); // PCFactorSetLevels(pc, 0); // KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); KSPSetType(ksp, KSPPREONLY); PCSetType(pc, PCCHOLESKY); PCFactorSetMatOrderingType(pc, MATORDERING_ND); KSPSetFromOptions(ksp); KSPSetOperators(ksp, m, m, SAME_PRECONDITIONER); PCSetUp(pc); ierr = IIMCreate(&iim, 12); CHKERRQ(ierr); // IIMSetForceComponents(iim, ForceComponentNormalSimple, ForceComponentTangentialSimple); CreateGrid2D(d1, d2, &rhsC); CreateGrid2D(d1, d2, &p); CreateGrid2D(d1, d2, &px); CreateGrid2D(d1, d2, &py); CreateGrid2D(d1, d2, &u); CreateGrid2D(d1, d2, &v); ierr = CreateLevelSet2D(d1,d2,&ls); CHKERRQ(ierr); CreateLevelSet2D(d1, d2, &lstemp); LevelSetInitializeToStar(ls,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE); // LevelSetInitializeToCircle(ls, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE); Vec u2,v2,magVel; VecCreateSeq(PETSC_COMM_SELF,d1*d2,&magVel); VecDuplicate(magVel, &u2); VecDuplicate(magVel, &v2); PetscReal maxVel, dt; WriteVec wv, wvC, wvP, wvPX, wvPY, wvU, wvV; WriteVecCreate(lstemp->g2d->v, "ls",&wv); WriteVecCreate(rhsC->v, "c", &wvC); WriteVecCreate(p->v, "p", &wvP); WriteVecCreate(px->v, "px", &wvPY); WriteVecCreate(py->v, "py", &wvPX); WriteVecCreate(u->v, "u", &wvU); WriteVecCreate(v->v, "v", &wvV); for( int t = 0; t < 100; t++ ) { printf("TIME: %d ",t); VecSet(rhsC->v, 0.); // TODO: this is important when reusing vecs, need to reset back to zero after time step ierr = IIMComputeCorrection2D(iim,JumpConditionPressure,ls,rhsC); CHKERRQ(ierr); IIMDiscreteCompatabilityCondition(ls, rhsC); ierr = KSPSolve(ksp, rhsC->v, p->v); CHKERRQ(ierr); IrregularNodeListWrite(ls,t); ierr = IIMPressureGradient(ls,p,px,py); CHKERRQ(ierr); ierr = IIMComputeCorrection2D(iim,JumpConditionXVelocity,ls,px); CHKERRQ(ierr); px->v2[d1/2][d2/2] = 0.; ierr = KSPSolve(ksp, px->v, u->v); CHKERRQ(ierr); ierr = IIMComputeCorrection2D(iim,JumpConditionYVelocity,ls,py); CHKERRQ(ierr); py->v2[d1/2][d2/2] = 0.; ierr = KSPSolve(ksp, py->v, v->v); CHKERRQ(ierr); VecPointwiseMult(u2, u->v, u->v); VecPointwiseMult(v2, u->v, u->v); VecWAXPY(magVel, 1, u2, v2); VecSqrt( magVel ); VecNorm(magVel, NORM_INFINITY, &maxVel); dt = PetscMin(1/(2*maxVel),.1); LevelSetAdvect2D(dt, u, v, ls, lstemp); printf("\tmaxVel: %f",maxVel); printf("\tdt: %f\n",dt); WriteVecToDisk(wv); WriteVecToDisk(wvC); WriteVecToDisk(wvP); WriteVecToDisk(wvPX); WriteVecToDisk(wvPY); WriteVecToDisk(wvU); WriteVecToDisk(wvV); VecCopy( lstemp->g2d->v, ls->g2d->v); UpdateIrregularNodeList(ls); // ReinitializeLevelSet(ls); } ierr = VecDestroy(magVel); CHKERRQ(ierr); ierr = VecDestroy(u2); CHKERRQ(ierr); ierr = VecDestroy(v2); CHKERRQ(ierr); ierr = DestroyLevelSet2D(ls); CHKERRQ(ierr); ierr = DestroyGrid2D(rhsC); CHKERRQ(ierr); ierr = DestroyGrid2D(p); CHKERRQ(ierr); ierr = DestroyGrid2D(px); CHKERRQ(ierr); ierr = DestroyGrid2D(py); CHKERRQ(ierr); ierr = DestroyGrid2D(u); CHKERRQ(ierr); ierr = DestroyGrid2D(v); CHKERRQ(ierr); ierr = IIMDestroy(iim); CHKERRQ(ierr); }
PetscErrorCode bsscr_buildK2(KSP ksp){ KSP_BSSCR *bsscr = (KSP_BSSCR *)ksp->data; Mat K,G,M, K2=0; Vec f2=0; //MatStokesBlockScaling BA=bsscr->BA; Stokes_SLE* stokesSLE = (Stokes_SLE*)bsscr->st_sle; //PetscErrorCode ierr; PetscFunctionBegin; K = stokesSLE->kStiffMat->matrix; G = stokesSLE->gStiffMat->matrix; if(bsscr->K2){ Stg_MatDestroy(&bsscr->K2); bsscr->K2 = PETSC_NULL; } switch (bsscr->k2type) { case (K2_DGMGD): {/* Should only do this one if scaling is turned off. */ Vec D; AugLagStokes_SLE * stokesSLE = (AugLagStokes_SLE*)bsscr->st_sle; if (stokesSLE->mStiffMat){ MatGetVecs( K, PETSC_NULL, &D ); //MatGetRowMax( K, D, PETSC_NULL ); MatGetDiagonal( K, D ); VecSqrt( D ); M = stokesSLE->mStiffMat->matrix; bsscr_GMiGt(&K2,K,G,M); /* K2 created */ MatDiagonalScale(K2, D, D ); /* K2 = D*G*inv(M)*Gt*D */ Stg_VecDestroy(&D); PetscPrintf( PETSC_COMM_WORLD, "\n\n----- K2_DGMGD ------\n\n"); } else{ PetscPrintf( PETSC_COMM_WORLD,"The Augmented Lagrangian Method DGMGD was specified but the SLE has no mStiffMat on it.\n"); PetscPrintf( PETSC_COMM_WORLD,"You need to use the AugLagStokes_SLE class.\n"); } } break; case (K2_GMG): { AugLagStokes_SLE * stokesSLE = (AugLagStokes_SLE*)bsscr->st_sle; if (stokesSLE->mStiffMat){ M = stokesSLE->mStiffMat->matrix; bsscr_GMiGt(&K2,K,G,M); /* K2 created */ PetscPrintf( PETSC_COMM_WORLD, "\n\n----- K2_GMG ------\n\n"); }else{ PetscPrintf( PETSC_COMM_WORLD,"The Augmented Lagrangian Method GMG was specified but the SLE has no mStiffMat on it.\n"); PetscPrintf( PETSC_COMM_WORLD,"You need to use the AugLagStokes_SLE class.\n"); } } break; case (K2_GG): { bsscr_GGt(&K2,K,G); /* K2 created */ PetscPrintf( PETSC_COMM_WORLD, "\n\n----- K2_GG ------\n\n"); } break; case (K2_SLE): { AugLagStokes_SLE * stokesSLE = (AugLagStokes_SLE*)bsscr->st_sle; if (stokesSLE->mStiffMat){ K2 = stokesSLE->mStiffMat->matrix; if(stokesSLE->jForceVec){ f2 = stokesSLE->jForceVec->vector; } PetscPrintf( PETSC_COMM_WORLD, "\n\n----- K2_SLE ------\n\n"); }else{ PetscPrintf( PETSC_COMM_WORLD,"The Augmented Lagrangian Method SLE was specified but the SLE has no Matrix on it.\n"); PetscPrintf( PETSC_COMM_WORLD,"You need to use the AugLagStokes_SLE class.\n"); } } break; default: PetscPrintf( PETSC_COMM_WORLD, "\n\n----- NO K2 ------\n\n"); } bsscr->f2 = f2; bsscr->K2 = K2; PetscFunctionReturn(0); }
static PetscErrorCode KSPSolve_LSQR(KSP ksp) { PetscErrorCode ierr; PetscInt i,size1,size2; PetscScalar rho,rhobar,phi,phibar,theta,c,s,tmp,tau,alphac; PetscReal beta,alpha,rnorm; Vec X,B,V,V1,U,U1,TMP,W,W2,SE,Z = PETSC_NULL; Mat Amat,Pmat; MatStructure pflag; KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data; PetscTruth diagonalscale,nopreconditioner; PetscFunctionBegin; ierr = PCDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr); if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name); ierr = PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);CHKERRQ(ierr); ierr = PetscTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);CHKERRQ(ierr); /* nopreconditioner =PETSC_FALSE; */ /* Calculate norm of right hand side */ ierr = VecNorm(ksp->vec_rhs,NORM_2,&lsqr->rhs_norm);CHKERRQ(ierr); /* Calculate norm of matrix*/ ierr = MatNorm( Amat, NORM_FROBENIUS, &lsqr->anorm);CHKERRQ(ierr); /* vectors of length m, where system size is mxn */ B = ksp->vec_rhs; U = lsqr->vwork_m[0]; U1 = lsqr->vwork_m[1]; /* vectors of length n */ X = ksp->vec_sol; W = lsqr->vwork_n[0]; V = lsqr->vwork_n[1]; V1 = lsqr->vwork_n[2]; W2 = lsqr->vwork_n[3]; if (!nopreconditioner) { Z = lsqr->vwork_n[4]; } /* standard error vector */ SE = lsqr->se; if (SE){ ierr = VecGetSize(SE,&size1);CHKERRQ(ierr); ierr = VecGetSize(X ,&size2);CHKERRQ(ierr); if (size1 != size2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Standard error vector (size %d) does not match solution vector (size %d)",size1,size2); ierr = VecSet(SE,0.0);CHKERRQ(ierr); } /* Compute initial residual, temporarily use work vector u */ if (!ksp->guess_zero) { ierr = KSP_MatMult(ksp,Amat,X,U);CHKERRQ(ierr); /* u <- b - Ax */ ierr = VecAYPX(U,-1.0,B);CHKERRQ(ierr); } else { ierr = VecCopy(B,U);CHKERRQ(ierr); /* u <- b (x is 0) */ } /* Test for nothing to do */ ierr = VecNorm(U,NORM_2,&rnorm);CHKERRQ(ierr); ierr = PetscObjectTakeAccess(ksp);CHKERRQ(ierr); ksp->its = 0; ksp->rnorm = rnorm; ierr = PetscObjectGrantAccess(ksp);CHKERRQ(ierr); KSPLogResidualHistory(ksp,rnorm); KSPMonitor(ksp,0,rnorm); ierr = (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) PetscFunctionReturn(0); ierr = VecCopy(B,U);CHKERRQ(ierr); ierr = VecNorm(U,NORM_2,&beta);CHKERRQ(ierr); ierr = VecScale(U,1.0/beta);CHKERRQ(ierr); ierr = KSP_MatMultTranspose(ksp,Amat,U,V); CHKERRQ(ierr); if (nopreconditioner) { ierr = VecNorm(V,NORM_2,&alpha); CHKERRQ(ierr); } else { ierr = PCApply(ksp->pc,V,Z);CHKERRQ(ierr); ierr = VecDot(V,Z,&alphac);CHKERRQ(ierr); if (PetscRealPart(alphac) <= 0.0) { ksp->reason = KSP_DIVERGED_BREAKDOWN; PetscFunctionReturn(0); } alpha = sqrt(PetscRealPart(alphac)); ierr = VecScale(Z,1.0/alpha); CHKERRQ(ierr); } ierr = VecScale(V,1.0/alpha);CHKERRQ(ierr); if (nopreconditioner){ ierr = VecCopy(V,W);CHKERRQ(ierr); } else { ierr = VecCopy(Z,W);CHKERRQ(ierr); } ierr = VecSet(X,0.0);CHKERRQ(ierr); lsqr->arnorm = alpha * beta; phibar = beta; rhobar = alpha; tau = -beta; i = 0; do { if (nopreconditioner) { ierr = KSP_MatMult(ksp,Amat,V,U1);CHKERRQ(ierr); } else { ierr = KSP_MatMult(ksp,Amat,Z,U1);CHKERRQ(ierr); } ierr = VecAXPY(U1,-alpha,U);CHKERRQ(ierr); ierr = VecNorm(U1,NORM_2,&beta);CHKERRQ(ierr); if (beta == 0.0){ ksp->reason = KSP_DIVERGED_BREAKDOWN; break; } ierr = VecScale(U1,1.0/beta);CHKERRQ(ierr); ierr = KSP_MatMultTranspose(ksp,Amat,U1,V1);CHKERRQ(ierr); ierr = VecAXPY(V1,-beta,V);CHKERRQ(ierr); if (nopreconditioner) { ierr = VecNorm(V1,NORM_2,&alpha); CHKERRQ(ierr); } else { ierr = PCApply(ksp->pc,V1,Z);CHKERRQ(ierr); ierr = VecDot(V1,Z,&alphac);CHKERRQ(ierr); if (PetscRealPart(alphac) <= 0.0) { ksp->reason = KSP_DIVERGED_BREAKDOWN; break; } alpha = sqrt(PetscRealPart(alphac)); ierr = VecScale(Z,1.0/alpha);CHKERRQ(ierr); } ierr = VecScale(V1,1.0/alpha);CHKERRQ(ierr); rho = PetscSqrtScalar(rhobar*rhobar + beta*beta); c = rhobar / rho; s = beta / rho; theta = s * alpha; rhobar = - c * alpha; phi = c * phibar; phibar = s * phibar; tau = s * phi; ierr = VecAXPY(X,phi/rho,W);CHKERRQ(ierr); /* x <- x + (phi/rho) w */ if (SE) { ierr = VecCopy(W,W2);CHKERRQ(ierr); ierr = VecSquare(W2);CHKERRQ(ierr); ierr = VecScale(W2,1.0/(rho*rho));CHKERRQ(ierr); ierr = VecAXPY(SE, 1.0, W2);CHKERRQ(ierr); /* SE <- SE + (w^2/rho^2) */ } if (nopreconditioner) { ierr = VecAYPX(W,-theta/rho,V1);CHKERRQ(ierr); /* w <- v - (theta/rho) w */ } else { ierr = VecAYPX(W,-theta/rho,Z);CHKERRQ(ierr); /* w <- z - (theta/rho) w */ } lsqr->arnorm = alpha*PetscAbsScalar(tau); rnorm = PetscRealPart(phibar); ierr = PetscObjectTakeAccess(ksp);CHKERRQ(ierr); ksp->its++; ksp->rnorm = rnorm; ierr = PetscObjectGrantAccess(ksp);CHKERRQ(ierr); KSPLogResidualHistory(ksp,rnorm); KSPMonitor(ksp,i+1,rnorm); ierr = (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; SWAP(U1,U,TMP); SWAP(V1,V,TMP); i++; } while (i<ksp->max_it); if (i >= ksp->max_it && !ksp->reason) { ksp->reason = KSP_DIVERGED_ITS; } /* Finish off the standard error estimates */ if (SE) { tmp = 1.0; ierr = MatGetSize(Amat,&size1,&size2);CHKERRQ(ierr); if ( size1 > size2 ) tmp = size1 - size2; tmp = rnorm / PetscSqrtScalar(tmp); ierr = VecSqrt(SE);CHKERRQ(ierr); ierr = VecScale(SE,tmp);CHKERRQ(ierr); } PetscFunctionReturn(0); }