Add HYPRE_MGRSetFSolverAtLevel (#983)
				
					
				
			Add hypre_MGRSetFSolverAtLevel and its public interface (HYPRE_). This function is used to set the F-relaxation solver at a specific level in MGR.
This commit is contained in:
		
							parent
							
								
									1593fbc474
								
							
						
					
					
						commit
						3e4087bac3
					
				| @ -4245,6 +4245,17 @@ HYPRE_Int HYPRE_MGRSetFSolver(HYPRE_Solver             solver, | ||||
|                               HYPRE_PtrToParSolverFcn  fine_grid_solver_setup, | ||||
|                               HYPRE_Solver             fsolver ); | ||||
| 
 | ||||
| /**
 | ||||
|  * (Optional) Set the F-relaxation solver at a given level. | ||||
|  * | ||||
|  * @param level [IN] MGR solver level | ||||
|  * @param solver [IN] MGR solver/preconditioner object | ||||
|  * @param fsolver [IN] F-relaxation solver object | ||||
|  **/ | ||||
| HYPRE_Int HYPRE_MGRSetFSolverAtLevel(HYPRE_Int     level, | ||||
|                                      HYPRE_Solver  solver, | ||||
|                                      HYPRE_Solver  fsolver ); | ||||
| 
 | ||||
| /**
 | ||||
|  * (Optional) Extract A_FF block from matrix A. | ||||
|  * | ||||
|  | ||||
| @ -222,6 +222,20 @@ HYPRE_MGRSetFSolver(HYPRE_Solver          solver, | ||||
|                                  (void *) fsolver ) ); | ||||
| } | ||||
| 
 | ||||
| /*--------------------------------------------------------------------------
 | ||||
|  * HYPRE_MGRSetFSolverAtLevel | ||||
|  *--------------------------------------------------------------------------*/ | ||||
| 
 | ||||
| HYPRE_Int | ||||
| HYPRE_MGRSetFSolverAtLevel(HYPRE_Int     level, | ||||
|                            HYPRE_Solver  solver, | ||||
|                            HYPRE_Solver  fsolver ) | ||||
| { | ||||
|    return ( hypre_MGRSetFSolverAtLevel( level, | ||||
|                                         (void *) solver, | ||||
|                                         (void *) fsolver ) ); | ||||
| } | ||||
| 
 | ||||
| /*--------------------------------------------------------------------------
 | ||||
|  * HYPRE_MGRBuildAff | ||||
|  *--------------------------------------------------------------------------*/ | ||||
|  | ||||
| @ -3601,6 +3601,7 @@ HYPRE_Int hypre_MGRSetFSolver( void *mgr_vdata, | ||||
|                                HYPRE_Int (*fine_grid_solver_solve)(void*, void*, void*, void*), | ||||
|                                HYPRE_Int (*fine_grid_solver_setup)(void*, void*, void*, void*), | ||||
|                                void *fsolver ); | ||||
| HYPRE_Int hypre_MGRSetFSolverAtLevel( HYPRE_Int level, void *mgr_vdata, void *fsolver ); | ||||
| HYPRE_Int hypre_MGRSetup( void *mgr_vdata, hypre_ParCSRMatrix *A, | ||||
|                           hypre_ParVector *f, hypre_ParVector *u ); | ||||
| HYPRE_Int hypre_MGRSolve( void *mgr_vdata, hypre_ParCSRMatrix *A, | ||||
|  | ||||
| @ -169,7 +169,8 @@ hypre_MGRCreate(void) | ||||
| HYPRE_Int | ||||
| hypre_MGRDestroy( void *data ) | ||||
| { | ||||
|    hypre_ParMGRData * mgr_data = (hypre_ParMGRData*) data; | ||||
|    hypre_ParMGRData  *mgr_data = (hypre_ParMGRData*) data; | ||||
|    hypre_Solver      *aff_base; | ||||
| 
 | ||||
|    HYPRE_Int i; | ||||
|    HYPRE_Int num_coarse_levels = (mgr_data -> num_coarse_levels); | ||||
| @ -350,15 +351,13 @@ hypre_MGRDestroy( void *data ) | ||||
|       { | ||||
|          if ((mgr_data -> aff_solver)[i]) | ||||
|          { | ||||
|             hypre_BoomerAMGDestroy((mgr_data -> aff_solver)[i]); | ||||
|             aff_base = (hypre_Solver*) (mgr_data -> aff_solver)[i]; | ||||
|             hypre_SolverDestroy(aff_base)((HYPRE_Solver) (aff_base)); | ||||
|          } | ||||
|       } | ||||
|       if (mgr_data -> fsolver_mode == 2) | ||||
|       { | ||||
|          if ((mgr_data -> aff_solver)[0]) | ||||
|          { | ||||
|             hypre_BoomerAMGDestroy((mgr_data -> aff_solver)[0]); | ||||
|          } | ||||
|          hypre_BoomerAMGDestroy((mgr_data -> aff_solver)[0]); | ||||
|       } | ||||
|       hypre_TFree(mgr_data -> aff_solver, HYPRE_MEMORY_HOST); | ||||
|       (mgr_data -> aff_solver) = NULL; | ||||
| @ -5267,6 +5266,51 @@ hypre_MGRSetFSolver( void  *mgr_vdata, | ||||
|    return hypre_error_flag; | ||||
| } | ||||
| 
 | ||||
| /*--------------------------------------------------------------------------
 | ||||
|  * hypre_MGRSetFSolverAtLevel | ||||
|  * | ||||
|  * set F-relaxation solver for a given MGR level. | ||||
|  * | ||||
|  * Note this function asks for a level identifier and doesn't expect an array | ||||
|  * of function pointers for each level (as done by SetLevel functions). | ||||
|  *--------------------------------------------------------------------------*/ | ||||
| 
 | ||||
| HYPRE_Int | ||||
| hypre_MGRSetFSolverAtLevel( HYPRE_Int   level, | ||||
|                             void       *mgr_vdata, | ||||
|                             void       *fsolver ) | ||||
| { | ||||
|    hypre_ParMGRData *mgr_data = (hypre_ParMGRData*) mgr_vdata; | ||||
| 
 | ||||
|    if (!mgr_data) | ||||
|    { | ||||
|       hypre_error_in_arg(1); | ||||
|       return hypre_error_flag; | ||||
|    } | ||||
|    HYPRE_Int        max_num_coarse_levels = (mgr_data -> max_num_coarse_levels); | ||||
|    HYPRE_Solver   **aff_solver = (mgr_data -> aff_solver); | ||||
| 
 | ||||
|    /* Check if the requested level makes sense */ | ||||
|    if (level < 0 || level >= max_num_coarse_levels) | ||||
|    { | ||||
|       hypre_error_in_arg(2); | ||||
|       return hypre_error_flag; | ||||
|    } | ||||
| 
 | ||||
|    /* Allocate aff_solver if needed */ | ||||
|    if (!aff_solver) | ||||
|    { | ||||
|       (mgr_data -> aff_solver) = aff_solver = hypre_CTAlloc(HYPRE_Solver*, | ||||
|                                                             max_num_coarse_levels, | ||||
|                                                             HYPRE_MEMORY_HOST); | ||||
|    } | ||||
| 
 | ||||
|    aff_solver[level] = (HYPRE_Solver *) fsolver; | ||||
|    (mgr_data -> fsolver_mode)  = 0; | ||||
| 
 | ||||
|    return hypre_error_flag; | ||||
| } | ||||
| 
 | ||||
| /* set coarse grid solver */ | ||||
| HYPRE_Int | ||||
| hypre_MGRSetCoarseSolver( void  *mgr_vdata, | ||||
|  | ||||
| @ -86,7 +86,8 @@ hypre_MGRSetup( void               *mgr_vdata, | ||||
|    hypre_ParCSRMatrix  *A_CF = NULL; | ||||
|    hypre_ParCSRMatrix  *A_CC = NULL; | ||||
| #endif | ||||
|    HYPRE_Solver **aff_solver = (mgr_data -> aff_solver); | ||||
|    hypre_Solver         *aff_base; | ||||
|    HYPRE_Solver        **aff_solver = (mgr_data -> aff_solver); | ||||
|    hypre_ParCSRMatrix  **A_ff_array = (mgr_data -> A_ff_array); | ||||
|    hypre_ParVector    **F_fine_array = (mgr_data -> F_fine_array); | ||||
|    hypre_ParVector    **U_fine_array = (mgr_data -> U_fine_array); | ||||
| @ -872,15 +873,14 @@ hypre_MGRSetup( void               *mgr_vdata, | ||||
|       { | ||||
|          if (aff_solver[j]) | ||||
|          { | ||||
|             hypre_BoomerAMGDestroy(aff_solver[j]); | ||||
|             aff_base = (hypre_Solver*) aff_solver[j]; | ||||
|             hypre_SolverDestroy(aff_base)((HYPRE_Solver) (aff_base)); | ||||
|             aff_solver[j] = NULL; | ||||
|          } | ||||
|       } | ||||
|       if (mgr_data -> fsolver_mode == 2) | ||||
|       { | ||||
|          if (aff_solver[0]) | ||||
|          { | ||||
|             hypre_BoomerAMGDestroy(aff_solver[0]); | ||||
|          } | ||||
|          hypre_BoomerAMGDestroy(aff_solver[0]); | ||||
|       } | ||||
|    } | ||||
| 
 | ||||
| @ -1527,7 +1527,7 @@ hypre_MGRSetup( void               *mgr_vdata, | ||||
|                   (mgr_data -> fsolver_mode) = 1; | ||||
|                } | ||||
|             } | ||||
|             else | ||||
|             else if (aff_solver[lev]) | ||||
|             { | ||||
|                hypre_sprintf(msg, "Warning!! User-prescribed F-solver for the first level\n\
 | ||||
|                              reduction (set in HYPRE_MGRSetFSolver()) only supports AMG\n\ | ||||
| @ -1536,6 +1536,19 @@ hypre_MGRSetup( void               *mgr_vdata, | ||||
|                hypre_error_w_msg(0, msg); | ||||
|             } | ||||
|          } | ||||
|          else if (aff_solver[lev]) | ||||
|          { | ||||
|             aff_base = (hypre_Solver*) aff_solver[lev]; | ||||
| 
 | ||||
|             /* Save A_FF splitting */ | ||||
|             A_ff_array[lev] = A_FF; | ||||
| 
 | ||||
|             /* Call setup function */ | ||||
|             hypre_SolverSetup(aff_base)((HYPRE_Solver) aff_solver[lev], | ||||
|                                         (HYPRE_Matrix) A_ff_array[lev], | ||||
|                                         (HYPRE_Vector) F_fine_array[lev + 1], | ||||
|                                         (HYPRE_Vector) U_fine_array[lev + 1]); | ||||
|          } | ||||
|          else if (Frelax_type[lev] == 2) /* Construct default AMG solver */ | ||||
|          { | ||||
|             /* Save A_FF splitting */ | ||||
|  | ||||
| @ -25,9 +25,9 @@ hypre_MGRSolve( void               *mgr_vdata, | ||||
| { | ||||
| 
 | ||||
|    MPI_Comm              comm = hypre_ParCSRMatrixComm(A); | ||||
|    hypre_ParMGRData   *mgr_data = (hypre_ParMGRData*) mgr_vdata; | ||||
|    hypre_ParMGRData     *mgr_data = (hypre_ParMGRData*) mgr_vdata; | ||||
| 
 | ||||
|    hypre_ParCSRMatrix  **A_array = (mgr_data -> A_array); | ||||
|    hypre_ParCSRMatrix **A_array = (mgr_data -> A_array); | ||||
|    hypre_ParVector    **F_array = (mgr_data -> F_array); | ||||
|    hypre_ParVector    **U_array = (mgr_data -> U_array); | ||||
| 
 | ||||
| @ -35,10 +35,10 @@ hypre_MGRSolve( void               *mgr_vdata, | ||||
|    HYPRE_Int            logging = (mgr_data -> logging); | ||||
|    HYPRE_Int            print_level = (mgr_data -> print_level); | ||||
|    HYPRE_Int            max_iter = (mgr_data -> max_iter); | ||||
|    HYPRE_Real           *norms = (mgr_data -> rel_res_norms); | ||||
|    hypre_ParVector      *Vtemp = (mgr_data -> Vtemp); | ||||
|    HYPRE_Real          *norms = (mgr_data -> rel_res_norms); | ||||
|    hypre_ParVector     *Vtemp = (mgr_data -> Vtemp); | ||||
|    //   hypre_ParVector      *Utemp = (mgr_data -> Utemp);
 | ||||
|    hypre_ParVector      *residual; | ||||
|    hypre_ParVector     *residual; | ||||
| 
 | ||||
|    HYPRE_Complex        fp_zero = 0.0; | ||||
|    HYPRE_Complex        fp_one = 1.0; | ||||
| @ -507,6 +507,7 @@ hypre_MGRCycle( void              *mgr_vdata, | ||||
| { | ||||
|    MPI_Comm               comm; | ||||
|    hypre_ParMGRData      *mgr_data = (hypre_ParMGRData*) mgr_vdata; | ||||
|    hypre_Solver          *aff_base; | ||||
| 
 | ||||
|    HYPRE_Int              local_size; | ||||
|    HYPRE_Int              level; | ||||
| @ -984,10 +985,23 @@ hypre_MGRCycle( void              *mgr_vdata, | ||||
|             if (Frelax_type[level] == 2) | ||||
|             { | ||||
|                /* Do F-relaxation using AMG */ | ||||
|                fine_grid_solver_solve((mgr_data -> aff_solver)[fine_grid], | ||||
|                                       A_ff_array[fine_grid], | ||||
|                                       F_fine_array[coarse_grid], | ||||
|                                       U_fine_array[coarse_grid]); | ||||
|                if (level == 0) | ||||
|                { | ||||
|                   /* TODO (VPM): unify with the next block */ | ||||
|                   fine_grid_solver_solve((mgr_data -> aff_solver)[fine_grid], | ||||
|                                          A_ff_array[fine_grid], | ||||
|                                          F_fine_array[coarse_grid], | ||||
|                                          U_fine_array[coarse_grid]); | ||||
|                } | ||||
|                else | ||||
|                { | ||||
|                   aff_base = (hypre_Solver*) (mgr_data -> aff_solver)[level]; | ||||
| 
 | ||||
|                   hypre_SolverSolve(aff_base)((HYPRE_Solver) (mgr_data -> aff_solver)[level], | ||||
|                                               (HYPRE_Matrix) A_ff_array[level], | ||||
|                                               (HYPRE_Vector) F_fine_array[level + 1], | ||||
|                                               (HYPRE_Vector) U_fine_array[level + 1]); | ||||
|                } | ||||
|             } | ||||
|             else | ||||
|             { | ||||
|  | ||||
| @ -2181,6 +2181,7 @@ HYPRE_Int hypre_MGRSetFSolver( void *mgr_vdata, | ||||
|                                HYPRE_Int (*fine_grid_solver_solve)(void*, void*, void*, void*), | ||||
|                                HYPRE_Int (*fine_grid_solver_setup)(void*, void*, void*, void*), | ||||
|                                void *fsolver ); | ||||
| HYPRE_Int hypre_MGRSetFSolverAtLevel( HYPRE_Int level, void *mgr_vdata, void *fsolver ); | ||||
| HYPRE_Int hypre_MGRSetup( void *mgr_vdata, hypre_ParCSRMatrix *A, | ||||
|                           hypre_ParVector *f, hypre_ParVector *u ); | ||||
| HYPRE_Int hypre_MGRSolve( void *mgr_vdata, hypre_ParCSRMatrix *A, | ||||
|  | ||||
		Loading…
	
		Reference in New Issue
	
	Block a user
	 Victor A. P. Magri
						Victor A. P. Magri