PetscErrorCode KSPBuildPressure_CB_Nullspace_BSSCR(KSP ksp) { KSP_BSSCR *bsscr = (KSP_BSSCR *)ksp->data; FeEquationNumber *eq_num = bsscr->solver->st_sle->pSolnVec->feVariable->eqNum; FeMesh *feMesh = bsscr->solver->st_sle->pSolnVec->feVariable->feMesh; /* is the pressure mesh */ unsigned ijk[3]; Vec t, v; int numLocalNodes, globalNodeNumber, i, j, eq; MatStokesBlockScaling BA = bsscr->BA; PetscErrorCode ierr; Mat Amat,Pmat, G; MatStructure pflag; PetscFunctionBegin; /* get G matrix from Amat matrix operator on ksp */ ierr=Stg_PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);CHKERRQ(ierr); MatNestGetSubMat( Amat, 0,1, &G );/* G should always exist */ /* now create Vecs t and v to match size of G: i.e. pressure */ /* NOTE: not using "h" vector from ksp->vec_rhs because this part of the block vector doesn't always exist */ MatGetVecs( G, &t, PETSC_NULL );/* t and v are destroyed in KSPDestroy_BSSCR */ MatGetVecs( G, &v, PETSC_NULL );/* t and v such that can do G*t */ numLocalNodes = Mesh_GetLocalSize( feMesh, MT_VERTEX); /* number of nodes on current proc not counting any shadow nodes */ for(j=0;j<numLocalNodes;j++){ i = globalNodeNumber = Mesh_DomainToGlobal( feMesh, MT_VERTEX, j); RegularMeshUtils_Element_1DTo3D(feMesh, i, ijk); eq = eq_num->destinationArray[j][0];/* get global equation number -- 2nd arg is always 0 because pressure has only one dof */ if(eq != -1){ if( (ijk[0]+ijk[1]+ijk[2])%2 ==0 ){ VecSetValue(t,eq,1.0,INSERT_VALUES); } else{ VecSetValue(v,eq,1.0,INSERT_VALUES); }}} VecAssemblyBegin( t ); VecAssemblyEnd( t ); VecAssemblyBegin( v ); VecAssemblyEnd( v ); /* Scaling the null vectors here because it easier at the moment *//* maybe should do this in the original scaling function */ if( BA->scaling_exists == PETSC_TRUE ){ Vec R2; /* Get the scalings out the block mat data */ VecNestGetSubVec( BA->Rz, 1, &R2 ); VecPointwiseDivide( t, t, R2); /* x <- x * 1/R2 */ VecPointwiseDivide( v, v, R2); } bsscr_writeVec( t, "t", "Writing t vector"); bsscr_writeVec( v, "v", "Writing v vector"); bsscr->t=t; bsscr->v=v; PetscFunctionReturn(0); }
PetscErrorCode KSPRemovePressureNullspace_BSSCR(KSP ksp, Vec h_hat) { KSP_BSSCR *bsscr = (KSP_BSSCR *)ksp->data; MatStokesBlockScaling BA = bsscr->BA; PetscErrorCode ierr; PetscScalar norm, a, a1, a2, hnorm, pnorm, gnorm; Vec t2, t,v; Mat Amat,Pmat, D; MatStructure pflag; double nstol = bsscr->nstol; /* set in KSPCreate_BSSCR */ PetscFunctionBegin; t=bsscr->t; v=bsscr->v; /* get G matrix from Amat matrix operator on ksp */ ierr=Stg_PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);CHKERRQ(ierr); //MatNestGetSubMat( Amat, 0,1, &G );/* G should always exist */ MatNestGetSubMat( Amat, 1,0, &D );/* D should always exist */ /* now create Vec t2 to match left hand size of G: i.e. velocity */ MatGetVecs( D, &t2, PETSC_NULL ); MatNorm(D,NORM_INFINITY,&gnorm); /* seems like not a bad estimate of the largest eigenvalue for this matrix */ if(t){/* assumes that v and t are initially set to PETSC_NULL (in KSPCreate_BSSCR ) */ MatMultTranspose( D, t, t2); VecNorm(t2, NORM_2, &norm); VecNorm(t, NORM_2, &a); norm=norm/a;/* so we are using unit null vector */ if(norm < nstol*gnorm){/* then t in NS of G */ VecDot(t,h_hat, &a1); VecDot(t,t, &a2); a=-a1/a2; VecAXPY(h_hat, a, t); VecDot(t,h_hat, &a1); PetscPrintf( PETSC_COMM_WORLD, "\n\t* Found t NS norm= %g : Dot = %g\n\n", norm, a1); } } if(v){ MatMultTranspose( D, v, t2); VecNorm(t2, NORM_2, &norm); VecNorm(v, NORM_2, &a); norm=norm/a;/* so we are using unit null vector */ if(norm < nstol*gnorm){/* then v in NS of G */ VecDot(v,h_hat, &a1); VecDot(v,v, &a2); a=-a1/a2; VecAXPY(h_hat, a, v); VecDot(v,h_hat, &a1); PetscPrintf( PETSC_COMM_WORLD, "\n\t* Found v NS norm= %g : Dot = %g\n\n", norm, a1); } } bsscr_writeVec( h_hat, "h", "Writing h_hat Vector in Solver"); Stg_VecDestroy(&t2); PetscFunctionReturn(0); }
PetscErrorCode KSPBuildPressure_Const_Nullspace_BSSCR(KSP ksp) { KSP_BSSCR *bsscr = (KSP_BSSCR *)ksp->data; FeEquationNumber *eq_num = bsscr->solver->st_sle->pSolnVec->feVariable->eqNum; FeMesh *feMesh = bsscr->solver->st_sle->pSolnVec->feVariable->feMesh; /* is the pressure mesh */ unsigned ijk[3]; Vec t; int numLocalNodes, globalNodeNumber, i, eq; MatStokesBlockScaling BA = bsscr->BA; PetscErrorCode ierr; Mat Amat,Pmat, G; MatStructure pflag; PetscFunctionBegin; /* get G matrix from Amat matrix operator on ksp */ ierr=Stg_PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);CHKERRQ(ierr); MatNestGetSubMat( Amat, 0,1, &G );/* G should always exist */ /* now create Vec t to match size of G: i.e. pressure */ /* NOTE: not using "h" vector from ksp->vec_rhs because this part of the block vector doesn't always exist */ MatGetVecs( G, &t, PETSC_NULL );/* t is destroyed in KSPDestroy_BSSCR */ VecSet(t, 1.0); VecAssemblyBegin( t ); VecAssemblyEnd( t ); /* Scaling the null vectors here because it easier at the moment */ if( BA->scaling_exists == PETSC_TRUE ){ Vec R2; /* Get the scalings out the block mat data */ VecNestGetSubVec( BA->Rz, 1, &R2 ); VecPointwiseDivide( t, t, R2); /* x <- x * 1/R2 */ } bsscr_writeVec( t, "t", "Writing t vector"); bsscr->t=t; PetscFunctionReturn(0); }
PetscErrorCode KSPSolve_BSSCR(KSP ksp) { Mat Amat,Pmat; /* Stokes Matrix and it's Preconditioner matrix: Both 2x2 PetscExt block matrices */ Vec B, X; /* rhs and solution vectors */ //MatStructure pflag; PetscErrorCode ierr; KSP_BSSCR * bsscr; Stokes_SLE * SLE; //PETScMGSolver * MG; Mat K,D,ApproxS; MatStokesBlockScaling BA; PetscTruth flg, sym, augment; double TotalSolveTime; PetscFunctionBegin; TotalSolveTime = MPI_Wtime(); PetscPrintf( PETSC_COMM_WORLD, "\nBSSCR -- Block Stokes Schur Compliment Reduction Solver \n"); /** Get the stokes Block matrix and its preconditioner matrix */ ierr = Stg_PCGetOperators(ksp->pc,&Amat,&Pmat,PETSC_NULL);CHKERRQ(ierr); /** In Petsc proper, KSP's ksp->data is usually set in KSPCreate_XXX function. Here it is set in the _StokesBlockKSPInterface_Solve function instead so that we can ensure that the solver has everything it needs */ bsscr = (KSP_BSSCR*)ksp->data; //MG = (PETScMGSolver*)bsscr->mg; SLE = (Stokes_SLE*)bsscr->st_sle; X = ksp->vec_sol; B = ksp->vec_rhs; if( bsscr->do_scaling ){ (*bsscr->scale)(ksp); /* scales everything including the UW preconditioner */ BA = bsscr->BA; } if( (bsscr->k2type != 0) ){ if(bsscr->buildK2 != PETSC_NULL) { (*bsscr->buildK2)(ksp); /* building K2 from scaled version of stokes operators: K2 lives on bsscr struct = ksp->data */ } } /* get sub matrix / vector objects */ MatNestGetSubMat( Amat, 0,0, &K ); /* Underworld preconditioner matrix*/ ApproxS = PETSC_NULL; if( ((StokesBlockKSPInterface*)SLE->solver)->preconditioner ) { /* SLE->solver->st_sle == SLE here, by the way */ StiffnessMatrix *preconditioner; preconditioner = ((StokesBlockKSPInterface*)SLE->solver)->preconditioner; ApproxS = BSSCR_GetPetscMatrix( preconditioner->matrix ); } sym = bsscr->DIsSym; MatNestGetSubMat( Amat, 1,0, &D );if(!D){ PetscPrintf( PETSC_COMM_WORLD, "D does not exist but should!!\n"); exit(1); } /**********************************************************/ /******* SOLVE!! ******************************************/ /**********************************************************/ flg = PETSC_FALSE; augment = PETSC_TRUE; PetscOptionsGetTruth(PETSC_NULL, "-augmented_lagrangian", &augment, &flg); BSSCR_DRIVER_auglag( ksp, Amat, X, B, ApproxS, BA, sym, bsscr ); /**********************************************************/ /***** END SOLVE!! ****************************************/ /**********************************************************/ if( bsscr->do_scaling ){ (*bsscr->unscale)(ksp); } if( (bsscr->k2type != 0) && bsscr->K2 != PETSC_NULL ){ if(bsscr->k2type != K2_SLE){/* don't destroy here, as in this case, K2 is just pointing to an existing matrix on the SLE */ Stg_MatDestroy(&bsscr->K2 ); } bsscr->K2built = PETSC_FALSE; bsscr->K2 = PETSC_NULL; } ksp->reason = KSP_CONVERGED_RTOL; TotalSolveTime = MPI_Wtime() - TotalSolveTime; PetscPrintf( PETSC_COMM_WORLD, " Total BSSCR Linear solve time: %lf seconds\n\n", TotalSolveTime); bsscr->solver->stats.total_time=TotalSolveTime; PetscFunctionReturn(0); }
PetscErrorCode BSSCR_DRIVER_auglag( KSP ksp, Mat stokes_A, Vec stokes_x, Vec stokes_b, Mat approxS, MatStokesBlockScaling BA, PetscTruth sym, KSP_BSSCR * bsscrp_self ) { AugLagStokes_SLE * stokesSLE = (AugLagStokes_SLE*)bsscrp_self->solver->st_sle; PetscTruth uzawastyle, KisJustK=PETSC_TRUE, restorek, change_A11rhspresolve; PetscTruth usePreviousGuess, useNormInfStoppingConditions, useNormInfMonitor, found, forcecorrection; PetscTruth change_backsolve, mg_active; PetscErrorCode ierr; //PetscInt monitor_index; PetscInt max_it,min_it; KSP ksp_inner, ksp_S, ksp_new_inner; PC pc_S, pcInner; Mat K,G,D,C, S, K2;// Korig; Vec u,p,f,f2=0,f3=0,h, h_hat,t; MGContext mgCtx; double mgSetupTime, scrSolveTime, a11SingleSolveTime, penaltyNumber;// hFactor; static int been_here = 0; /* Ha Ha Ha !! */ char name[PETSC_MAX_PATH_LEN]; char matname[PETSC_MAX_PATH_LEN]; char suffix[PETSC_MAX_PATH_LEN]; char str[PETSC_MAX_PATH_LEN]; PetscTruth flg, extractMats; PetscLogDouble flopsA,flopsB; /***************************************************************************************************************/ /***************************************************************************************************************/ //if( bsscrp_self->solver->st_sle->context->loadFromCheckPoint ){ // been_here=1; //} /* get sub matrix / vector objects */ /* note that here, the matrix D should always exist. It is set up in _StokesBlockKSPInterface_Solve in StokesBlockKSPInterface.c */ /* now extract K,G etc from a MatNest object */ MatNestGetSubMat( stokes_A, 0,0, &K ); MatNestGetSubMat( stokes_A, 0,1, &G ); MatNestGetSubMat( stokes_A, 1,0, &D );if(!D){ PetscPrintf( PETSC_COMM_WORLD, "D does not exist but should!!\n"); exit(1); } MatNestGetSubMat( stokes_A, 1,1, &C ); VecNestGetSubVec( stokes_x, 0, &u ); VecNestGetSubVec( stokes_x, 1, &p ); VecNestGetSubVec( stokes_b, 0, &f ); VecNestGetSubVec( stokes_b, 1, &h ); PetscPrintf( PETSC_COMM_WORLD, "\n\n---------- AUGMENTED LAGRANGIAN K2 METHOD ---------\n\n" ); PetscPrintf( PETSC_COMM_WORLD, "----------- Penalty = %f\n\n", bsscrp_self->solver->penaltyNumber ); sprintf(suffix,"%s","x"); PetscOptionsGetString( PETSC_NULL, "-matsuffix", suffix, PETSC_MAX_PATH_LEN-1, &extractMats ); flg=0; PetscOptionsGetString( PETSC_NULL, "-matdumpdir", name, PETSC_MAX_PATH_LEN-1, &flg ); if(flg){ sprintf(str,"%s/",name); sprintf(matname,"K%s",suffix); bsscr_dirwriteMat( K, matname,str, "Writing K matrix in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"G%s",suffix); bsscr_dirwriteMat( G, matname,str, "Writing G matrix in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"D%s",suffix); bsscr_dirwriteMat( D, matname,str, "Writing D matrix in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"f%s",suffix); bsscr_dirwriteVec( f, matname,str, "Writing f vector in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"h%s",suffix); bsscr_dirwriteVec( h, matname,str, "Writing h vector in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"Shat%s",suffix); bsscr_dirwriteMat( approxS, matname,str, "Writing Shat matrix in al Solver"); if(C){ sprintf(str,"%s/",name); sprintf(matname,"C%s",suffix); bsscr_dirwriteMat( C, matname,str, "Writing C matrix in al Solver"); } } mg_active=PETSC_TRUE; PetscOptionsGetTruth( PETSC_NULL ,"-A11_mg_active", &mg_active, &found ); bsscrp_self->solver->mg_active=mg_active; penaltyNumber = bsscrp_self->solver->penaltyNumber; //hFactor = stokesSLE->hFactor; /***************************************************************************************************************/ /***************************************************************************************************************/ /****** GET K2 ****************************************************************************************/ if(penaltyNumber > 1e-10 && bsscrp_self->k2type){ flg=0; PetscOptionsGetString( PETSC_NULL, "-matdumpdir", name, PETSC_MAX_PATH_LEN-1, &flg ); if(flg){ sprintf(str,"%s/",name); sprintf(matname,"K2%s",suffix); bsscr_dirwriteMat( bsscrp_self->K2, matname,str, "Writing K2 matrix in al Solver"); } K2=bsscrp_self->K2; scrSolveTime = MPI_Wtime(); ierr=MatAXPY(K,penaltyNumber,K2,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);/* Computes K = penaltyNumber*K2 + K */ scrSolveTime = MPI_Wtime() - scrSolveTime; PetscPrintf( PETSC_COMM_WORLD, "\n\t* K+p*K2 in time: %lf seconds\n\n", scrSolveTime); KisJustK=PETSC_FALSE; forcecorrection=PETSC_TRUE; PetscOptionsGetTruth( PETSC_NULL ,"-force_correction", &forcecorrection, &found ); if(forcecorrection){ if(bsscrp_self->f2 && forcecorrection){ f2=bsscrp_self->f2; ierr=VecAXPY(f,penaltyNumber,f2);/* f <- f +a*f2 */ }else{ switch (bsscrp_self->k2type) { case (K2_GG): { VecDuplicate( u, &f3 ); MatMult( G, h, f3); ierr=VecAXPY(f,penaltyNumber,f3);/* f <- f +a*f2 */ } break; case (K2_GMG): { Mat M; Vec Mdiag; VecDuplicate( u, &f3 ); M = bsscrp_self->solver->mStiffMat->matrix; MatGetVecs( M, &Mdiag, PETSC_NULL ); MatGetDiagonal( M, Mdiag ); VecReciprocal(Mdiag); VecPointwiseMult(Mdiag,Mdiag,h); MatMult( G, Mdiag, f3); ierr=VecAXPY(f,penaltyNumber,f3);/* f <- f +a*f2 */ Stg_VecDestroy(&Mdiag); } break; case (K2_DGMGD): { Mat M; Vec Mdiag; VecDuplicate( u, &f3 ); M = bsscrp_self->solver->mStiffMat->matrix; MatGetVecs( M, &Mdiag, PETSC_NULL ); MatGetDiagonal( M, Mdiag ); VecReciprocal(Mdiag); VecPointwiseMult(Mdiag,Mdiag,h); MatMult( G, Mdiag, f3); ierr=VecAXPY(f,penaltyNumber,f3);/* f <- f +a*f2 */ Stg_VecDestroy(&Mdiag); } break; case (K2_NULL): { ; } break; case (K2_SLE): { ; } break; } } } } /* Create Schur complement matrix */ MatCreateSchurComplement(K,K,G,D,C, &S); MatSchurComplementGetKSP( S, &ksp_inner); KSPGetPC( ksp_inner, &pcInner ); /***************************************************************************************************************/ /***************************************************************************************************************/ /********* SET PREFIX FOR INNER/VELOCITY KSP *************************************************************/ KSPSetOptionsPrefix( ksp_inner, "A11_" ); KSPSetFromOptions( ksp_inner ); Stg_KSPSetOperators(ksp_inner, K, K, DIFFERENT_NONZERO_PATTERN); useNormInfStoppingConditions = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL ,"-A11_use_norm_inf_stopping_condition", &useNormInfStoppingConditions, &found ); if(useNormInfStoppingConditions) BSSCR_KSPSetNormInfConvergenceTest( ksp_inner ); useNormInfMonitor = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-A11_ksp_norm_inf_monitor", &useNormInfMonitor, &found ); if(useNormInfMonitor) KSPMonitorSet( ksp_inner, BSSCR_KSPNormInfMonitor, PETSC_NULL, PETSC_NULL ); useNormInfMonitor = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-A11_ksp_norm_inf_to_norm_2_monitor", &useNormInfMonitor, &found ); if(useNormInfMonitor) KSPMonitorSet( ksp_inner, BSSCR_KSPNormInfToNorm2Monitor, PETSC_NULL, PETSC_NULL ); /***************************************************************************************************************/ /***************************************************************************************************************/ /* If multigrid is enabled, set it now. */ change_A11rhspresolve = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-change_A11rhspresolve", &change_A11rhspresolve, &found ); if(bsscrp_self->solver->mg_active && !change_A11rhspresolve) { mgSetupTime=setupMG( bsscrp_self, ksp_inner, pcInner, K, &mgCtx ); } /***************************************************************************************************************/ /***************************************************************************************************************/ /* create right hand side */ if(change_A11rhspresolve){ //Stg_KSPDestroy(&ksp_inner ); KSPCreate(PETSC_COMM_WORLD, &ksp_new_inner); Stg_KSPSetOperators(ksp_new_inner, K, K, DIFFERENT_NONZERO_PATTERN); KSPSetOptionsPrefix(ksp_new_inner, "rhsA11_"); MatSchurComplementSetKSP( S, ksp_new_inner );/* this call destroys the ksp_inner that is already set on S */ ksp_inner=ksp_new_inner; KSPGetPC( ksp_inner, &pcInner ); KSPSetFromOptions(ksp_inner); /* make sure we are setting up our solver how we want it */ } MatGetVecs( S, PETSC_NULL, &h_hat ); Vec f_tmp; /* It may be the case that the current velocity solution might not be bad guess for f_tmp? <-- maybe not */ MatGetVecs( K, PETSC_NULL, &f_tmp ); scrSolveTime = MPI_Wtime(); KSPSolve(ksp_inner, f, f_tmp); scrSolveTime = MPI_Wtime() - scrSolveTime; PetscPrintf( PETSC_COMM_WORLD, "\n\t* KSPSolve for RHS setup Finished in time: %lf seconds\n\n", scrSolveTime); MatMult(D, f_tmp, h_hat); VecAYPX(h_hat, -1.0, h); /* Computes y = x + alpha y. h_hat -> h - Gt*K^(-1)*f*/ Stg_VecDestroy(&f_tmp); if(bsscrp_self->solver->mg_active && change_A11rhspresolve) { //Stg_KSPDestroy(&ksp_inner ); KSPCreate(PETSC_COMM_WORLD, &ksp_new_inner); Stg_KSPSetOperators(ksp_new_inner, K, K, DIFFERENT_NONZERO_PATTERN); KSPSetOptionsPrefix( ksp_new_inner, "A11_" ); MatSchurComplementSetKSP( S, ksp_new_inner ); ksp_inner=ksp_new_inner; //MatSchurSetKSP( S, ksp_inner ); KSPGetPC( ksp_inner, &pcInner ); KSPSetFromOptions( ksp_inner ); mgSetupTime=setupMG( bsscrp_self, ksp_inner, pcInner, K, &mgCtx ); } /* create solver for S p = h_hat */ KSPCreate( PETSC_COMM_WORLD, &ksp_S ); KSPSetOptionsPrefix( ksp_S, "scr_"); /* By default use the UW approxS Schur preconditioner -- same as the one used by the Uzawa solver */ /* Note that if scaling is activated then the approxS matrix has been scaled already */ /* so no need to rebuild in the case of scaling as we have been doing */ if(!approxS){ PetscPrintf( PETSC_COMM_WORLD, "WARNING approxS is NULL\n"); } Stg_KSPSetOperators( ksp_S, S, S, SAME_NONZERO_PATTERN ); KSPSetType( ksp_S, "cg" ); KSPGetPC( ksp_S, &pc_S ); BSSCR_BSSCR_StokesCreatePCSchur2( K,G,D,C,approxS, pc_S, sym, bsscrp_self ); flg=0; PetscOptionsGetString( PETSC_NULL, "-NN", name, PETSC_MAX_PATH_LEN-1, &flg ); if(flg){ Mat Smat, Pmat; //MatStructure mstruct; Stg_PCGetOperators( pc_S, &Smat, &Pmat, NULL ); sprintf(str,"%s/",name); sprintf(matname,"Pmat%s",suffix); bsscr_dirwriteMat( Pmat, matname,str, "Writing Pmat matrix in al Solver"); } uzawastyle=PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-uzawa_style", &uzawastyle, &found ); if(uzawastyle){ /* now want to set up the ksp_S->pc to be of type ksp (gmres) by default to match Uzawa */ KSP pc_ksp; KSPGetPC( ksp_S, &pc_S ); PCSetType(pc_S,PCKSP); PCKSPGetKSP( pc_S, &pc_ksp); KSPSetType(pc_ksp, "gmres" ); KSPSetOptionsPrefix( pc_ksp, "scrPCKSP_"); KSPSetFromOptions( pc_ksp ); } KSPSetFromOptions( ksp_S ); /* Set specific monitor test */ KSPGetTolerances( ksp_S, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_it ); // Weirdness with petsc 3.2 here...look at it later //BSSCR_KSPLogSetMonitor( ksp_S, max_it, &monitor_index ); /***************************************************************************************************************/ /* Pressure / Velocity Solve */ /***************************************************************************************************************/ PetscPrintf( PETSC_COMM_WORLD, "\t* Pressure / Velocity Solve \n"); /***************************************************************************************************************/ /***************************************************************************************************************/ usePreviousGuess = PETSC_FALSE; if(been_here) PetscOptionsGetTruth( PETSC_NULL, "-scr_use_previous_guess", &usePreviousGuess, &found ); if(usePreviousGuess) { /* Note this should actually look at checkpoint information */ KSPSetInitialGuessNonzero( ksp_S, PETSC_TRUE ); } else { KSPSetInitialGuessNonzero( ksp_S, PETSC_FALSE ); } /***************************************************************************************************************/ /***************************************************************************************************************/ /******* SET CONVERGENCE TESTS *************************************************************************/ useNormInfStoppingConditions = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL ,"-scr_use_norm_inf_stopping_condition", &useNormInfStoppingConditions, &found ); if(useNormInfStoppingConditions) BSSCR_KSPSetNormInfConvergenceTest(ksp_S); useNormInfMonitor = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-scr_ksp_norm_inf_monitor", &useNormInfMonitor, &found ); if(useNormInfMonitor) KSPMonitorSet( ksp_S, BSSCR_KSPNormInfToNorm2Monitor, PETSC_NULL, PETSC_NULL ); /***************************************************************************************************************/ /***************************************************************************************************************/ /******* PRESSURE SOLVE ************************************************************************************/ PetscPrintf( PETSC_COMM_WORLD, "\t* KSPSolve( ksp_S, h_hat, p )\n"); /* if h_hat needs to be fixed up ..take out any nullspace vectors here */ /* we want to check that there is no "noise" in the null-space in the h vector */ /* this causes problems when we are trying to solve a Jacobian system when the Residual is almost converged */ if(bsscrp_self->check_pressureNS){ bsscrp_self->buildPNS(ksp);/* build and set nullspace vectors on bsscr - which is on ksp (function pointer is set in KSPSetUp_BSSCR) */ } PetscScalar hnorm, gnorm; MatNorm(G,NORM_INFINITY,&gnorm); VecNorm(h_hat, NORM_2, &hnorm); hnorm=hnorm/gnorm; if((hnorm < 1e-6) && (hnorm > 1e-20)){ VecScale(h_hat,1.0/hnorm); } /* test to see if v or t are in nullspace of G and orthogonalize wrt h_hat if needed */ KSPRemovePressureNullspace_BSSCR(ksp, h_hat); /* set convergence test to use min_it */ found = PETSC_FALSE; min_it = 0; PetscOptionsGetInt( PETSC_NULL,"-scr_ksp_set_min_it_converge", &min_it, &found); if(found && min_it > 0){ BSSCR_KSPSetConvergenceMinIts(ksp_S, min_it, bsscrp_self); } /** Pressure Solve **/ PetscGetFlops(&flopsA); scrSolveTime = MPI_Wtime(); KSPSolve( ksp_S, h_hat, p ); scrSolveTime = MPI_Wtime() - scrSolveTime; PetscGetFlops(&flopsB); PetscPrintf( PETSC_COMM_WORLD, "\n\t* KSPSolve( ksp_S, h_hat, p ) Solve Finished in time: %lf seconds\n\n", scrSolveTime); bsscrp_self->solver->stats.pressure_time=scrSolveTime; bsscrp_self->solver->stats.pressure_flops=(double)(flopsB-flopsA); /***************************************/ if((hnorm < 1e-6) && (hnorm > 1e-20)){ VecScale(h_hat,hnorm); VecScale(p,hnorm); } KSPRemovePressureNullspace_BSSCR(ksp, p); /***************************************************************************************************************/ /***************************************************************************************************************/ /* restore K and f for the Velocity back solve */ found = PETSC_FALSE; restorek = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-restore_K", &restorek, &found); //PetscOptionsGetString( PETSC_NULL, "-restore_K", name, PETSC_MAX_PATH_LEN-1, &flg ); if(penaltyNumber > 1e-10 && bsscrp_self->k2type){ if(restorek){ penaltyNumber = -penaltyNumber; if(f2) { ierr=VecAXPY(f,penaltyNumber,f2); }/* f <- f +a*f2 */ if(f3) { ierr=VecAXPY(f,penaltyNumber,f3); }/* f <- f +a*f3 */ ierr=MatAXPY(K,penaltyNumber,K2,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);/* Computes K = penaltyNumber*K2 + K */ //K=Korig; Stg_KSPSetOperators(ksp_inner, K, K, DIFFERENT_NONZERO_PATTERN); KisJustK=PETSC_TRUE; } } if(f3){ Stg_VecDestroy(&f3 ); } /* always destroy this local vector if was created */ /* obtain solution for u */ VecDuplicate( u, &t ); MatMult( G, p, t); VecAYPX( t, -1.0, f ); /*** t <- -t + f = f - G*p ***/ MatSchurComplementGetKSP( S, &ksp_inner ); a11SingleSolveTime = MPI_Wtime(); /* ---------------------------------- Final V Solve */ if(usePreviousGuess) KSPSetInitialGuessNonzero( ksp_inner, PETSC_TRUE ); /***************************************************************************************************************/ /***************************************************************************************************************/ /******* VELOCITY SOLVE ************************************************************************************/ /** Easier to just create a new KSP here if we want to do backsolve diffferently. (getting petsc errors now when switching from fgmres) */ change_backsolve=PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-change_backsolve", &change_backsolve, &found ); if(change_backsolve){ //Stg_KSPDestroy(&ksp_inner ); KSPCreate(PETSC_COMM_WORLD, &ksp_new_inner); Stg_KSPSetOperators(ksp_new_inner, K, K, DIFFERENT_NONZERO_PATTERN); KSPSetOptionsPrefix(ksp_new_inner, "backsolveA11_"); KSPSetFromOptions(ksp_new_inner); /* make sure we are setting up our solver how we want it */ MatSchurComplementSetKSP( S, ksp_new_inner );/* need to give the Schur it's inner ksp back for when we destroy it at end */ ksp_inner=ksp_new_inner; } PetscGetFlops(&flopsA); KSPSolve(ksp_inner, t, u); /* Solve, then restore default tolerance and initial guess */ PetscGetFlops(&flopsB); bsscrp_self->solver->stats.velocity_backsolve_flops=(double)(flopsB-flopsA); a11SingleSolveTime = MPI_Wtime() - a11SingleSolveTime; /* ------------------ Final V Solve */ bsscrp_self->solver->stats.velocity_backsolve_time=a11SingleSolveTime; flg=0; PetscOptionsGetString( PETSC_NULL, "-solutiondumpdir", name, PETSC_MAX_PATH_LEN-1, &flg ); if(flg){ sprintf(str,"%s/",name); sprintf(matname,"p%s",suffix); bsscr_dirwriteVec( p, matname,str, "Writing p vector in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"u%s",suffix); bsscr_dirwriteVec( u, matname,str, "Writing u vector in al Solver"); sprintf(str,"%s/",name); sprintf(matname,"h_hat%s",suffix); bsscr_dirwriteVec( h_hat, matname,str, "Writing h_hat vector in al Solver"); } found = PETSC_FALSE; restorek = PETSC_FALSE; PetscOptionsGetTruth( PETSC_NULL, "-restore_K_after_solve", &restorek, &found); if(penaltyNumber > 1e-10 && bsscrp_self->k2type){ if(restorek){ penaltyNumber = -penaltyNumber; ierr=MatAXPY(K,penaltyNumber,K2,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);/* Computes K = penaltyNumber*K2 + K */ KisJustK=PETSC_TRUE; } } /***************************************************************************************************************/ /***************************************************************************************************************/ /****** SOLUTION SUMMARY ******************************************************************************/ bsscr_summary(bsscrp_self,ksp_S,ksp_inner,K,K2,D,G,C,u,p,f,h,t,penaltyNumber,KisJustK, mgSetupTime, scrSolveTime, a11SingleSolveTime); //bsscr_summary(bsscrp_self,ksp_S,ksp_inner,K,Korig,K2,D,G,C,u,p,f,h,t,penaltyNumber,KisJustK, mgSetupTime, scrSolveTime, a11SingleSolveTime); /***************************************************************************************************************/ /***************************************************************************************************************/ Stg_VecDestroy(&t ); Stg_KSPDestroy(&ksp_S ); Stg_VecDestroy(&h_hat ); Stg_MatDestroy(&S );//This will destroy ksp_inner: also.. pcInner == pc_MG and is destroyed when ksp_inner is been_here = 1; PetscFunctionReturn(0); }