Arcane  v3.14.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
AlephParams.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2022 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* AlephParams.cc (C) 2010 */
9/* */
10/*---------------------------------------------------------------------------*/
11/*---------------------------------------------------------------------------*/
12#include "AlephArcane.h"
13
14/*---------------------------------------------------------------------------*/
15/*---------------------------------------------------------------------------*/
16
17ARCANE_BEGIN_NAMESPACE
18
19/*---------------------------------------------------------------------------*/
20/*---------------------------------------------------------------------------*/
21
22/*---------------------------------------------------------------------------*/
23/*---------------------------------------------------------------------------*/
24
25AlephParams::
26AlephParams()
27: TraceAccessor(nullptr)
28, m_param_epsilon(1.0e-10)
29, m_param_max_iteration(1024)
30, m_param_preconditioner_method(TypesSolver::DIAGONAL)
31, m_param_solver_method(TypesSolver::PCG)
32, m_param_gamma(-1)
33, m_param_alpha(-1.0)
34, m_param_xo_user(false)
35, m_param_check_real_residue(false)
36, m_param_print_real_residue(false)
37, m_param_debug_info(false)
38, m_param_min_rhs_norm(1.e-20)
39, m_param_convergence_analyse(false)
40, m_param_stop_error_strategy(true)
41, m_param_write_matrix_to_file_error_strategy(false)
42, m_param_write_matrix_name_error_strategy("SolveErrorAlephMatrix.dbg")
43, m_param_listing_output(false)
44, m_param_threshold(0.0)
45, m_param_print_cpu_time_resolution(false)
46, m_param_amg_coarsening_method(0)
47, m_param_output_level(0)
48, m_param_amg_cycle(1)
49, m_param_amg_solver_iterations(1)
50, m_param_amg_smoother_iterations(1)
51, m_param_amg_smootherOption(TypesSolver::SymHybGSJ_smoother)
52, m_param_amg_coarseningOption(TypesSolver::ParallelRugeStuben)
53, m_param_amg_coarseSolverOption(TypesSolver::CG_coarse_solver)
54, m_param_keep_solver_structure(false)
55, m_param_sequential_solver(false)
56, m_param_criteria_stop(TypesSolver::RB)
57{
58 //debug() << "\33[1;4;33m\t[AlephParams] New"<<"\33[0m";
59}
60
61AlephParams::
62AlephParams(ITraceMng* tm,
63 Real epsilon, // epsilon de convergence
64 Integer max_iteration, // nb max iterations
65 TypesSolver::ePreconditionerMethod preconditioner_method, // préconditionnement utilisé (defaut DIAG)
66 TypesSolver::eSolverMethod solver_method, // méthode de résolution par defaut PCG
67 Integer gamma, // destine au parametrage des préconditionnements
68 Real alpha, // destine au parametrage des préconditionnements
69 bool xo_user, // permet a l'utilisateur d'initialiser le PGC avec un Xo different de zero
70 bool check_real_residue,
71 bool print_real_residue,
72 bool debug_info,
73 Real min_rhs_norm,
74 bool convergence_analyse,
75 bool stop_error_strategy,
76 bool write_matrix_to_file_error_strategy,
77 String write_matrix_name_error_strategy,
78 bool listing_output,
79 Real threshold,
80 bool print_cpu_time_resolution,
81 Integer amg_coarsening_method,
82 Integer output_level,
83 Integer amg_cycle,
84 Integer amg_solver_iterations,
85 Integer amg_smoother_iterations,
86 TypesSolver::eAmgSmootherOption amg_smootherOption,
87 TypesSolver::eAmgCoarseningOption amg_coarseningOption,
88 TypesSolver::eAmgCoarseSolverOption amg_coarseSolverOption,
89 bool keep_solver_structure,
90 bool sequential_solver,
91 TypesSolver::eCriteriaStop param_criteria_stop)
92: TraceAccessor(tm)
93, m_param_epsilon(epsilon)
94, m_param_max_iteration(max_iteration)
95, m_param_preconditioner_method(preconditioner_method)
96, m_param_solver_method(solver_method)
97, m_param_gamma(gamma)
98, m_param_alpha(alpha)
99, m_param_xo_user(xo_user)
100, m_param_check_real_residue(check_real_residue)
101, m_param_print_real_residue(print_real_residue)
102, m_param_debug_info(debug_info)
103, m_param_min_rhs_norm(min_rhs_norm)
104, m_param_convergence_analyse(convergence_analyse)
105, m_param_stop_error_strategy(stop_error_strategy)
106, m_param_write_matrix_to_file_error_strategy(write_matrix_to_file_error_strategy)
107, m_param_write_matrix_name_error_strategy(write_matrix_name_error_strategy)
108, m_param_listing_output(listing_output)
109, m_param_threshold(threshold)
110, m_param_print_cpu_time_resolution(print_cpu_time_resolution)
111, m_param_amg_coarsening_method(amg_coarsening_method)
112, m_param_output_level(output_level)
113, m_param_amg_cycle(amg_cycle)
114, m_param_amg_solver_iterations(amg_solver_iterations)
115, m_param_amg_smoother_iterations(amg_smoother_iterations)
116, m_param_amg_smootherOption(amg_smootherOption)
117, m_param_amg_coarseningOption(amg_coarseningOption)
118, m_param_amg_coarseSolverOption(amg_coarseSolverOption)
119, m_param_keep_solver_structure(keep_solver_structure)
120, m_param_sequential_solver(sequential_solver)
121, m_param_criteria_stop(param_criteria_stop)
122{
123 //debug() << "\33[1;4;33m\t[AlephParams] New"<<"\33[0m";
124}
125
126AlephParams::
127~AlephParams()
128{
129 //debug() << "\33[1;4;33m\t[~AlephParams]"<<"\33[0m";
130}
131
132// set
133void AlephParams::setEpsilon(const Real epsilon)
134{
135 m_param_epsilon = epsilon;
136}
137void AlephParams::setMaxIter(const Integer max_iteration)
138{
139 m_param_max_iteration = max_iteration;
140}
141void AlephParams::setPrecond(const TypesSolver::ePreconditionerMethod preconditioner_method)
142{
143 m_param_preconditioner_method = preconditioner_method;
144}
145void AlephParams::setMethod(const TypesSolver::eSolverMethod solver_method)
146{
147 m_param_solver_method = solver_method;
148}
149void AlephParams::setAlpha(const Real alpha)
150{
151 m_param_alpha = alpha;
152}
153void AlephParams::setGamma(const Integer gamma)
154{
155 m_param_gamma = gamma;
156}
157void AlephParams::setXoUser(const bool xo_user)
158{
159 m_param_xo_user = xo_user;
160}
161void AlephParams::setCheckRealResidue(const bool check_real_residue)
162{
163 m_param_check_real_residue = check_real_residue;
164}
165void AlephParams::setPrintRealResidue(const bool print_real_residue)
166{
167 m_param_print_real_residue = print_real_residue;
168}
169void AlephParams::setDebugInfo(const bool debug_info)
170{
171 m_param_debug_info = debug_info;
172}
173void AlephParams::setMinRHSNorm(const Real min_rhs_norm)
174{
175 m_param_min_rhs_norm = min_rhs_norm;
176}
177void AlephParams::setConvergenceAnalyse(const bool convergence_analyse)
178{
179 m_param_convergence_analyse = convergence_analyse;
180}
181void AlephParams::setStopErrorStrategy(const bool stop_error_strategy)
182{
183 m_param_stop_error_strategy = stop_error_strategy;
184}
185void AlephParams::setWriteMatrixToFileErrorStrategy(const bool write_matrix_to_file_error_strategy)
186{
187 m_param_write_matrix_to_file_error_strategy = write_matrix_to_file_error_strategy;
188}
189void AlephParams::setWriteMatrixNameErrorStrategy(const String& write_matrix_name_error_strategy)
190{
191 m_param_write_matrix_name_error_strategy = write_matrix_name_error_strategy;
192}
193void AlephParams::setDDMCParameterListingOutput(const bool listing_output)
194{
195 m_param_listing_output = listing_output;
196}
197void AlephParams::setDDMCParameterAmgDiagonalThreshold(const Real threshold)
198{
199 m_param_threshold = threshold;
200}
201void AlephParams::setPrintCpuTimeResolution(const bool print_cpu_time_resolution)
202{
203 m_param_print_cpu_time_resolution = print_cpu_time_resolution;
204}
205
206void AlephParams::
207setAmgCoarseningMethod(const TypesSolver::eAmgCoarseningMethod method)
208{
209 switch (method) {
210 case TypesSolver::AMG_COARSENING_AUTO:
211 m_param_amg_coarsening_method = 6;
212 break;
213 case TypesSolver::AMG_COARSENING_HYPRE_0:
214 m_param_amg_coarsening_method = 0;
215 break;
216 case TypesSolver::AMG_COARSENING_HYPRE_1:
217 m_param_amg_coarsening_method = 1;
218 break;
219 case TypesSolver::AMG_COARSENING_HYPRE_3:
220 m_param_amg_coarsening_method = 3;
221 break;
222 case TypesSolver::AMG_COARSENING_HYPRE_6:
223 m_param_amg_coarsening_method = 6;
224 break;
225 case TypesSolver::AMG_COARSENING_HYPRE_7:
226 m_param_amg_coarsening_method = 7;
227 break;
228 case TypesSolver::AMG_COARSENING_HYPRE_8:
229 m_param_amg_coarsening_method = 8;
230 break;
231 case TypesSolver::AMG_COARSENING_HYPRE_9:
232 m_param_amg_coarsening_method = 9;
233 break;
234 case TypesSolver::AMG_COARSENING_HYPRE_10:
235 m_param_amg_coarsening_method = 10;
236 break;
237 case TypesSolver::AMG_COARSENING_HYPRE_11:
238 m_param_amg_coarsening_method = 11;
239 break;
240 case TypesSolver::AMG_COARSENING_HYPRE_21:
241 m_param_amg_coarsening_method = 21;
242 break;
243 case TypesSolver::AMG_COARSENING_HYPRE_22:
244 m_param_amg_coarsening_method = 22;
245 break;
246 default:
247 throw NotImplementedException(A_FUNCINFO);
248 }
249}
250void AlephParams::setOutputLevel(const Integer output_level)
251{
252 m_param_output_level = output_level;
253}
254void AlephParams::setAmgCycle(const Integer amg_cycle)
255{
256 m_param_amg_cycle = amg_cycle;
257}
258void AlephParams::setAmgSolverIter(const Integer amg_solver_iterations)
259{
260 m_param_amg_solver_iterations = amg_solver_iterations;
261}
262void AlephParams::setAmgSmootherIter(const Integer amg_smoother_iterations)
263{
264 m_param_amg_smoother_iterations = amg_smoother_iterations;
265}
266void AlephParams::setAmgSmootherOption(const TypesSolver::eAmgSmootherOption amg_smootherOption)
267{
268 m_param_amg_smootherOption = amg_smootherOption;
269}
270void AlephParams::setAmgCoarseningOption(const TypesSolver::eAmgCoarseningOption amg_coarseningOption)
271{
272 m_param_amg_coarseningOption = amg_coarseningOption;
273}
274void AlephParams::setAmgCoarseSolverOption(const TypesSolver::eAmgCoarseSolverOption amg_coarseSolverOption)
275{
276 m_param_amg_coarseSolverOption = amg_coarseSolverOption;
277}
278void AlephParams::setKeepSolverStructure(const bool keep_solver_structure)
279{
280 m_param_keep_solver_structure = keep_solver_structure;
281}
282void AlephParams::setSequentialSolver(const bool sequential_solver)
283{
284 m_param_sequential_solver = sequential_solver;
285}
286void AlephParams::setCriteriaStop(const TypesSolver::eCriteriaStop criteria_stop)
287{
288 m_param_criteria_stop = criteria_stop;
289}
290
291// get
292Real AlephParams::epsilon() const
293{
294 return m_param_epsilon;
295}
296int AlephParams::maxIter() const
297{
298 return m_param_max_iteration;
299}
300Real AlephParams::alpha() const
301{
302 return m_param_alpha;
303}
304int AlephParams::gamma() const
305{
306 return m_param_gamma;
307}
308TypesSolver::ePreconditionerMethod AlephParams::precond()
309{
310 return m_param_preconditioner_method;
311}
312TypesSolver::eSolverMethod AlephParams::method()
313{
314 return m_param_solver_method;
315}
316bool AlephParams::xoUser() const
317{
318 return m_param_xo_user;
319}
320bool AlephParams::checkRealResidue() const
321{
322 return m_param_check_real_residue;
323}
324bool AlephParams::printRealResidue() const
325{
326 return m_param_print_real_residue;
327}
328bool AlephParams::debugInfo()
329{
330 return m_param_debug_info;
331}
332Real AlephParams::minRHSNorm()
333{
334 return m_param_min_rhs_norm;
335}
336bool AlephParams::convergenceAnalyse()
337{
338 return m_param_convergence_analyse;
339}
340bool AlephParams::stopErrorStrategy()
341{
342 return m_param_stop_error_strategy;
343}
344bool AlephParams::writeMatrixToFileErrorStrategy()
345{
346 return m_param_write_matrix_to_file_error_strategy;
347}
348String AlephParams::writeMatrixNameErrorStrategy()
349{
350 return m_param_write_matrix_name_error_strategy;
351}
352bool AlephParams::DDMCParameterListingOutput() const
353{
354 return m_param_listing_output;
355}
356Real AlephParams::DDMCParameterAmgDiagonalThreshold() const
357{
358 return m_param_threshold;
359}
360bool AlephParams::printCpuTimeResolution() const
361{
362 return m_param_print_cpu_time_resolution;
363}
364int AlephParams::amgCoarseningMethod() const
365{
366 return m_param_amg_coarsening_method;
367} // -1 pour Sloop
368int AlephParams::getOutputLevel() const
369{
370 return m_param_output_level;
371}
372int AlephParams::getAmgCycle() const
373{
374 return m_param_amg_cycle;
375}
376int AlephParams::getAmgSolverIter() const
377{
378 return m_param_amg_solver_iterations;
379}
380int AlephParams::getAmgSmootherIter() const
381{
382 return m_param_amg_smoother_iterations;
383}
384TypesSolver::eAmgSmootherOption AlephParams::getAmgSmootherOption() const
385{
386 return m_param_amg_smootherOption;
387}
388TypesSolver::eAmgCoarseningOption AlephParams::getAmgCoarseningOption() const
389{
390 return m_param_amg_coarseningOption;
391}
392TypesSolver::eAmgCoarseSolverOption AlephParams::getAmgCoarseSolverOption() const
393{
394 return m_param_amg_coarseSolverOption;
395}
396bool AlephParams::getKeepSolverStructure() const
397{
398 return m_param_keep_solver_structure;
399}
400bool AlephParams::getSequentialSolver() const
401{
402 return m_param_sequential_solver;
403}
404TypesSolver::eCriteriaStop AlephParams::getCriteriaStop() const
405{
406 return m_param_criteria_stop;
407}
408
409/*---------------------------------------------------------------------------*/
410/*---------------------------------------------------------------------------*/
411
412ARCANE_END_NAMESPACE
413
414/*---------------------------------------------------------------------------*/
415/*---------------------------------------------------------------------------*/
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120