#include <schurcompl.h>
Public Member Functions | |
void | gather_bound_vect (double *lv, double *gv, long *domproc) |
void | globcnnum_pdd (gtopology *top, long *ltg, long *domproc, FILE *out) |
void | initiate (partopjk *ptop, FILE *out) |
void | initiate (partop *ptop, selnodes *selnodschur) |
double | pss_gather_bound_vect (double *lv1, double *gv1, double *lv2, double *gv2, long *domproc) |
double | pss_gather_bound_vect_old (double *lv, double *gv, long *domproc) |
void | scatter_bound_vect (double *lv, double *gv, long *domproc) |
schurcompl (int np, int mr, long nd) | |
void | schurordering_new (gtopology *top, long *nodeidentif) |
long | selected_norm_calculation (long cid, double err, double thresh, partopjk *ptop, gtopology *top, double *lv1, double *gv1, double *rhs, double *grhs, long *domproc, long n, FILE *out) |
void | solve_red_sys (long *domproc, double *condmat, double *condvect, long decomp, FILE *out) |
void | solve_red_sys_fin (long *domproc, double *condmat, double *condvect, long decomp, FILE *out) |
void | solve_red_sys_iter (long *domproc, double *condmat, double *condvect, FILE *out) |
void | solve_system (gtopology *top, gmatrix *gm, long *domproc, double *lhs, double *rhs, FILE *out) |
double | unbalanced_values (double *lv, long *domproc, FILE *out) |
~schurcompl () | |
Public Attributes | |
double | aerrcgsch |
long | anicgsch |
gmatrix * | arr |
long | destr_nbdofmas |
storagetype | dmstor |
type of storage of domain matrix | |
double | errcgsch |
gtopology * | gtop |
long | maxnbdof |
int | myrank |
rank of processor | |
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 | ndom |
number of subdomain | |
long | nicgsch |
long | nidof |
number of internal degrees of freedom (unknowns) on subdomain | |
int | nproc |
number of processors | |
storagetype | rsmstor |
storage of reduced system matrix | |
slesolv * | ssle |
type of linear system solver | |
redsystsolver | trssol |
type of reduced system solver | |
double | zero |
class paral contains information about local-global ordering
PRIMAL DOMAIN DECOMPOSITION METHOD
nproc, myrank, ndom number of processors, myrank and number of domain are established after constructor execution
ndof, nidof total number of local unknowns and number of internal unknowns on subdomains are established after execution of function schurordering
ngdof (M) total number of global unknowns is established after execution of function globcnnum_pdd
maxnbn, maxnrdof, totmaxndofn maximum number of boundary nodes, maximum number of boundary unknowns on subdomain and maximum number of degrees of freedom on node are established after execution of globcnnum_pdd
gcn array of local to global code numbers correspondence is established after execution of function globcnnum_pdd
nbndom (M) array of numbers of boundary nodes on subdomains is established after execution of function globcnnum_pdd
nrdofdom (M) array of numbers of reduced unknowns on subdomains is established after execution of function globcnnum_pdd
masgcn (M) array of local to global code numbers of all subdomains is established after execution of function globcnnum_pdd
Definition at line 55 of file schurcompl.h.
schurcompl | ( | int | np, | |
int | mr, | |||
long | nd | |||
) |
constructor
np | - number of processors | |
mr | - my rank (process id) | |
nd | - number of subdomains |
JK, revised 24.7.2008
Definition at line 19 of file schurcompl.cpp.
References aerrcgsch, anicgsch, arr, destr_nbdofmas, errcgsch, gtop, maxnbdof, myrank, nbdof, nbdofmas, ndof, ndom, nicgsch, nidof, nproc, rsmstor, ssle, trssol, and zero.
~schurcompl | ( | ) |
Definition at line 72 of file schurcompl.cpp.
References arr, destr_nbdofmas, gtop, myrank, nbdofmas, and ssle.
void gather_bound_vect | ( | double * | lv, | |
double * | gv, | |||
long * | domproc | |||
) |
function gathers boundary/interface contributions from subdomains into one global vector the global vector must be deleted outside of the subroutine (the global vector contains all components defined on the boundaries/interfaces) (the local vector contains all components defined on a subdomain)
lv | - local vector | |
gv | - global vector, it is allocated only on the master processor the array must be deleted outside of the subroutine | |
domproc | - domain-processor correspondence |
JK, 6.11.2004
Definition at line 1242 of file schurcompl.cpp.
References gtopology::give_code_numbers(), gtopology::give_ndofe(), gtop, locglob(), maxnbdof, MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, myrank, ndof, nidof, nproc, and stat.
Referenced by unbalanced_values().
void globcnnum_pdd | ( | gtopology * | top, | |
long * | ltg, | |||
long * | domproc, | |||
FILE * | out | |||
) |
void initiate | ( | partopjk * | ptop, | |
FILE * | out | |||
) |
Definition at line 392 of file schurcompl.cpp.
References partopjk::bdofd, gtopology::cnstate, gtop, gtopology::initiate(), maxnbdof, partopjk::maxnbdofd, myrank, partopjk::nbdof, nbdof, partopjk::nbdofd, nbdofmas, partopjk::ndof, ndof, partopjk::ndofc, ndofcp, partopjk::nidof, nidof, and nproc.
Function search prescribed common code numbers on boundary nodes in the given subdomain.
top | - pointer to the sequential general topology | |
ltg | - local to global array | |
idgnn | - array with node numbers of boundary nodes | |
nbn | - number of boundary nodes | |
gnnc | - array with boundary node numbers at positions of nodes with prescribed common code numbers | |
id | - searched node id | |
ccn | - searched common code number id function generates ordering of unknowns on subdomains | |
top | - pointer to the sequential general topology | |
ltg | - local to global array |
JK, 7.8.2007 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, 8.7.2005 JK, 7.8.2007
Definition at line 366 of file schurcompl.cpp.
References selnodes::cnmas, gtopology::cnstate, gtop, gtopology::initiate(), maxnbdof, selnodes::maxsnndof, myrank, partop::nbdof, nbdof, nbdofmas, partop::ndof, ndof, ndofcp, partop::nidof, nidof, nproc, selnodes::snndofmas, and selnodes::tndofsn.
double pss_gather_bound_vect | ( | double * | lv1, | |
double * | gv1, | |||
double * | lv2, | |||
double * | gv2, | |||
long * | domproc | |||
) |
function computes dot product of two vectors which are split on several processors it gathers boundary/interface contributions from subdomains into global vectors it gathers results of dot product of internal components of the vectors the global vectors must be deleted outside of the subroutine
lv1,lv2 | - local vectors | |
gv1,gv2 | - global vectors, they are allocated only on the master processor the arrays must be deleted outside of the subroutine |
JK, 6.11.2004, revised 2.6.2011
Definition at line 1307 of file schurcompl.cpp.
References gtopology::give_code_numbers(), gtopology::give_ndofe(), gtop, locglob(), maxnbdof, MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, myrank, ndof, nidof, nproc, and stat.
double pss_gather_bound_vect_old | ( | double * | lv, | |
double * | gv, | |||
long * | domproc | |||
) |
function gathers contributions from subdomains into one global vector the global vector must be deleted outside of the subroutine
lv | - local vector | |
gv | - global vector, it is allocated only on the master processor the array must be deleted outside of the subroutine |
JK, 6.11.2004
Definition at line 1617 of file schurcompl.cpp.
References gtopology::give_code_numbers(), gtopology::give_ndofe(), gtop, locglob(), maxnbdof, MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, myrank, ndof, nidof, nproc, and stat.
void scatter_bound_vect | ( | double * | lv, | |
double * | gv, | |||
long * | domproc | |||
) |
function scatters contributions from one global vector into local vectors function rewrites part of the local vector by contributions from boundaries
lv | - local vector | |
gv | - global vector, it is allocated only on the master processor |
JK, 6.11.2004
Definition at line 1699 of file schurcompl.cpp.
References gtopology::give_code_numbers(), gtopology::give_ndofe(), globloc(), gtop, maxnbdof, MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), myrank, nbdofmas, ndof, nidof, nproc, nullv(), and stat.
Referenced by unbalanced_values().
void schurordering_new | ( | gtopology * | top, | |
long * | nodeidentif | |||
) |
long selected_norm_calculation | ( | long | cid, | |
double | err, | |||
double | thresh, | |||
partopjk * | ptop, | |||
gtopology * | top, | |||
double * | lv1, | |||
double * | gv1, | |||
double * | rhs, | |||
double * | grhs, | |||
long * | domproc, | |||
long | n, | |||
FILE * | out | |||
) |
function computes norm of a vector where only selected components contribute the vector is split on several processors
function assumes that no boundary/interface node is equipped with the Dirichlet boundary condition function assumes that no boundary/interface node is coupled with any other node
cid | - component id (id of required component) | |
ntm | - number of transported media | |
err | - required residual | |
ptop | - pointer to the parallel topology | |
top | - pointer to the subdomain topology | |
lv1,rhs | - local vectors | |
gv1,grhs | - global vectors, they are allocated only on the master processor the arrays must be deleted outside of the subroutine | |
domproc | - domain-processor correspondence |
if the vectors gv1 and grhs are not used outside of this function, they will be removed
JK, 13.7.2011
Definition at line 1405 of file schurcompl.cpp.
References partopjk::bnid, gtopology::give_code_numbers(), gtopology::give_dof(), gtopology::give_ndofe(), gtopology::give_ndofn(), gtop, partopjk::inid, locglob(), maxnbdof, MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_LONG, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, myrank, partopjk::nbn, ndofcp, partopjk::nin, nproc, nullv(), and stat.
void solve_red_sys | ( | long * | domproc, | |
double * | condmat, | |||
double * | condvect, | |||
long | decomp, | |||
FILE * | out | |||
) |
function solves reduced system of linear algebraic equations
domproc | - domain-processor correspondence | |
condmat | - array containing reduced matrices (Schur complements) | |
condvect | - array containing reduced vectors | |
decomp | - indicator for valid condensed matrix of reduced system used for elimination of reduced system only. (0=condensation have to be performed/1=only backward substitution is necessary) | |
out | - output file (for auxiliary output) |
JK
Definition at line 998 of file schurcompl.cpp.
References master_sol, myrank, par_print_err(), paral_sol, solve_red_sys_fin(), solve_red_sys_iter(), and trssol.
Referenced by solve_system().
void solve_red_sys_fin | ( | long * | domproc, | |
double * | condmat, | |||
double * | condvect, | |||
long | decomp, | |||
FILE * | out | |||
) |
function solves reduced system of linear algebraic equations on the master processor by a direct method
domproc | - domain-processor correspondence | |
condmat | - array containing reduced matrices (Schur complements) | |
condvect | - array containing reduced vectors | |
decomp | - indicator for valid condensed matrix of reduced system. (0=condensation have to be performed/1=only backward substitution is necessary) | |
out | - output file (for auxiliary output) |
JK
Definition at line 828 of file schurcompl.cpp.
References gmatrix::alloc(), arr, gtopology::give_code_numbers(), gtopology::give_ndofe(), globloc(), gtop, gmatrix::initiate(), gmatrix::localized(), locglob(), maxnbdof, memset(), MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, myrank, nbdof, nbdofmas, ndofcp, ndom, nproc, nullv(), gmatrix::prepmat(), rsmstor, gmatrix::setval(), slesolv::solve_system(), ssle, stat, and gmatrix::ts.
Referenced by solve_red_sys().
void solve_red_sys_iter | ( | long * | domproc, | |
double * | condmat, | |||
double * | condvect, | |||
FILE * | out | |||
) |
function orders unknowns in coarse problem
function initializates object gtop of the class gtopology it represents subdomains as superelements
ptop | - pointer to parallel topology |
JK, 8.7.2005 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
domproc | - domain-processor correspondence | |
condmat | - array containing reduced matrices (Schur complements) | |
condvect | - array containing reduced vectors | |
out | - output file (for auxiliary output) |
JK, 2.12.2001
Definition at line 599 of file schurcompl.cpp.
References aerrcgsch, anicgsch, copyv(), errcgsch, gtopology::give_code_numbers(), gtopology::give_ndofe(), globloc(), gtop, locglob(), maxnbdof, memset(), MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, mxv(), myrank, nbdof, nbdofmas, ndof, ndofcp, nicgsch, nidof, nproc, nullv(), p, ss(), stat, and zero.
Referenced by solve_red_sys().
void solve_system | ( | gtopology * | top, | |
gmatrix * | gm, | |||
long * | domproc, | |||
double * | lhs, | |||
double * | rhs, | |||
FILE * | out | |||
) |
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 | |
domproc | - domain-processor correspondence | |
lhs | - array containing solution | |
rhs | - array containing right hand side | |
out | - output file (for auxiliary output) |
JK
Definition at line 1029 of file schurcompl.cpp.
References gmatrix::condense(), condmat(), cp, gmatrix::decomp(), maxnbdof, memset(), MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_LONG, MPI_Recv(), MPI_Send(), MPI_Status::MPI_TAG, myrank, nbdof, nproc, solve_red_sys(), stat, and time.
double unbalanced_values | ( | double * | lv, | |
long * | domproc, | |||
FILE * | out | |||
) |
function computes norm of vector of unbalanced values function gathers vector of unbalanced values on boundaries function scatters correct values into subdomains function returns norm of vector of unbalanced values
lv | - local vector | |
domproc | - domain-processor correspondence | |
out | - output file |
JK, 6.11.2004
Definition at line 1755 of file schurcompl.cpp.
References gather_bound_vect(), MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_Barrier(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_Recv(), MPI_Send(), myrank, ndofcp, nidof, nproc, nullv(), scatter_bound_vect(), and stat.
double aerrcgsch |
Definition at line 113 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_iter().
long anicgsch |
Definition at line 109 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_iter().
Definition at line 149 of file schurcompl.h.
Referenced by schurcompl(), solve_red_sys_fin(), and ~schurcompl().
long destr_nbdofmas |
Definition at line 146 of file schurcompl.h.
Referenced by schurcompl(), and ~schurcompl().
type of storage of domain matrix
Definition at line 138 of file schurcompl.h.
double errcgsch |
Definition at line 111 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_iter().
Definition at line 148 of file schurcompl.h.
Referenced by gather_bound_vect(), initiate(), pss_gather_bound_vect(), pss_gather_bound_vect_old(), scatter_bound_vect(), schurcompl(), selected_norm_calculation(), solve_red_sys_fin(), solve_red_sys_iter(), and ~schurcompl().
long maxnbdof |
Definition at line 104 of file schurcompl.h.
Referenced by gather_bound_vect(), initiate(), pss_gather_bound_vect(), pss_gather_bound_vect_old(), scatter_bound_vect(), schurcompl(), selected_norm_calculation(), solve_red_sys_fin(), solve_red_sys_iter(), and solve_system().
int myrank |
rank of processor
Definition at line 86 of file schurcompl.h.
Referenced by gather_bound_vect(), initiate(), pss_gather_bound_vect(), pss_gather_bound_vect_old(), scatter_bound_vect(), schurcompl(), selected_norm_calculation(), solve_red_sys(), solve_red_sys_fin(), solve_red_sys_iter(), solve_system(), unbalanced_values(), and ~schurcompl().
long nbdof |
number of boundary degrees of freedom (unknowns) on subdomain
Definition at line 95 of file schurcompl.h.
Referenced by initiate(), schurcompl(), solve_red_sys_fin(), solve_red_sys_iter(), and solve_system().
long* nbdofmas |
Definition at line 144 of file schurcompl.h.
Referenced by initiate(), scatter_bound_vect(), schurcompl(), solve_red_sys_fin(), solve_red_sys_iter(), and ~schurcompl().
long ndof |
number of degrees of freedom on subdomain
Definition at line 91 of file schurcompl.h.
Referenced by gather_bound_vect(), initiate(), pss_gather_bound_vect(), pss_gather_bound_vect_old(), scatter_bound_vect(), schurcompl(), and solve_red_sys_iter().
long ndofcp |
number of global DOFs = total number of boundary DOFs number of DOFs of coarse problem
Definition at line 126 of file schurcompl.h.
Referenced by initiate(), selected_norm_calculation(), solve_red_sys_fin(), solve_red_sys_iter(), and unbalanced_values().
long ndom |
number of subdomain
Definition at line 88 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_fin().
long nicgsch |
Definition at line 107 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_iter().
long nidof |
number of internal degrees of freedom (unknowns) on subdomain
Definition at line 93 of file schurcompl.h.
Referenced by gather_bound_vect(), initiate(), pss_gather_bound_vect(), pss_gather_bound_vect_old(), scatter_bound_vect(), schurcompl(), solve_red_sys_iter(), and unbalanced_values().
int nproc |
number of processors
Definition at line 84 of file schurcompl.h.
Referenced by gather_bound_vect(), initiate(), pss_gather_bound_vect(), pss_gather_bound_vect_old(), scatter_bound_vect(), schurcompl(), selected_norm_calculation(), solve_red_sys_fin(), solve_red_sys_iter(), solve_system(), and unbalanced_values().
storage of reduced system matrix
Definition at line 132 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_fin().
type of linear system solver
Definition at line 135 of file schurcompl.h.
Referenced by schurcompl(), solve_red_sys_fin(), and ~schurcompl().
type of reduced system solver
Definition at line 129 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys().
double zero |
Definition at line 115 of file schurcompl.h.
Referenced by schurcompl(), and solve_red_sys_iter().