8 #ifndef OPENFPM_NUMERICS_SRC_SOLVERS_PETSC_SOLVER_AMG_REPORT_HPP_
9 #define OPENFPM_NUMERICS_SRC_SOLVERS_PETSC_SOLVER_AMG_REPORT_HPP_
16 #include "Solvers/petsc_solver.hpp"
17 #include "Vector/map_vector.hpp"
86 auto x_ = solver.
solve(A,b);
98 double time1 = tm_solve.
getwct();
99 double time2 = tm_solve2.
getwct();
115 std::cout <<
"Time1: " << time1 <<
" Time2: " << time2 << std::endl;
147 A.getMatrixTriplets();
148 if (v_cl.
getProcessUnitID() == 0) {std::cout <<
"Benchmarking HMIS coarsener" << std::endl;}
154 A.getMatrixTriplets();
155 if (v_cl.
getProcessUnitID() == 0) {std::cout <<
"Benchmarking PMIS coarsener" << std::endl;}
161 A.getMatrixTriplets();
162 if (v_cl.
getProcessUnitID() == 0) {std::cout <<
"Benchmarking Falgout coarsener" << std::endl;}
168 A.getMatrixTriplets();
169 if (v_cl.
getProcessUnitID() == 0) {std::cout <<
"Benchmarking modifiedRuge-Stueben coarsener" << std::endl;}
172 method.add(
"modifiedRuge-Stueben");
175 A.getMatrixTriplets();
176 if (v_cl.
getProcessUnitID() == 0) {std::cout <<
"Benchmarking Ruge-Stueben coarsener" << std::endl;}
179 method.add(
"Ruge-Stueben");
182 A.getMatrixTriplets();
183 if (v_cl.
getProcessUnitID() == 0) {std::cout <<
"Benchmarking CLJP coarsener" << std::endl;}
209 std::cout <<
"Composing coarsening report" << std::endl;
221 yn.add(
"Norm infinity");
222 yn.add(
"time to setup");
223 yn.add(
"time to solve");
227 y.add({perf_amg_coars.get(i).residual,perf_amg_coars.get(i).time_setup,perf_amg_coars.get(i).time_solve});
232 options.
title = std::string(
"Overview of different coarsener on a V cycle with 6 relaxation SOR/Jacobi sweeps for each level up and down");
233 options.
yAxis = std::string(
"Error/time");
234 options.
xAxis = std::string(
"Method");
235 options.
stype = std::string(
"bars");
236 options.
more = std::string(
"vAxis: {scaleType: 'log'}");
238 std::stringstream ss;
240 cg.
addHTML(
"<h2>Coarsener</h2>");
252 cg.
addHTML(
"<h2>V cycle sweep</h2>");
280 options.
title = std::string(
"Coarsener: " + coars.get(k)) + std::string(
": Overview of V cycle with n sweeps symmetrically up and down");
281 options.
yAxis = std::string(
"error/time/score");
282 options.
xAxis = std::string(
"sweeps");
283 options.
stype = std::string(
"bars");
284 options.
more = std::string(
"vAxis: {scaleType: 'log'}");
286 std::stringstream ss;
300 cg.
addHTML(
"<h2>V cycle sweep asymmetric test</h2>");
301 for(
size_t k = 0 ; k < coars.
size() ; k++)
331 options.
title = std::string(
"Coarsener: " + coars.get(k)) + std::string(
": Overview of V cycle with n sweeps asymmetrically distributed");
332 options.
yAxis = std::string(
"error/time/score");
333 options.
xAxis = std::string(
"sweeps");
334 options.
stype = std::string(
"bars");
335 options.
more = std::string(
"vAxis: {scaleType: 'log'}");
337 std::stringstream ss;
366 for (
size_t k = 0 ; k < coarsener_to_test.
size(); k++)
380 A.getMatrixTriplets();
388 std::cout <<
"Coarsener: " << coarsener_to_test.get(k) <<
" Sweep: " << i <<
" Score: " << score << std::endl;
421 for (
size_t k = 0 ; k < coarsener_to_test.
size(); k++)
424 size_t num_tot_sweep = 2*ids_ts.get(k);
425 for (
size_t i = 0 ; i <= num_tot_sweep ; i++)
437 A.getMatrixTriplets();
445 std::cout <<
"Coarsener: " << coarsener_to_test.get(k) <<
" Sweep: (" << i <<
"," << num_tot_sweep-i <<
") Score: " << score << std::endl;
462 if (create_vcluster().getProcessUnitID() != 0)
473 cg.
write(
"gc_AMG_preconditioners.html");
503 coars.add(
method.get(fastest));
504 coars.add(
method.get(accurate));
518 size_t page = perf.
size() / nm;
520 for (
size_t k = 0 ; k < nm ; k++)
522 int score_id = page*k;
525 for (
size_t i = page*k ; i < page*k+page ; i++)
Small class to produce graph with Google chart in HTML.
void write(std::string file)
It write the graphs on file in html format using Google charts.
void AddHistGraph(openfpm::vector< Y > &y)
Add an histogram graph.
void addHTML(const std::string &html)
Add HTML text.
Sparse Matrix implementation.
void execute()
Execute all the requests.
size_t getProcessUnitID()
Get the process unit id.
void max(T &num)
Get the maximum number across all processors (or reduction with infinity norm)
Implementation of VCluster class.
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support.
Class to test AMG solvers.
void test_cycle_asym(SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b, openfpm::vector< size_t > &ids_ts, openfpm::vector< std::string > &coarsener_to_test)
test best asymmetric combination of sweeps
void test_coarsener(SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b)
test the corasener for this problem
double target_accuracy
Target accuracy for the preconditioner.
void setMaxSweep(int max_sw)
Set the target accuracy to score the AMG solver.
double score_solver(double t_solve, double t_m)
Score the solver.
void benchmark(SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b, petsc_solver< double > &solver, openfpm::vector< AMG_time_err_coars > &perf_amg)
benchmark the selected setting for the preconditioner
openfpm::vector< AMG_time_err_coars > perf_amg_sweep_sym
It contain the performance of several AMG methods sweep configuration.
openfpm::vector< std::string > method
List of methods.
openfpm::vector< std::pair< size_t, size_t > > v_cycle_tested_asym
List of tested cycles in case of asymmetric.
openfpm::vector< size_t > v_cycle_tested_sym
List of tested cycles.
void test_cycle_type(SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b, openfpm::vector< std::string > &coarsener_to_test)
test the corasener for this problem
openfpm::vector< size_t > asym_tests
starting point of each coarsening test for the asymmetric test
void write_report(openfpm::vector< std::string > &coars)
void write_report_cycle(GoogleChart &cg, openfpm::vector< std::string > &coars)
Write the report for coarsening.
size_t num_sweep_test
Number of sweeps for test.
void best_coarsener(openfpm::vector< std::string > &coars)
take the most accurate and the fastest AMG
void best_score(openfpm::vector< AMG_time_err_coars > &perf, openfpm::vector< size_t > &sw_optimal, size_t nm)
Return the best scoring solver.
openfpm::vector< AMG_time_err_coars > perf_amg_coars
It contains the performance of several AMG methods coarsener.
void try_solve(SparseMatrix< double, int, PETSC_BASE > &A, const Vector< double, PETSC_BASE > &b)
Try to use AMG pre-conditioner and check how they they perform.
void setTergetAMGAccuracy(double t_a)
Set the target accuracy to score the AMG solver.
void write_report_cycle_asym(GoogleChart &cg, openfpm::vector< std::string > &coars)
Write the report for coarsening.
void write_report_coars(GoogleChart &cg)
Write the report for coarsening.
openfpm::vector< AMG_time_err_coars > perf_amg_sweep_asym
It contain the performance of several AMG methods sweep configuration.
This class is able to do Matrix inversion in parallel with PETSC solvers.
Vector< T, PETSC_BASE > solve(SparseMatrix< T, int, PETSC_BASE > &A, const Vector< T, PETSC_BASE > &b, bool initial_guess=false)
Here we invert the matrix and solve the system.
void setPreconditionerAMG_cycleType(const std::string &cycle_type, int sweep_up=-1, int sweep_dw=-1, int sweep_crs=-1)
It set the type of cycle and optionally the number of sweep.
solError get_residual_error(SparseMatrix< T, int, PETSC_BASE > &A, const Vector< T, PETSC_BASE > &x, const Vector< T, PETSC_BASE > &b)
It return the resiual error.
void setPreconditionerAMG_nl(int nl)
Set the number of levels for the algebraic-multigrid preconditioner.
void log_monitor()
Set the Petsc solver.
void setPreconditionerAMG_coarsenNodalType(int norm)
Set the block coarsening norm type.
void setPreconditioner(PCType type)
Set the preconditioner of the linear solver.
void setPreconditionerAMG_coarsen(const std::string &type)
Set the coarsening method for the algebraic-multi-grid preconditioner.
void setSolver(KSPType type)
Set the Petsc solver.
void setPreconditionerAMG_relax(const std::string &type, int k=REL_ALL)
Set the relaxation method for the algebraic-multi-grid preconditioner.
void setMaxIter(PetscInt n)
Set the maximum number of iteration for Krylov solvers.
void setAbsTol(PetscReal abstol_)
Set the absolute tolerance as stop criteria.
void setPreconditionerAMG_maxit(int nit)
Set the maximum number of V or W cycle for algebraic-multi-grid.
Class for cpu time benchmarking.
void stop()
Stop the timer.
void start()
Start the timer.
double getwct()
Return the elapsed real time.
It contain information about the performance of the AMG.
double score
assign a score to the solver
double time_solve
time to solve
double residual
residual norm
double time_setup
time to setup
std::string xAxis
X axis name.
std::string title
Title of the chart.
std::string yAxis
Y axis name.
It contain statistic of the error of the calculated solution.
PetscReal err_inf
infinity norm of the error