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 BSSCR_DRIVER_flex( KSP ksp, Mat stokes_A, Vec stokes_x, Vec stokes_b, Mat approxS, KSP ksp_K, MatStokesBlockScaling BA, PetscTruth sym, KSP_BSSCR * bsscrp_self ) { char name[PETSC_MAX_PATH_LEN]; char ubefore[100]; char uafter[100]; char pbefore[100]; char pafter[100]; PetscTruth flg, flg2, truth, useAcceleratingSmoothingMG, useFancySmoothingMG; PetscTruth usePreviousGuess, useNormInfStoppingConditions, useNormInfMonitor, found, extractMats; Mat K,G,D,C; Vec u,p,f,h; Mat S; Vec h_hat,t,t2,q,v; KSP ksp_inner; KSP ksp_S; KSP ksp_cleaner; KSPType ksp_inner_type; PetscTruth has_cnst_nullspace; PC pc_S, pc_MG, pcInner; PetscInt monitor_index,max_it,min_it; Vec nsp_vec = PETSC_NULL; PetscReal scr_rtol; PetscReal inner_rtol; PetscReal vSolver_rtol; PetscScalar uNormInf, pNormInf; PetscScalar uNorm, pNorm, rNorm, fNorm; PetscInt uSize, pSize; PetscInt lmin,lmax; PetscInt iterations; PetscReal min,max; PetscReal p_sum; MGContext mgCtx; PC shellPC; double t0, t1; double mgSetupTime, scrSolveTime, a11SingleSolveTime, solutionAnalysisTime; Index nx,ny,nz; PetscInt j,start,end; static int been_here = 0; /* Ha Ha Ha !! */ /* get sub matrix / vector objects */ MatNestGetSubMat( stokes_A, 0,0, &K ); MatNestGetSubMat( stokes_A, 0,1, &G ); MatNestGetSubMat( stokes_A, 1,0, &D ); 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, "\t Adress of stokes_x is %p\n", stokes_x); */ /* VecNorm( u, NORM_2, &uNorm ); */ /* PetscPrintf( PETSC_COMM_WORLD, "\t u Norm is %.6e in %s: address is %p\n",uNorm,__func__,u); */ /* VecNorm( p, NORM_2, &pNorm ); */ /* PetscPrintf( PETSC_COMM_WORLD, "\t p Norm is %.6e in %s: addres is %p\n",pNorm,__func__,p); */ /* Create Schur complement matrix */ //MatCreateSchurFromBlock( stokes_A, 0.0, "MatSchur_A11", &S ); MatCreateSchurComplement(K,K,G,D,C, &S); /* configure inner solver */ if (ksp_K!=PETSC_NULL) { MatSchurComplementSetKSP( S, ksp_K ); MatSchurComplementGetKSP( S, &ksp_inner ); } else { abort(); MatSchurComplementGetKSP( S, &ksp_inner ); KSPSetType( ksp_inner, "cg" ); } KSPGetPC( ksp_inner, &pcInner ); /* If we're using multigrid, replace the preconditioner here so we get the same options prefix. */ if(bsscrp_self->mg) { mgSetupTime=setupMG( bsscrp_self, ksp_inner, pcInner, K, &mgCtx ); } /* SETFROMOPTIONS MIGHT F**K MG UP */ KSPSetOptionsPrefix( ksp_inner, "A11_" ); KSPSetFromOptions( ksp_inner ); 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 ); /* create right hand side */ /* h_hat = G'*inv(K)*f - h */ MatGetVecs(K,PETSC_NULL,&t); MatGetVecs( S, PETSC_NULL, &h_hat ); KSPSetOptionsPrefix( ksp_inner, "A11_" ); KSPSetFromOptions( ksp_inner ); KSPSolve(ksp_inner,f,t);/* t=f/K */ //bsscr_writeVec( t, "ts", "Writing t vector"); MatMult(D,t,h_hat);/* G'*t */ VecAXPY(h_hat, -1.0, h);/* h_hat = h_hat - h */ Stg_VecDestroy(&t); //bsscr_writeVec( h_hat, "h_hat", "Writing h_hat Vector in Solver"); //MatSchurApplyReductionToVecFromBlock( S, stokes_b, h_hat ); /* create solver for S p = h_hat */ KSPCreate( PETSC_COMM_WORLD, &ksp_S ); KSPSetOptionsPrefix( ksp_S, "scr_"); Stg_KSPSetOperators( ksp_S, S,S, SAME_NONZERO_PATTERN ); KSPSetType( ksp_S, "cg" ); /* Build preconditioner for S */ KSPGetPC( ksp_S, &pc_S ); BSSCR_BSSCR_StokesCreatePCSchur2( K,G,D,C,approxS,pc_S,sym, bsscrp_self ); KSPSetFromOptions(ksp_S); /* Set specific monitor test */ KSPGetTolerances( ksp_S, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_it ); //BSSCR_KSPLogSetMonitor( ksp_S, max_it, &monitor_index ); /* Pressure / Velocity Solve */ scrSolveTime = MPI_Wtime(); 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 ); } //KSPSetRelativeRhsConvergenceTest( ksp_S ); 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 ); 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 norm, a, a1, a2, hnorm, pnorm, 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); } KSPSolve( ksp_S, h_hat, p ); sprintf(pafter,"psafter_%d",been_here); bsscr_writeVec( p, pafter, "Writing p Vector in Solver"); /***************************************/ if((hnorm < 1e-6) && (hnorm > 1e-20)){ VecScale(h_hat,hnorm); VecScale(p,hnorm); } KSPRemovePressureNullspace_BSSCR(ksp, p); scrSolveTime = MPI_Wtime() - scrSolveTime; PetscPrintf( PETSC_COMM_WORLD, "\n\t* KSPSolve( ksp_S, h_hat, p ) Solve Finished in time: %lf seconds\n\n", scrSolveTime); /* Resolve with this pressure to obtain solution for u */ /* obtain solution for u */ VecDuplicate( u, &t ); MatMult( G, p, t); VecAYPX( t, -1.0, f ); /* t <- -t + f */ MatSchurComplementGetKSP( S, &ksp_inner ); a11SingleSolveTime = MPI_Wtime(); /* ---------------------------------- Final V Solve */ if(usePreviousGuess) KSPSetInitialGuessNonzero( ksp_inner, PETSC_TRUE ); KSPSetOptionsPrefix( ksp_inner, "backsolveA11_" ); KSPSetFromOptions( ksp_inner ); KSPSolve( ksp_inner, t, u ); /* Solve, then restore default tolerance and initial guess */ a11SingleSolveTime = MPI_Wtime() - a11SingleSolveTime; /* ------------------ Final V Solve */ PetscPrintf( PETSC_COMM_WORLD, "\n\nSCR Solver Summary:\n\n"); if(bsscrp_self->mg) PetscPrintf( PETSC_COMM_WORLD, " Multigrid setup: = %.4g secs \n", mgSetupTime); KSPGetIterationNumber( ksp_S, &iterations); PetscPrintf( PETSC_COMM_WORLD, " Pressure Solve: = %.4g secs / %d its\n", scrSolveTime, iterations); KSPGetIterationNumber( ksp_inner, &iterations); PetscPrintf( PETSC_COMM_WORLD, " Final V Solve: = %.4g secs / %d its\n\n", a11SingleSolveTime, iterations); /* Analysis of solution: This can be somewhat time consuming as it requires allocation / de-allocation, computing vector norms etc. So we make it optional.. This should be put into a proper KSP monitor now? */ flg = PETSC_TRUE; PetscOptionsGetTruth( PETSC_NULL, "-scr_ksp_solution_summary", &flg, &found ); if(flg) { solutionAnalysisTime = MPI_Wtime(); VecGetSize( u, &uSize ); VecGetSize( p, &pSize ); VecDuplicate( u, &t2 ); MatMult( K, u, t2); VecAYPX( t2, -1.0, t ); /* t2 <- -t2 + t ... should be the formal residual vector */ VecNorm( t2, NORM_2, &rNorm ); VecNorm( f, NORM_2, &fNorm ); PetscPrintf( PETSC_COMM_WORLD, " |f - K u - G p|/|f| = %.6e\n", rNorm/fNorm ); VecDuplicate( p, &q ); MatMult( D, u, q ); VecNorm( u, NORM_2, &uNorm ); VecNorm( q, NORM_2, &rNorm ); PetscPrintf( PETSC_COMM_WORLD, " |G^T u|_2/|u|_2 = %.6e\n", sqrt( (double) uSize / (double) pSize ) * rNorm / uNorm); VecNorm( q, NORM_INFINITY, &rNorm ); PetscPrintf( PETSC_COMM_WORLD, " |G^T u|_infty/|u|_2 = %.6e\n", sqrt( (double) uSize ) * rNorm / uNorm); VecNorm( u, NORM_INFINITY, &uNormInf ); VecNorm( u, NORM_2, &uNorm ); VecGetSize( u, &uSize ); VecNorm( p, NORM_INFINITY, &pNormInf ); VecNorm( p, NORM_2, &pNorm ); PetscPrintf( PETSC_COMM_WORLD, " |u|_{\\infty} = %.6e , u_rms = %.6e\n", uNormInf, uNorm / sqrt( (double) uSize ) ); PetscPrintf( PETSC_COMM_WORLD, " |p|_{\\infty} = %.6e , p_rms = %.6e\n", pNormInf, pNorm / sqrt( (double) pSize ) ); VecMax( u, &lmax, &max ); VecMin( u, &lmin, &min ); PetscPrintf( PETSC_COMM_WORLD, " min/max(u) = %.6e [%d] / %.6e [%d]\n",min,lmin,max,lmax); VecMax( p, &lmax, &max ); VecMin( p, &lmin, &min ); PetscPrintf( PETSC_COMM_WORLD, " min/max(p) = %.6e [%d] / %.6e [%d]\n",min,lmin,max,lmax); VecSum( p, &p_sum ); PetscPrintf( PETSC_COMM_WORLD, " \\sum_i p_i = %.6e \n", p_sum ); solutionAnalysisTime = MPI_Wtime() - solutionAnalysisTime; PetscPrintf( PETSC_COMM_WORLD, "\n Time for this analysis = %.4g secs\n\n",solutionAnalysisTime); Stg_VecDestroy(&t2 ); Stg_VecDestroy(&q ); } if(bsscrp_self->mg) { //MG_inner_solver_pcmg_shutdown( pcInner ); } Stg_VecDestroy(&t ); // KSPLogDestroyMonitor( ksp_S ); Stg_KSPDestroy(&ksp_S ); //Stg_KSPDestroy(&ksp_inner ); Stg_VecDestroy(&h_hat ); Stg_MatDestroy(&S ); /* Destroy nullspace vector if it exists. */ if(nsp_vec) Stg_VecDestroy(&nsp_vec); //been_here = 1; been_here++; PetscFunctionReturn(0); }