#include <seqschur.h>
Public Member Functions | |
void | assemble_subdom_unknowns (gtopology *top, FILE *out) |
void | globcnnum_pdd (gtopology *top, long *ltg, long *domproc, FILE *out) |
void | initiate (seqselnodes *selnodschur) |
void | print (FILE *out) |
void | read (gtopology *top, long mespr, XFILE *in) |
void | schurordering_new (gtopology *top, long *nodeidentif) |
seqschur () | |
void | solve_red_sys (double **condmat, double **condvect, FILE *out) |
void | solve_red_sys_fin (double **condmat, double **condvect, FILE *out) |
void | solve_red_sys_iter (double **condmat, double **condvect, FILE *out) |
void | solve_system (gtopology *top, gmatrix *gm, double *lhs, double *rhs, FILE *out) |
void | subdomain_matrices (gmatrix *gm, FILE *out) |
~seqschur () | |
Public Attributes | |
double | aerrcg |
attained error | |
long | anicg |
number of performed iterations in conjugate gradient method | |
gmatrix * | arr |
long ** | cndom |
storagetype | dmstor |
type of storage of domain matrix | |
double | errcg |
required error | |
gtopology * | gtop |
generalized topology, each subdomain is assumed as generalized element (superelement) | |
long | maxnbdof |
long | nbdof |
number of boundary degrees of freedom (unknowns) on subdomain | |
long * | nbdofmas |
long | ndof |
number of degrees of freedom on subdomain | |
long | ndofcp |
long * | ndofdom |
long | nicg |
maximum number of iterations in conjugate gradient method | |
long | nidof |
number of internal degrees of freedom (unknowns) on subdomain | |
long | ns |
number of subdomains | |
long * | nsid |
storagetype | rsmstor |
storage of reduced system matrix | |
skyline * | smsky |
subdomain matrices stored in the skyline storage | |
slesolv * | ssle |
type of linear system solver | |
redsystsolver | trssol |
type of reduced system solver | |
double | zero |
computer zero |
class seqschur deals with the Schur complement method this is single-processor implementation, it is different from schur in PARGEF which is multi-processor implementation
JK, 13.5.2009
Definition at line 20 of file seqschur.h.
seqschur | ( | ) |
~seqschur | ( | ) |
Definition at line 66 of file seqschur.cpp.
References ssle.
void assemble_subdom_unknowns | ( | gtopology * | top, | |
FILE * | out | |||
) |
function orders unknowns in coarse problem
function initializates object gtop of the class gtopology it represents subdomains as superelements
out | - output file |
JK, 8.7.2005 function assembles list of unknowns belonging to the subdomains
ndofdom - numbers of DOFs on subdomains cndom - list of DOFs on subdomains
top | - pointer to topology | |
out | - output file |
JK, 21.5.2009
Definition at line 575 of file seqschur.cpp.
References cndom, g, seqtop::ggnbn, seqtop::ggnin, gtopology::give_dof(), gtopology::give_ndofn(), seqtop::nbnd, ndofdom, seqtop::nind, ns, and gtopology::stop.
Referenced by slesolv::prepare_data().
void globcnnum_pdd | ( | gtopology * | top, | |
long * | ltg, | |||
long * | domproc, | |||
FILE * | out | |||
) |
void initiate | ( | seqselnodes * | selnodschur | ) |
function generates ordering of unknowns on subdomains
top | - pointer to the sequential general topology | |
ltg | - local to global array | |
out | - output stream |
JK, 13.5.2009 function orders unknowns on subdomains with respect of Schur complement method inner variables are ordered at the beginning boundary variables are ordered at he end
top | - pointer to domain topology | |
ptop | - pointer to topology for paralle computation |
JK, 13.5.2009 function initializes variables and arrays
namely:
ndofcp - number of coarse DOFs maxnbdof - maximum number of boundary/interface DOFs on subdomain nbdofmas - array of numbers of boundary/interface DOFs on subdomains
selnodschur | - pointer to the object of the class seqselnodes |
JK, 13.5.2009
Definition at line 404 of file seqschur.cpp.
References seqselnodes::cndofmas, gtopology::cnstate, gtop, gtopology::initiate(), maxnbdof, nbdofmas, ndofcp, ns, seqselnodes::snndofmas, and seqselnodes::tndofsn.
Referenced by slesolv::prepare_data().
void print | ( | FILE * | out | ) |
function prints basic data
out | - output file |
JK, 13.5.2009
Definition at line 120 of file seqschur.cpp.
References errcg, master_sol, nicg, ns, paral_sol, slesolv::print(), print_err(), rsmstor, ssle, and trssol.
function reads basic data
i | - type of solver of system of equations | |
in | - input file |
JK, 13.5.2009
Definition at line 80 of file seqschur.cpp.
References errcg, master_sol, nicg, ns, paral_sol, print_err(), slesolv::read(), rsmstor, ssle, trssol, and xfscanf().
Referenced by slesolv::read().
void schurordering_new | ( | gtopology * | top, | |
long * | nodeidentif | |||
) |
void solve_red_sys | ( | double ** | condmat, | |
double ** | condvect, | |||
FILE * | out | |||
) |
function solves reduced system of linear algebraic equations
condmat | - array containing reduced matrices (Schur complements) | |
condvect | - array containing reduced vectors | |
out | - output file (for auxiliary output) |
JK, 13.5.2009
Definition at line 1108 of file seqschur.cpp.
References master_sol, paral_sol, print_err(), solve_red_sys_fin(), solve_red_sys_iter(), and trssol.
Referenced by solve_system().
void solve_red_sys_fin | ( | double ** | condmat, | |
double ** | condvect, | |||
FILE * | out | |||
) |
function solves reduced system of linear algebraic equations on the master processor by a direct method
condmat | - array containing reduced matrices (Schur complements) | |
condvect | - array containing reduced vectors | |
out | - output file (for auxiliary output) |
JK, 13.5.2009
Definition at line 1041 of file seqschur.cpp.
References arr, gtopology::give_code_numbers(), gtopology::give_ndofe(), globloc(), gtop, gmatrix::initiate(), gmatrix::localized(), locglob(), memset(), nbdofmas, ndofcp, ns, nullv(), gmatrix::prepmat(), rsmstor, gmatrix::setval(), gmatrix::solve_system(), and ssle.
Referenced by solve_red_sys().
void solve_red_sys_iter | ( | double ** | condmat, | |
double ** | condvect, | |||
FILE * | out | |||
) |
function solves reduced system of equations from primal domain decomposition by conjugate gradient method matrix of reduced system is not assembled iteration process is directed by the master processor
condmat | - array containing reduced matrices (Schur complements) | |
condvect | - array containing reduced vectors | |
out | - output file (for auxiliary output) |
JK, 13.5.2009
Definition at line 813 of file seqschur.cpp.
Referenced by solve_red_sys().
function solves system of linear algebraic equations by the Schur complment method
top | - pointer to the sequential general topology | |
gm | - pointer to the subdomain matrix stored in the form gmatrix | |
lhs | - array containing solution | |
rhs | - array containing right hand side | |
out | - output file (for auxiliary output) |
JK, 13.5.2009
Definition at line 1138 of file seqschur.cpp.
References cndom, condmat(), globloc(), skyline::ldlkon_sky(), nbdofmas, ndofdom, ns, smsky, solve_red_sys(), subdomain_matrices(), and time.
Referenced by slesolv::solve_system().
void subdomain_matrices | ( | gmatrix * | gm, | |
FILE * | out | |||
) |
function assembles subdomain matrices from the matrix of the whole system
top | - general topology | |
gm | - pointer to the matrix of the system | |
out | - output file |
JK, 13.5.2009
Definition at line 743 of file seqschur.cpp.
References symcomprow::a, symcomprow::adr, skyline::assemble_from_scr(), symcomprow::ci, cndom, skyline::n, ndofdom, ns, gmatrix::scr, and smsky.
Referenced by solve_system().
double aerrcg |
long anicg |
number of performed iterations in conjugate gradient method
Definition at line 69 of file seqschur.h.
Referenced by seqschur().
Definition at line 139 of file seqschur.h.
Referenced by seqschur(), and solve_red_sys_fin().
long** cndom |
list of DOFs on subdomains cndom[i][j]=k - the j-th DOF on the i-th subdomain has number k
Definition at line 123 of file seqschur.h.
Referenced by assemble_subdom_unknowns(), solve_system(), and subdomain_matrices().
type of storage of domain matrix
Definition at line 131 of file seqschur.h.
double errcg |
required error
Definition at line 71 of file seqschur.h.
Referenced by print(), read(), and seqschur().
generalized topology, each subdomain is assumed as generalized element (superelement)
Definition at line 96 of file seqschur.h.
Referenced by initiate(), seqschur(), and solve_red_sys_fin().
long maxnbdof |
Definition at line 86 of file seqschur.h.
Referenced by initiate(), and seqschur().
long nbdof |
number of boundary degrees of freedom (unknowns) on subdomain
Definition at line 83 of file seqschur.h.
Referenced by seqschur().
long* nbdofmas |
Definition at line 93 of file seqschur.h.
Referenced by initiate(), seqschur(), solve_red_sys_fin(), and solve_system().
long ndof |
number of degrees of freedom on subdomain
Definition at line 79 of file seqschur.h.
Referenced by seqschur().
long ndofcp |
number of coarse DOFs = total number of boundary/interface DOFs number of DOFs of the coarse problem
Definition at line 90 of file seqschur.h.
Referenced by initiate(), and solve_red_sys_fin().
long* ndofdom |
numbers of DOFs on subdomains ndofdom[i]=j - the i-th subdomain contains j DOFs
Definition at line 119 of file seqschur.h.
Referenced by assemble_subdom_unknowns(), solve_system(), and subdomain_matrices().
long nicg |
maximum number of iterations in conjugate gradient method
Definition at line 67 of file seqschur.h.
Referenced by print(), read(), and seqschur().
long nidof |
number of internal degrees of freedom (unknowns) on subdomain
Definition at line 81 of file seqschur.h.
Referenced by seqschur().
long ns |
number of subdomains
Definition at line 55 of file seqschur.h.
Referenced by assemble_subdom_unknowns(), initiate(), print(), read(), seqschur(), solve_red_sys_fin(), solve_system(), and subdomain_matrices().
long* nsid |
node-subdomain correspondence nsid[i]=j - the i-th node belongs to the j-th subdomain
Definition at line 114 of file seqschur.h.
storage of reduced system matrix
Definition at line 64 of file seqschur.h.
Referenced by print(), read(), seqschur(), and solve_red_sys_fin().
subdomain matrices stored in the skyline storage
Definition at line 99 of file seqschur.h.
Referenced by solve_system(), and subdomain_matrices().
type of linear system solver
Definition at line 61 of file seqschur.h.
Referenced by print(), read(), seqschur(), solve_red_sys_fin(), and ~seqschur().
type of reduced system solver
Definition at line 58 of file seqschur.h.
Referenced by print(), read(), seqschur(), and solve_red_sys().
double zero |