My Project
programmer's documentation
Functions
cs_multigrid.c File Reference
#include "cs_defs.h"
#include <stdarg.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <assert.h>
#include <math.h>
#include "bft_mem.h"
#include "bft_error.h"
#include "bft_printf.h"
#include "cs_base.h"
#include "cs_blas.h"
#include "cs_file.h"
#include "cs_grid.h"
#include "cs_halo.h"
#include "cs_log.h"
#include "cs_matrix.h"
#include "cs_matrix_default.h"
#include "cs_mesh.h"
#include "cs_mesh_quantities.h"
#include "cs_post.h"
#include "cs_sles.h"
#include "cs_sles_it.h"
#include "cs_sles_pc.h"
#include "cs_timer.h"
#include "cs_time_plot.h"
#include "cs_time_step.h"
#include "cs_multigrid.h"
Include dependency graph for cs_multigrid.c:

Functions

void cs_multigrid_initialize (void)
 Initialize multigrid solver API. More...
 
void cs_multigrid_finalize (void)
 Finalize multigrid solver API. More...
 
bool cs_multigrid_needed (void)
 Indicate if multigrid solver API is used for at least one system. More...
 
cs_multigrid_tcs_multigrid_define (int f_id, const char *name, cs_multigrid_type_t mg_type)
 Define and associate a multigrid sparse linear system solver for a given field or equation name. More...
 
cs_multigrid_tcs_multigrid_create (cs_multigrid_type_t mg_type)
 Create multigrid linear system solver info and context. More...
 
void cs_multigrid_destroy (void **context)
 Destroy multigrid linear system solver info and context. More...
 
void * cs_multigrid_copy (const void *context)
 Create multigrid sparse linear system solver info and context based on existing info and context. More...
 
void cs_multigrid_log (const void *context, cs_log_t log_type)
 Log multigrid solver info. More...
 
void cs_multigrid_set_coarsening_options (cs_multigrid_t *mg, int aggregation_limit, int coarsening_type, int n_max_levels, cs_gnum_t min_g_rows, double p0p1_relax, int postprocess)
 Set multigrid coarsening parameters. More...
 
void cs_multigrid_set_solver_options (cs_multigrid_t *mg, cs_sles_it_type_t descent_smoother_type, cs_sles_it_type_t ascent_smoother_type, cs_sles_it_type_t coarse_solver_type, int n_max_cycles, int n_max_iter_descent, int n_max_iter_ascent, int n_max_iter_coarse, int poly_degree_descent, int poly_degree_ascent, int poly_degree_coarse, double precision_mult_descent, double precision_mult_ascent, double precision_mult_coarse)
 Set multigrid parameters for associated iterative solvers. More...
 
cs_sles_it_type_t cs_multigrid_get_fine_solver_type (const cs_multigrid_t *mg)
 Return solver type used on fine mesh. More...
 
void cs_multigrid_setup (void *context, const char *name, const cs_matrix_t *a, int verbosity)
 Setup multigrid sparse linear equation solver. More...
 
void cs_multigrid_setup_conv_diff (void *context, const char *name, const cs_matrix_t *a, const cs_matrix_t *a_conv, const cs_matrix_t *a_diff, int verbosity)
 Setup multigrid sparse linear equation solver. More...
 
cs_sles_convergence_state_t cs_multigrid_solve (void *context, const char *name, const cs_matrix_t *a, int verbosity, cs_halo_rotation_t rotation_mode, double precision, double r_norm, int *n_iter, double *residue, const cs_real_t *rhs, cs_real_t *vx, size_t aux_size, void *aux_vectors)
 Call multigrid sparse linear equation solver. More...
 
void cs_multigrid_free (void *context)
 Free multigrid sparse linear equation solver setup context. More...
 
cs_sles_pc_tcs_multigrid_pc_create (cs_multigrid_type_t mg_type)
 Create a multigrid preconditioner. More...
 
bool cs_multigrid_error_post_and_abort (cs_sles_t *sles, cs_sles_convergence_state_t state, const cs_matrix_t *a, cs_halo_rotation_t rotation_mode, const cs_real_t rhs[], cs_real_t vx[])
 Error handler for multigrid sparse linear equation solver. More...
 
void cs_multigrid_set_plot_options (cs_multigrid_t *mg, const char *base_name, bool use_iteration)
 Set plotting options for multigrid. More...
 

Function Documentation

◆ cs_multigrid_copy()

void* cs_multigrid_copy ( const void *  context)

Create multigrid sparse linear system solver info and context based on existing info and context.

Parameters
[in]contextpointer to reference info and context (actual type: cs_multigrid_t *)
Returns
pointer to newly created solver info object (actual type: cs_multigrid_t *)

◆ cs_multigrid_create()

cs_multigrid_t* cs_multigrid_create ( cs_multigrid_type_t  mg_type)

Create multigrid linear system solver info and context.

The multigrid variant is an ACM (Additive Corrective Multigrid) method.

Parameters
[in]mg_typetype of multigrid algorithm to use
Returns
pointer to new multigrid info and context

◆ cs_multigrid_define()

cs_multigrid_t* cs_multigrid_define ( int  f_id,
const char *  name,
cs_multigrid_type_t  mg_type 
)

Define and associate a multigrid sparse linear system solver for a given field or equation name.

If this system did not previously exist, it is added to the list of "known" systems. Otherwise, its definition is replaced by the one defined here.

This is a utility function: if finer control is needed, see cs_sles_define and cs_multigrid_create.

Note that this function returns a pointer directly to the multigrid solver management structure. This may be used to set further options, for example calling cs_multigrid_set_coarsening_options and cs_multigrid_set_solver_options. If needed, cs_sles_find may be used to obtain a pointer to the matching cs_sles_t container.

Parameters
[in]f_idassociated field id, or < 0
[in]nameassociated name if f_id < 0, or NULL
[in]mg_typetype of multigrid algorithm to use
Returns
pointer to new multigrid info and context

◆ cs_multigrid_destroy()

void cs_multigrid_destroy ( void **  context)

Destroy multigrid linear system solver info and context.

Parameters
[in,out]contextpointer to multigrid linear solver info (actual type: cs_multigrid_t **)

◆ cs_multigrid_error_post_and_abort()

bool cs_multigrid_error_post_and_abort ( cs_sles_t sles,
cs_sles_convergence_state_t  state,
const cs_matrix_t a,
cs_halo_rotation_t  rotation_mode,
const cs_real_t  rhs[],
cs_real_t  vx[] 
)

Error handler for multigrid sparse linear equation solver.

In case of divergence or breakdown, this error handler outputs postprocessing data to assist debugging, then aborts the run. It does nothing in case the maximum iteration count is reached.

Parameters
[in,out]slespointer to solver object
[in]stateconvergence state
[in]amatrix
[in]rotation_modehalo update option for rotational periodicity
[in]rhsright hand side
[in,out]vxsystem solution
Returns
false (do not attempt new solve)

◆ cs_multigrid_finalize()

void cs_multigrid_finalize ( void  )

Finalize multigrid solver API.

◆ cs_multigrid_free()

void cs_multigrid_free ( void *  context)

Free multigrid sparse linear equation solver setup context.

This function frees resolution-related data, incuding the current grid hierarchy, but does not free the whole context, as info used for logging (especially performance data) is maintained.

Parameters
[in,out]contextpointer to multigrid solver info and context (actual type: cs_multigrid_t *)

◆ cs_multigrid_get_fine_solver_type()

cs_sles_it_type_t cs_multigrid_get_fine_solver_type ( const cs_multigrid_t mg)

Return solver type used on fine mesh.

Parameters
[in]mgpointer to multigrid info and context
Returns
type of smoother for descent (used for fine mesh)

◆ cs_multigrid_initialize()

void cs_multigrid_initialize ( void  )

Initialize multigrid solver API.

◆ cs_multigrid_log()

void cs_multigrid_log ( const void *  context,
cs_log_t  log_type 
)

Log multigrid solver info.

Parameters
[in]contextpointer to iterative solver info and context (actual type: cs_multigrid_t *)
[in]log_typelog type

◆ cs_multigrid_needed()

bool cs_multigrid_needed ( void  )

Indicate if multigrid solver API is used for at least one system.

Returns
true if at least one system uses a multigrid solver, false otherwise

◆ cs_multigrid_pc_create()

cs_sles_pc_t* cs_multigrid_pc_create ( cs_multigrid_type_t  mg_type)

Create a multigrid preconditioner.

Parameters
[in]mg_typetype of multigrid algorithm to use
Returns
pointer to newly created preconditioner object.

◆ cs_multigrid_set_coarsening_options()

void cs_multigrid_set_coarsening_options ( cs_multigrid_t mg,
int  aggregation_limit,
int  coarsening_type,
int  n_max_levels,
cs_gnum_t  min_g_rows,
double  p0p1_relax,
int  postprocess 
)

Set multigrid coarsening parameters.

Parameters
[in,out]mgpointer to multigrid info and context
[in]aggregation_limitmaximum allowed fine rows per coarse row
[in]coarsening_typecoarsening type; see cs_grid_coarsening_t
[in]n_max_levelsmaximum number of grid levels
[in]min_g_rowsglobal number of rows on coarse grids under which no coarsening occurs
[in]p0p1_relaxp0/p1 relaxation_parameter
[in]postprocessif > 0, postprocess coarsening (uses coarse row numbers modulo this value)

◆ cs_multigrid_set_plot_options()

void cs_multigrid_set_plot_options ( cs_multigrid_t mg,
const char *  base_name,
bool  use_iteration 
)

Set plotting options for multigrid.

Parameters
[in,out]mgpointer to multigrid info and context
[in]base_namebase plot name to activate, NULL otherwise
[in]use_iterationif true, use iteration as time stamp otherwise, use wall clock time

◆ cs_multigrid_set_solver_options()

void cs_multigrid_set_solver_options ( cs_multigrid_t mg,
cs_sles_it_type_t  descent_smoother_type,
cs_sles_it_type_t  ascent_smoother_type,
cs_sles_it_type_t  coarse_solver_type,
int  n_max_cycles,
int  n_max_iter_descent,
int  n_max_iter_ascent,
int  n_max_iter_coarse,
int  poly_degree_descent,
int  poly_degree_ascent,
int  poly_degree_coarse,
double  precision_mult_descent,
double  precision_mult_ascent,
double  precision_mult_coarse 
)

Set multigrid parameters for associated iterative solvers.

Parameters
[in,out]mgpointer to multigrid info and context
[in]descent_smoother_typetype of smoother for descent
[in]ascent_smoother_typetype of smoother for ascent
[in]coarse_solver_typetype of solver for coarsest grid
[in]n_max_cyclesmaximum number of cycles
[in]n_max_iter_descentmaximum iterations per descent smoothing
[in]n_max_iter_ascentmaximum iterations per ascent smmothing
[in]n_max_iter_coarsemaximum iterations per coarsest solution
[in]poly_degree_descentpreconditioning polynomial degree for descent phases (0: diagonal)
[in]poly_degree_ascentpreconditioning polynomial degree for ascent phases (0: diagonal)
[in]poly_degree_coarsepreconditioning polynomial degree for coarse solver (0: diagonal)
[in]precision_mult_descentprecision multiplier for descent smoothers (levels >= 1)
[in]precision_mult_ascentprecision multiplier for ascent smoothers
[in]precision_mult_coarseprecision multiplier for coarsest grid

◆ cs_multigrid_setup()

void cs_multigrid_setup ( void *  context,
const char *  name,
const cs_matrix_t a,
int  verbosity 
)

Setup multigrid sparse linear equation solver.

Parameters
[in,out]contextpointer to multigrid solver info and context (actual type: cs_multigrid_t *)
[in]namepointer to name of linear system
[in]aassociated matrix
[in]verbosityassociated verbosity

◆ cs_multigrid_setup_conv_diff()

void cs_multigrid_setup_conv_diff ( void *  context,
const char *  name,
const cs_matrix_t a,
const cs_matrix_t a_conv,
const cs_matrix_t a_diff,
int  verbosity 
)

Setup multigrid sparse linear equation solver.

Parameters
[in,out]contextpointer to multigrid solver info and context (actual type: cs_multigrid_t *)
[in]namepointer to name of linear system
[in]aassociated matrix
[in]a_convassociated matrix (convection)
[in]a_diffassociated matrix (diffusion)
[in]verbosityassociated verbosity

◆ cs_multigrid_solve()

cs_sles_convergence_state_t cs_multigrid_solve ( void *  context,
const char *  name,
const cs_matrix_t a,
int  verbosity,
cs_halo_rotation_t  rotation_mode,
double  precision,
double  r_norm,
int *  n_iter,
double *  residue,
const cs_real_t rhs,
cs_real_t vx,
size_t  aux_size,
void *  aux_vectors 
)

Call multigrid sparse linear equation solver.

Parameters
[in,out]contextpointer to multigrid solver info and context (actual type: cs_multigrid_t *)
[in]namepointer to name of linear system
[in]amatrix
[in]verbosityassociated verbosity
[in]rotation_modehalo update option for rotational periodicity
[in]precisionsolver precision
[in]r_normresidue normalization
[out]n_iternumber of "equivalent" iterations
[out]residueresidue
[in]rhsright hand side
[in,out]vxsystem solution
[in]aux_sizesize of aux_vectors (in bytes)
aux_vectorsoptional working area (internal allocation if NULL)
Returns
convergence state