#include <slesolv.h>
Public Member Functions | |
void | initiate (slesolv *ssle) |
void | prepare_data (long *ndof, gtopology *top, FILE *out) |
void | print (FILE *out) |
void | read (gtopology *top, XFILE *in, long mespr) |
slesolv () | |
void | solve_system (gtopology *gt, gmatrix *gm, double *lhs, double *rhs, FILE *out) |
~slesolv () | |
Public Attributes | |
double | aerr |
attained norm of residual vector | |
long | ani |
number of performed iterations (atual number of iterations) | |
unsigned char | bsize |
block size in direct sparse solver | |
double | err |
required norm of residual vector | |
seqfeti | feti |
FETI solver. | |
long | iv |
long | nbdof |
number of rows which are not condensed | |
long | ni |
maximum number of iterations in iterative methods | |
precond | prec |
preconditioner | |
seqschur | schur |
solver based on the Schur complement method | |
saddpoint * | sp |
linsolvertype | tlinsol |
type of solver of system of linear algebraic equations | |
linsolvertype | tsol |
if tlinsol is conden, additional information is needed | |
double | zero |
computer zero |
class slesolv defines type of solver of systems of linear algebraic equations
JK
Definition at line 27 of file slesolv.h.
slesolv | ( | ) |
~slesolv | ( | ) |
Definition at line 45 of file slesolv.cpp.
void initiate | ( | slesolv * | ssle | ) |
void prepare_data | ( | long * | ndof, | |
gtopology * | top, | |||
FILE * | out | |||
) |
function prepares all necessary data for some solvers, namely for the Schur complement method and for the FETI method
top | - pointer to the general topology | |
out | - output file |
JK, 8.8.2010
Definition at line 320 of file slesolv.cpp.
References seqfeti::assemble_subdom_unknowns(), seqschur::assemble_subdom_unknowns(), gtopology::auto_subdom_detection(), seqtop::coarse_local_map(), feti, seqfeti::fetiimpl, seqtop::ggnbn, seqtop::ggnbncn, glob_glued, seqtop::icmultip, seqtop::icnbnmas, seqfeti::initiate(), seqschur::initiate(), seqtop::lnbncn, seqtop::md, seqtop::nnsd, seqtop::nodmultip, seqtop::ns, seqselnodes::prepare_feti(), seqselnodes::prepare_schur(), schur, seqtop::schur_ordering(), sfeti, seqtop::sid, sschur, gtopology::stop, tlinsol, and seqtop::tnbn.
Referenced by mefel_init(), and trfel_init().
void print | ( | FILE * | out | ) |
function prints
out | - output stream |
6.4.2003, JK
Definition at line 215 of file slesolv.cpp.
References bicg, cg, conden, err, feti, gauss_elim, ill, iv, lapack_sol, ldl, ll, lu, nbdof, ni, pardisolver, prec, precond::print(), seqfeti::print(), print_err(), saddle_point, sfeti, spdirldl, spdirll, spdirlu, and tlinsol.
Referenced by probdesct::print(), probdescc::print(), probdesc::print(), seqschur::print(), and aggregator::print().
function reads information about solver of linear equations
top | - pointer to the general topology | |
in | - input stream | |
mespr | - type of message printing |
6.4.2003, JK
Definition at line 59 of file slesolv.cpp.
References precond::agg, bicg, boss, cg, gtopology::cngen, conden, err, feti, gauss_elim, ill, aggregator::impl, iv, lapack_sol, ldl, linsolvertype_kwdset(), ll, lu, nbdof, ni, pardisolver, prec, print_err(), precond::pt, precond::read(), seqfeti::read(), seqschur::read(), gtopology::rst, saddle_point, schur, sfeti, spdirldl, spdirll, spdirlu, sschur, tlinsol, tsol, xfscanf(), and zero.
Referenced by probdesct::read(), probdescc::read(), probdesc::read(), seqschur::read(), and aggregator::read().
function solves system of linear algebraic equations
gt | - general topology | |
gm | - matrix of the system | |
lhs | - vector of unknowns | |
rhs | - vector of right hand side |
JK, 17.12.2006
Definition at line 430 of file slesolv.cpp.
References gnode::ai, bicg, cg, check_math_err(), conden, feti, gauss_elim, gtopology::gnodes, ill, precond::initiation(), lapack_sol, ldl, ll, long(), lu, gmatrix::n, gmatrix::nbdof, nbdof, gtopology::nidof, pardisolver, prec, print_err(), saddle_point, schur, sec(), sfeti, seqfeti::solve_system(), seqschur::solve_system(), gmatrix::solve_system(), spdirldl, spdirll, spdirlu, sschur, time, and tlinsol.
Referenced by arclengthrv(), arclengthrv1(), difference_method(), garclength(), gnewton_raphson(), gnewton_raphson_one_step(), inverse_iteration(), lin_nonstat_dform(), lin_nonstat_dform_resistance(), linear_nonstat_radiation_solv_dform(), linear_nonstat_solv_dform(), linear_nonstat_solv_dform_subcycl(), linear_nonstat_solv_vform(), newmark_method(), newton_raphson_coupl(), newton_raphson_coupl_new(), newton_raphson_gparcoupl_lin(), newton_raphson_gparcoupl_nonlin(), newton_raphson_parcoupl_lin(), newton_raphson_parcoupl_nonlin(), newton_raphson_parcoupl_nonlin_new(), newton_raphson_parcoupl_nonlin_old(), newton_raphson_restart(), newton_raphsonep(), newton_raphsonrv1(), nonlin_newmark_method(), nonlin_nonstat_dform(), nonlinear_nonstat_solv(), nonlinear_nonstat_solv_dform(), nonlinear_nonstat_solv_dform_dneska(), nonlinear_nonstat_solv_fnr_dform(), nonlinear_nonstat_solv_fnr_dform_old(), nonlinear_nonstat_solv_linesearch(), nonlinear_nonstat_solv_new(), nonlinear_nonstat_solv_nr_dform(), nonlinear_nonstat_solv_oldd(), nonlinear_nonstat_solv_pokus(), nonlinear_nonstat_solv_vform(), one_step_arcl(), one_step_linear(), one_step_nonlinear(), one_step_nonlinear_dform(), par_newton_raphson_parcoupl_lin(), paral_transport_homogenization(), solve_layered_linear_statics(), solve_linear_floating_subdomains(), solve_linear_floating_subdomains2(), solve_linear_statics(), solve_nonlinear_stationary_problem(), solve_nonlinear_stationary_problem_old(), solve_nonlinear_stationary_problem_pokus(), solve_nonstationary_growing_problem(), solve_nonstationary_growing_problem_nonlin(), solve_nonstationary_growing_vform(), solve_prob_constr_phases(), solve_radiation_stationary_problem(), schurcompl::solve_red_sys_fin(), solve_stationary_problem(), solve_var_stiff_method(), subspace_iter_jac(), subspace_iter_ortho(), subspace_shift_iter_ortho(), transport_homogenization(), and visco_solver().
double aerr |
attained norm of residual vector
Definition at line 59 of file slesolv.h.
Referenced by aggregator::boss(), initiate(), and slesolv().
long ani |
number of performed iterations (atual number of iterations)
Definition at line 57 of file slesolv.h.
Referenced by aggregator::boss(), initiate(), and slesolv().
unsigned char bsize |
block size in direct sparse solver
Definition at line 47 of file slesolv.h.
Referenced by initiate(), and slesolv().
double err |
required norm of residual vector
Definition at line 55 of file slesolv.h.
Referenced by aggregator::boss(), initiate(), print(), read(), gmatrix::setval(), and slesolv().
FETI solver.
Definition at line 74 of file slesolv.h.
Referenced by flsubdom::initiation(), prepare_data(), print(), transtop::read(), mechtop::read(), read(), solve_incremental_floating_subdomains(), flsubdom::solve_lin_alg_system(), and solve_system().
long iv |
initial values in iterative methods iv=0 - the input vector is multiplied by zero the iteration starts from zero vector iv=1 - the input vector is used as a start vector
Definition at line 65 of file slesolv.h.
Referenced by initiate(), print(), read(), gmatrix::setval(), and slesolv().
long nbdof |
long ni |
maximum number of iterations in iterative methods
Definition at line 53 of file slesolv.h.
Referenced by aggregator::boss(), initiate(), print(), read(), gmatrix::setval(), and slesolv().
preconditioner
Definition at line 71 of file slesolv.h.
Referenced by input(), inputt(), print(), read(), gmatrix::setval(), and solve_system().
solver based on the Schur complement method
Definition at line 77 of file slesolv.h.
Referenced by prepare_data(), read(), and solve_system().
type of solver of system of linear algebraic equations
Definition at line 41 of file slesolv.h.
Referenced by initiate(), input(), inputt(), prepare_data(), print(), read(), gmatrix::setval(), slesolv(), and solve_system().
if tlinsol is conden, additional information is needed
Definition at line 44 of file slesolv.h.
Referenced by initiate(), read(), gmatrix::setval(), and slesolv().
double zero |
computer zero
Definition at line 50 of file slesolv.h.
Referenced by aggregator::boss(), initiate(), read(), and slesolv().