Problem Definition

Collaboration diagram for Problem Definition:

Modules

 Problem Control

Defines

#define PDC_PI   3.141592654
#define PDC_SQRT_PI   1.772453851
#define PDC_SQRT_3   1.732050808
#define PDC_MAXEQ   30
#define PDC_MAXBCVAL   10
#define PDC_INTERIOR   -1
#define PDC_BC_DIRI   1
#define PDC_BC_NEUM   2
#define PDC_BC_MIXED   3
#define PDC_USE_CURRENT_PROBLEM_ID   -1
#define Current_solution_ID   1
#define Previous_iteration_sol_ID   2
#define Previous_time_step_sol_ID   3
#define Geometry_dofs_ID   999
#define PDC_UNIFORM_DEREF   -2
#define PDC_UNIFORM_REF   -1
#define PDC_NO_ADAPT   0
#define PDC_ADAPT_EXACT   1
#define PDC_ADAPT_ZZ   2
#define PDC_ADAPT_EXPL   3
#define PDC_MAX_DOF_PER_INT   70

Functions

int pdr_module_introduce (char *Problem_module_name)
int pdr_get_list_ent (int Problem_id, int *Nr_int_ent, int **List_int_ent_type, int **List_int_ent_id, int *Nr_dof_ent, int **List_dof_ent_type, int **List_dof_ent_id, int **List_dof_ent_nrdofs, int *Nrdofs_glob, int *Max_dofs_per_dof_ent)
int pdr_get_list_ent_coarse (int Problem_id, int Nr_int_ent_fine, int *List_int_ent_type_fine, int *List_int_ent_id_fine, int Nr_dof_ent_fine, int *List_dof_ent_type_fine, int *List_dof_ent_id_fine, int *List_dof_ent_nrdofs_fine, int Nrdofs_glob_fine, int Max_dof_per_ent_fine, int *Pdeg_coarse, int *Nr_int_ent, int **List_int_ent_type, int **List_int_ent_id, int *Nr_dof_ent, int **List_dof_ent_type, int **List_dof_ent_id, int **List_dof_ent_nrdofs, int *Nrdofs_glob, int *Max_dof_per_ent)
int pdr_get_max_num_grid_levels (int Problem_id)
int pdr_create_assemble_stiff_mat (int Problem_id, int Level_id, int Comp_type, int *Pdeg_coarse_p, int Nr_int_ent, int *L_int_ent_type, int *L_int_ent_id, int Nr_colors_elems, int *L_color_index_elems, int Nr_colors_faces, int *L_color_index_faces, int *Asse_pos_first_dof_int_ent, int *Assembly_table, int *Pos_first_dof_int_ent, int *Local_to_global, int Max_dofs_int_ent)
int pdr_assemble_local_stiff_mat_with_table (int Problem_id, int Level_id, int Comp_type, int Nr_dof_ent, int *Assembly_table_int_ent, int *Local_to_global_int_ent, double *Stiff_mat, double *Rhs_vect, char *Rewr_dofs)
int pdr_assemble_local_stiff_mat (int Problem_id, int Level_id, int Comp_type, int Nr_dof_ent, int *L_dof_ent_type, int *L_dof_ent_id, int *L_dof_ent_nrdofs, double *Stiff_mat, double *Rhs_vect, char *Rewr_dofs)
int pdr_comp_stiff_mat (int Problem_id, int Int_ent_type, int Int_ent_id, int Comp_sm, int *Pdeg, int *Nr_dof_ent, int *List_dof_ent_type, int *List_dof_ent_id, int *List_dof_ent_nrdofs, int *Nrdofs_loc, double *Stiff_mat, double *Rhs_vect, char *Rewr_dofs)
int pdr_comp_stiff_mat_uncon (int Problem_id, int Int_ent_type, int Int_ent_id, int Comp_sm, int *Pdeg, int *Nrdofs_loc, double *Stiff_mat, double *Rhs_vect, char *Rewr_dofs)
double pdr_err_indi (int Problem_id, int Mode, int El)
int pdr_read_sol_dofs (int Problem_id, int Dof_ent_type, int Dof_ent_id, int Nrdofs, double *Vect_dofs)
int pdr_write_sol_dofs (int Problem_id, int Dof_ent_type, int Dof_ent_id, int Nrdofs, double *Vect_dofs)
int pdr_L2_proj_sol (int Problem_id, int El, int *Pdeg_vec, double *Dofs, int *El_from, int *Pdeg_vec_from, double *Dofs_from)
int pdr_renum_coeff (int Problem_id, int Ent_type, int Ent_id, double *Ren_coeff)
int pdr_get_ent_pdeg (int Problem_id, int Ent_type, int Ent_id)
int pdr_dof_ent_sons (int Problem_id, int Ent_type, int Ent_id, int *Ent_sons)
int pdr_proj_sol_lev (int Problem_id, int Solver_id, int Ilev_from, double *Vec_from, int Ilev_to, double *Vec_to)
double pdr_vec_norm (int Problem_id, int Solver_id, int Level_id, int Nrdof, double *Vector)
double pdr_sc_prod (int Problem_id, int Solver_id, int Level_id, int Nrdof, double *Vector1, double *Vector2)
int pdr_create_exchange_tables (int Problem_id, int Solver_id, int Level_id, int Nr_dof_ent, int *L_dof_ent_type, int *L_dof_ent_id, int *L_bl_nrdof, int *L_bl_posg, int *L_elem_to_bl, int *L_face_to_bl, int *L_edge_to_bl, int *L_vert_to_bl)
int pdr_exchange_dofs (int Problem_id, int Solver_id, int Level_id, double *Vec_dofs)
int pdr_select_el_coeff_vect (int Problem_id, int *Coeff_vect_ind)
double * pdr_select_el_coeff (int Problem_id, double **Mval, double **Axx, double **Axy, double **Axz, double **Ayx, double **Ayy, double **Ayz, double **Azx, double **Azy, double **Azz, double **Bx, double **By, double **Bz, double **Tx, double **Ty, double **Tz, double **Cval, double **Lval, double **Qx, double **Qy, double **Qz, double **Sval)
int pdr_el_coeff (int Problem_id, int Elem, int Mat_num, double Hsize, int Pdeg, double *X_loc, double *Base_phi, double *Base_dphix, double *Base_dphiy, double *Base_dphiz, double *Xcoor, double *Uk_val, double *Uk_x, double *Uk_y, double *Uk_z, double *Un_val, double *Un_x, double *Un_y, double *Un_z, double *Mval, double *Axx, double *Axy, double *Axz, double *Ayx, double *Ayy, double *Ayz, double *Azx, double *Azy, double *Azz, double *Bx, double *By, double *Bz, double *Tx, double *Ty, double *Tz, double *Cval, double *Lval, double *Qx, double *Qy, double *Qz, double *Sval)

Variables

const int PDC_ELEMENT
const int PDC_FACE
const int PDC_EDGE
const int PDC_VERTEX
const int PDC_MIXED_ELEMENT
const int PDC_MIXED_FACE
const int PDC_MIXED_EDGE
const int PDC_MIXED_VERTEX
const int PDC_NO_COMP
const int PDC_COMP_SM
const int PDC_COMP_RHS
const int PDC_COMP_BOTH
const int PDC_COMP_MM
utt_patchespdv_patches

Define Documentation

#define Current_solution_ID   1

compute entries for sm and rhsv

#define Geometry_dofs_ID   999
#define PDC_ADAPT_EXACT   1
#define PDC_ADAPT_EXPL   3
#define PDC_ADAPT_ZZ   2
#define PDC_BC_DIRI   1
#define PDC_BC_MIXED   3
#define PDC_BC_NEUM   2
#define PDC_INTERIOR   -1
#define PDC_MAX_DOF_PER_INT   70
#define PDC_MAXBCVAL   10
#define PDC_MAXEQ   30
#define PDC_NO_ADAPT   0
#define PDC_PI   3.141592654
#define PDC_SQRT_3   1.732050808
#define PDC_SQRT_PI   1.772453851
#define PDC_UNIFORM_DEREF   -2
#define PDC_UNIFORM_REF   -1
#define PDC_USE_CURRENT_PROBLEM_ID   -1
#define Previous_iteration_sol_ID   2
#define Previous_time_step_sol_ID   3

Function Documentation

int pdr_assemble_local_stiff_mat ( int  Problem_id,
int  Level_id,
int  Comp_type,
int  Nr_dof_ent,
int *  L_dof_ent_type,
int *  L_dof_ent_id,
int *  L_dof_ent_nrdofs,
double *  Stiff_mat,
double *  Rhs_vect,
char *  Rewr_dofs 
)

----------------------------------------------------------- pdr_assemble_local_stiff_mat - to assemble an element stiffness matrix to the global SM ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Level_id in: solver ID (used to identify the subproblem)
Comp_type in: level ID
Rhs_vect in: stiffness matrix stored columnwise
Rewr_dofs in: rhs vector in: flag to rewrite or sum up entries 'T' - true, rewrite entries when assembling 'F' - false, sum up entries when assembling

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_assemble_local_stiff_mat_with_table ( int  Problem_id,
int  Level_id,
int  Comp_type,
int  Nr_dof_ent,
int *  Assembly_table_int_ent,
int *  Local_to_global_int_ent,
double *  Stiff_mat,
double *  Rhs_vect,
char *  Rewr_dofs 
)

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_comp_stiff_mat ( int  Problem_id,
int  Int_ent_type,
int  Int_ent_id,
int  Comp_sm,
int *  Pdeg,
int *  Nr_dof_ent,
int *  List_dof_ent_type,
int *  List_dof_ent_id,
int *  List_dof_ent_nrdofs,
int *  Nrdofs_loc,
double *  Stiff_mat,
double *  Rhs_vect,
char *  Rewr_dofs 
)

----------------------------------------------------------- pdr_comp_stiff_mat - to construct a stiffness matrix and a load vector for some given mesh entity ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Int_ent_type in: approximation field ID
Int_ent_id in: type of the integration entity
Comp_sm in: unique identifier of the integration entity
Pdeg in: indicator for the scope of computations: PDC_NO_COMP - do not compute anything PDC_COMP_SM - compute entries to stiff matrix only PDC_COMP_RHS - compute entries to rhs vector only PDC_COMP_BOTH - compute entries for sm and rhsv
Nr_dof_ent in: enforced degree of polynomial (if != NULL )
List_dof_ent_type in: size of List_dof_ent_... arrays, out: no of filled entries, i.e. number of mesh entities with which dofs and stiffness matrix blocks are associated
List_dof_ent_id out: list of types for 'dof' entity
List_dof_ent_nrdofs out: list of IDs for 'dof' entity
Nrdofs_loc out: list of no of dofs for 'dof' entity
Stiff_mat in(optional): size of Stiff_mat and Rhs_vect out(optional): actual number of dofs per integration entity
Rhs_vect out(optional): stiffness matrix stored columnwise
Rewr_dofs out(optional): rhs vector out(optional): flag to rewrite or sum up entries 'T' - true, rewrite entries when assembling 'F' - false, sum up entries when assembling

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_comp_stiff_mat_uncon ( int  Problem_id,
int  Int_ent_type,
int  Int_ent_id,
int  Comp_sm,
int *  Pdeg_vec,
int *  Nrdofs_loc,
double *  Stiff_mat,
double *  Rhs_vect,
char *  Rewr_dofs 
)

----------------------------------------------------------- pdr_comp_stiff_mat_uncon - to construct UNCONSTRAINED stiffness matrix and a load vector for some given mesh entity ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Int_ent_type in: approximation field ID
Int_ent_id in: type of the integration entity
Comp_sm in: unique identifier of the integration entity
Pdeg_vec in: indicator for the scope of computations: PDC_NO_COMP - do not compute anything PDC_COMP_SM - compute entries to stiff matrix only PDC_COMP_RHS - compute entries to rhs vector only PDC_COMP_BOTH - compute entries for sm and rhsv
Nrdofs_loc in: enforced degree of polynomial (if != NULL )
Stiff_mat in(optional): size of Stiff_mat and Rhs_vect out(optional): actual number of dofs per integration entity
Rhs_vect out(optional): stiffness matrix stored columnwise
Rewr_dofs out(optional): rhs vector out(optional): flag to rewrite or sum up entries 'T' - true, rewrite entries when assembling 'F' - false, sum up entries when assembling

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_create_assemble_stiff_mat ( int  Problem_id,
int  Level_id,
int  Comp_type,
int *  Pdeg_coarse_p,
int  Nr_int_ent,
int *  L_int_ent_type,
int *  L_int_ent_id,
int  Nr_colors_elems,
int *  L_color_index_elems,
int  Nr_colors_faces,
int *  L_color_index_faces,
int *  Asse_pos_first_dof_int_ent,
int *  Assembly_table,
int *  Pos_first_dof_int_ent,
int *  Local_to_global,
int  Max_dofs_int_ent 
)

----------------------------------------------------------- pdr_create_assemble_stiff_mat - to create element stiffness matrices and assemble them to the global SM ------------------------------------------------------------

Parameters:
Pdeg_coarse_p in: indicator for the scope of computations:

Here is the call graph for this function:

int pdr_create_exchange_tables ( int  Problem_id,
int  Solver_id,
int  Level_id,
int  Nr_dof_ent,
int *  L_dof_ent_type,
int *  L_dof_ent_id,
int *  L_bl_nrdof,
int *  L_bl_posg,
int *  L_elem_to_bl,
int *  L_face_to_bl,
int *  L_edge_to_bl,
int *  L_vert_to_bl 
)

-------------------------------------------------------- pdr_create_exchange_tables - to create tables to exchange dofs ---------------------------------------------------------

Parameters:
Problem_id returns: >=0 -success code, <0 -error code
Solver_id in: problem ID
Level_id in: solver data structure to be used
Nr_dof_ent in: level ID
L_dof_ent_type in: number of DOF entities in the level all four subsequent arrays are indexed by block IDs with 1(!!!) offset
L_dof_ent_id in: list of DOF entities associated with DOF blocks
L_bl_nrdof in: list of DOF entities associated with DOF blocks
L_bl_posg in: list of nrdofs for each dof block
L_elem_to_bl in: list of positions within the global vector of dofs for each dof block
L_face_to_bl in: list of DOF blocks associated with DOF entities
L_edge_to_bl in: list of DOF blocks associated with DOF entities
L_vert_to_bl in: list of DOF blocks associated with DOF entities in: list of DOF blocks associated with DOF entities

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_dof_ent_sons ( int  Problem_id,
int  Ent_type,
int  Ent_id,
int *  Ent_sons 
)

-------------------------------------------------------- pdr_dof_ent_sons - to return a list of dof entity sons ---------------------------------------------------------

Parameters:
Problem_id returns: success >=0 or <0 - error code
Ent_type in: problem ID
Ent_id in: type of mesh entity
Ent_sons in: mesh entity ID out: list of dof entity sons Ent_sons[0] - number of sons

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_el_coeff ( int  Problem_id,
int  Elem,
int  Mat_num,
double  Hsize,
int  Pdeg,
double *  X_loc,
double *  Base_phi,
double *  Base_dphix,
double *  Base_dphiy,
double *  Base_dphiz,
double *  Xcoor,
double *  Uk_val,
double *  Uk_x,
double *  Uk_y,
double *  Uk_z,
double *  Un_val,
double *  Un_x,
double *  Un_y,
double *  Un_z,
double *  Mval,
double *  Axx,
double *  Axy,
double *  Axz,
double *  Ayx,
double *  Ayy,
double *  Ayz,
double *  Azx,
double *  Azy,
double *  Azz,
double *  Bx,
double *  By,
double *  Bz,
double *  Tx,
double *  Ty,
double *  Tz,
double *  Cval,
double *  Lval,
double *  Qx,
double *  Qy,
double *  Qz,
double *  Sval 
)

!!!!!!!!!! THE FUNCTION BELOW IS USED BY APPROXIMATION MODULE !!!!!!!!!!! !!! DO NOT CHANGE ITS HEADER (OR CHANGE APPROXIMATION MODULE AS WELL) !!! -------------------------------------------------------- pdr_el_coeff - to return coefficients for element integrals ----------------------------------------------------------

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

---- CALCULATE FRICTIONAL FORCE IN THE MUSHY ZONE AT GAUSS POINT -----!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

---------------------- CALCULATE TOTAL SOURCE TERMS ------------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

-------------------- MATERIAL DATA AT GAUSS POINT ----------------- --!

----------------------------------------------------------------------!

----------------------------------------------------------------------!

------------ CALCULATE BOUSINESSQUE FORCE AT GAUSS POINT -------------!

----------------------------------------------------------------------!

Parameters:
Mat_num in: element number
Hsize in: material number
Pdeg in: size of an element
X_loc in: local degree of polynomial
Base_phi in: local coordinates of point within element
Base_dphix in: basis functions
Base_dphiy in: x-derivatives of basis functions
Base_dphiz in: y-derivatives of basis functions
Xcoor in: z-derivatives of basis functions
Uk_val in: global coordinates of a point
Uk_x in: computed solution from previous iteration
Uk_y in: gradient of computed solution Uk_val
Uk_z in: gradient of computed solution Uk_val
Un_val in: gradient of computed solution Uk_val
Un_x in: computed solution from previous time step
Un_y in: gradient of computed solution Un_val
Un_z in: gradient of computed solution Un_val
Mval in: gradient of computed solution Un_val
Axx out: mass matrix coefficient
Ayx out:diffusion coefficients
Azx e.g. Axy denotes scalar or matrix
Bx related to terms with dv/dx*du/dy second order derivatives in weak formulation (scalar for scalar problems matrix for vector problems) WARNING: if axy==NULL only diagonal (axx, ayy, azz) terms are considered in apr_num_int_el OPTIONS: azz!=NULL, axy!=NULL - all a.. matrices must be specified azz!=NULL, axy==NULL - axx, ayy, azz matrices must be specified azz==NULL - axx, axy, ayx, ayy matrices must be specified
Tx out: convection coefficients Bx denotes scalar or matrix related to terms with du/dx*v in weak form
Cval out: convection coefficients Tx denotes scalar or matrix related to terms with u*dv/dx in weak form
Lval out: reaction coefficients - for terms without derivatives in weak form (as usual: scalar for scalar problems, matrix for vectors)
Qx out: rhs coefficient for time term, Lval denotes scalar or matrix corresponding to time derivative - similar as mass matrix but with known solution at the previous time step (usually denoted by u_n)
Qy out: rhs coefficients for terms with derivatives
Qz Qy denotes scalar or matrix corresponding to terms with dv/dy
Sval derivatives in weak formulation out: rhs coefficients without derivatives (source terms)

Here is the call graph for this function:

double pdr_err_indi ( int  Problem_id,
int  Mode,
int  El 
)

----------------------------------------------------------- pdr_err_indi - to compute an error indicator for an element ------------------------------------------------------------

Parameters:
Problem_id returns error indicator for an element
Mode in: data structure to be used
El in: mode of operation in: element number

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_exchange_dofs ( int  Problem_id,
int  Solver_id,
int  Level_id,
double *  Vec_dofs 
)

-------------------------------------------------------- pdr_exchange_dofs - to exchange dofs between processors ---------------------------------------------------------

Parameters:
Solver_id in: problem ID
Level_id in: solver data structure to be used
Vec_dofs in: level number in: vector of dofs to be exchanged

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_get_ent_pdeg ( int  Problem_id,
int  Ent_type,
int  Ent_id 
)

----------------------------------------------------------- pdr_get_ent_pdeg - to return the degree of approximation index associated with a given mesh entity ------------------------------------------------------------

Parameters:
Problem_id returns: >0 - approximation index, <0 - error code
Ent_type in: approximation field ID
Ent_id in: type of mesh entity in: mesh entity ID

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_get_list_ent ( int  Problem_id,
int *  Nr_int_ent,
int **  List_int_ent_type,
int **  List_int_ent_id,
int *  Nr_dof_ent,
int **  List_dof_ent_type,
int **  List_dof_ent_id,
int **  List_dof_ent_nrdofs,
int *  Nrdofs_glob,
int *  Max_dofs_per_dof_ent 
)

----------------------------------------------------------- pdr_get_list_ent - to return the list of integration entities - entities for which stiffness matrices and load vectors are provided by the FEM code to the solver module ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Nr_int_ent in: problem (and solver) identification
List_int_ent_type out: number of integration entitites GHOST ENTITIES HAVE NEGATIVE TYPE !!!
List_int_ent_id out: list of types of integration entitites
Nr_dof_ent out: list of IDs of integration entitites
List_dof_ent_type out: number of dof entities (entities with which there are dofs associated by the given approximation)
List_dof_ent_id out: list of types of integration entitites
List_dof_ent_nrdofs out: list of IDs of integration entitites
Nrdofs_glob out: list of no of dofs for 'dof' entity
Max_dofs_per_dof_ent out: global number of degrees of freedom (unknowns) out: maximal number of dofs per dof entity

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_get_list_ent_coarse ( int  Problem_id,
int  Nr_int_ent_fine,
int *  List_int_ent_type_fine,
int *  List_int_ent_id_fine,
int  Nr_dof_ent_fine,
int *  List_dof_ent_type_fine,
int *  List_dof_ent_id_fine,
int *  List_dof_ent_nrdofs_fine,
int  Nrdofs_glob_fine,
int  Max_dof_per_ent_fine,
int *  Pdeg_coarse,
int *  Nr_int_ent,
int **  List_int_ent_type,
int **  List_int_ent_id,
int *  Nr_dof_ent,
int **  List_dof_ent_type,
int **  List_dof_ent_id,
int **  List_dof_ent_nrdofs,
int *  Nrdofs_glob,
int *  Max_dof_per_ent 
)

----------------------------------------------------------- pdr_get_list_ent_coarse - to return the list of integration entities - entities for which stiffness matrices and load vectors are provided by the FEM code to the solver module, and DOF entities - entities with which there are dofs associated by the given approximation for COARSE level given the corresponding lists from the fine level ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Nr_int_ent_fine in: problem (and solver) identification
List_int_ent_type_fine in: number of integration entitites
List_int_ent_id_fine in: list of types of integration entitites
Nr_dof_ent_fine in: list of IDs of integration entitites
List_dof_ent_type_fine in: number of dof entities (entities with which there are dofs associated by the given approximation)
List_dof_ent_id_fine in: list of types of integration entitites
List_dof_ent_nrdofs_fine in: list of IDs of integration entitites
Nrdofs_glob_fine in: list of no of dofs for 'dof' entity
Max_dof_per_ent_fine in: global number of degrees of freedom (unknowns)
Pdeg_coarse in: maximal number of dofs per dof entity
Nr_int_ent in: degree of approximation for coarse space
List_int_ent_type out: number of integration entitites
List_int_ent_id out: list of types of integration entitites
Nr_dof_ent out: list of IDs of integration entitites
List_dof_ent_type out: number of dof entities (entities with which there are dofs associated by the given approximation)
List_dof_ent_id out: list of types of integration entitites
List_dof_ent_nrdofs out: list of IDs of integration entitites
Nrdofs_glob out: list of no of dofs for 'dof' entity
Max_dof_per_ent out: global number of degrees of freedom (unknowns) out: maximal number of dofs per dof entity

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_get_max_num_grid_levels ( int  Problem_id  ) 

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_L2_proj_sol ( int  Problem_id,
int  El,
int *  Pdeg_vec,
double *  Dofs,
int *  El_from,
int *  Pdeg_vec_from,
double *  Dofs_from 
)

-------------------------------------------------------- pdr_L2_proj_sol - to project solution between elements of different generations ----------------------------------------------------------

Parameters:
El in: problem ID
Pdeg_vec in: element number
Dofs in: element degree of approximation
El_from out: workspace for degress of freedom of El NULL - write to data structure
Pdeg_vec_from in: list of elements to provide function
Dofs_from in: degree of polynomial for each El_from if element Pdeg is a vector - its components must be suitably placed in Pdeg_vec_from, according to conventions in approximation module in: Dofs of El_from or...

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_module_introduce ( char *  Problem_module_name  ) 

Here is the caller graph for this function:

int pdr_proj_sol_lev ( int  Problem_id,
int  Solver_id,
int  Ilev_from,
double *  Vec_from,
int  Ilev_to,
double *  Vec_to 
)

-------------------------------------------------------- pdr_proj_sol_lev - to L2 project solution dofs between mesh levels ---------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success; <0 - error code
Solver_id in: problem ID
Ilev_from in: solver data structure to be used
Vec_from in: level number to project from
Ilev_to in: vector of values to project
Vec_to in: level number to project to out: vector of projected values
int pdr_read_sol_dofs ( int  Problem_id,
int  Dof_ent_type,
int  Dof_ent_id,
int  Nrdofs,
double *  Vect_dofs 
)

----------------------------------------------------------- pdr_read_sol_dofs - to read a vector of dofs associated with a given mesh entity from approximation field data structure ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Dof_ent_type in: solver ID (used to identify the subproblem)
Vect_dofs in: dofs to be written

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_renum_coeff ( int  Problem_id,
int  Ent_type,
int  Ent_id,
double *  Ren_coeff 
)

-------------------------------------------------------- pdr_renum_coeff - to return a coefficient being a basis for renumbering ----------------------------------------------------------

Parameters:
Ent_type in: problem ID
Ent_id in: type of mesh entity
Ren_coeff in: mesh entity ID out: renumbering coefficient
double pdr_sc_prod ( int  Problem_id,
int  Solver_id,
int  Level_id,
int  Nrdof,
double *  Vector1,
double *  Vector2 
)

-------------------------------------------------------- pdr_sc_prod - to compute a scalar product of two global vectors ---------------------------------------------------------

Parameters:
Problem_id retruns: scalar product of Vector1 and Vector2
Solver_id in: problem ID
Level_id in: solver data structure to be used
Nrdof in: level number
Vector1 in: number of vector components
Vector2 in: local part of global Vector in: local part of global Vector

Here is the call graph for this function:

Here is the caller graph for this function:

double* pdr_select_el_coeff ( int  Problem_id,
double **  Mval,
double **  Axx,
double **  Axy,
double **  Axz,
double **  Ayx,
double **  Ayy,
double **  Ayz,
double **  Azx,
double **  Azy,
double **  Azz,
double **  Bx,
double **  By,
double **  Bz,
double **  Tx,
double **  Ty,
double **  Tz,
double **  Cval,
double **  Lval,
double **  Qx,
double **  Qy,
double **  Qz,
double **  Sval 
)

!!!!! OLD OBSOLETE VERSION !!!!!! ----------------------------------------------------------- pdr_select_el_coeff - to select coefficients returned to approximation routines for element integrals ------------------------------------------------------------

!!!!! OLD OBSOLETE VERSION !!!!!!

Parameters:
Problem_id returns: pointer !=NULL to indicate selection
Axx out: mass matrix coefficient
Ayx out:diffusion coefficients, e.g.
Azx Axy denotes scalar or matrix
Bx related to terms with dv/dx*du/dy second order derivatives in weak formulation (scalar for scalar problems matrix for vector problems) WARNING: if axy==NULL only diagonal (axx, ayy, azz) terms are considered in apr_num_int_el OPTIONS: azz!=NULL, axy!=NULL - all a.. matrices must be specified azz!=NULL, axy==NULL - axx, ayy, azz matrices must be specified azz==NULL - axx, axy, ayx, ayy matrices must be specified
Tx out: convection coefficients Bx denotes scalar or matrix related to terms with du/dx*v in weak form
Cval out: convection coefficients Tx denotes scalar or matrix related to terms with u*dv/dx in weak form
Lval out: reaction coefficients - for terms without derivatives in weak form (as usual: scalar for scalar problems, matrix for vectors)
Qx out: rhs coefficient for time term, Lval denotes scalar or matrix corresponding to time derivative - similar as mass matrix but with known solution at the previous time step (usually denoted by u_n)
Qy out: rhs coefficients for terms with derivatives
Qz Qy denotes scalar or matrix corresponding to terms with dv/dy
Sval derivatives in weak formulation out: rhs coefficients without derivatives (source terms)

Here is the call graph for this function:

int pdr_select_el_coeff_vect ( int  Problem_id,
int *  Coeff_vect_ind 
)

!!!!!!!!!! THE FUNCTION BELOW IS USED BY APPROXIMATION MODULE !!!!!!!!!!! !!! DO NOT CHANGE ITS HEADER (OR CHANGE APPROXIMATION MODULE AS WELL) !!! ----------------------------------------------------------- pdr_select_el_coeff_vect - to select coefficients returned to approximation routines for element integrals in weak formulation (the procedure indicates which terms are non-zero in weak form) ------------------------------------------------------------

Parameters:
Coeff_vect_ind out: coefficient indicator

Here is the call graph for this function:

double pdr_vec_norm ( int  Problem_id,
int  Solver_id,
int  Level_id,
int  Nrdof,
double *  Vector 
)

-------------------------------------------------------- pdr_vec_norm - to compute a norm of global vector (in parallel) ---------------------------------------------------------

Parameters:
Problem_id returns: L2 norm of global Vector
Solver_id in: problem ID
Level_id in: solver data structure to be used
Nrdof in: level number
Vector in: number of vector components in: local part of global Vector

Here is the call graph for this function:

Here is the caller graph for this function:

int pdr_write_sol_dofs ( int  Problem_id,
int  Dof_ent_type,
int  Dof_ent_id,
int  Nrdofs,
double *  Vect_dofs 
)

----------------------------------------------------------- pdr_write_sol_dofs - to write a vector of dofs associated with a given mesh entity to approximation field data structure ------------------------------------------------------------

Parameters:
Problem_id returns: >=0 - success code, <0 - error code
Dof_ent_type in: solver ID (used to identify the subproblem)
Vect_dofs in: dofs to be written

Here is the call graph for this function:

Here is the caller graph for this function:


Variable Documentation

const int PDC_COMP_BOTH

compute entries to rhs vector only

const int PDC_COMP_MM

compute entries for sm and rhsv

const int PDC_COMP_RHS

compute entries to stiff matrix only

const int PDC_COMP_SM

do not compute stiff matrix and rhs vector

const int PDC_EDGE
const int PDC_ELEMENT

types of integration and dof entities

const int PDC_FACE
const int PDC_MIXED_EDGE
const int PDC_MIXED_ELEMENT
const int PDC_MIXED_FACE
const int PDC_MIXED_VERTEX
const int PDC_NO_COMP

options for create_stiff_mat procedure

const int PDC_VERTEX
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Generated on 3 Jun 2020 for ModFEM by  doxygen 1.6.1