00001 /************************************************************************ 00002 File aph_std_quad.h - internal information for standard quadratic 00003 approximation module for tetrahedral and prismatic elements 00004 00005 Contains: 00006 - constants 00007 - data types 00008 - global variables (for the whole module) 00009 - headers for internal functions: 00010 apr_select_field - to select the proper field 00011 00012 ------------------------------ 00013 History: 00014 02.2002 - Krzysztof Banas, initial version 00015 10.2013 - Jan Bielański, implementation quadratic approximation without contrained nodes 00016 11.2013 - Jan Bielański, implementation quadratic approximation with contrained nodes 00017 02.2017 - Jan Bielański, implementation support of mixed approximation 00018 *************************************************************************/ 00019 00020 #ifndef _aph_std_quad_ 00021 #define _aph_std_quad_ 00022 00023 00024 #ifdef __cplusplus 00025 extern "C"{ 00026 #endif 00027 00028 /*** CONSTANTS ***/ 00029 // grab defines from aph interface 00030 // KM 06.2011 00031 #include "aph_intf.h" 00032 00033 // +++++++++++++++ Move to aph_intf.h +++++++++++++++++++++ 00035 //#define APC_BASE_TENSOR_DG 1 00036 //#define APC_BASE_COMPLETE_DG 2 00037 //#define APC_BASE_PRISM_STD 3 // for linear prismatic elements 00038 //#define APC_BASE_PRISM_QUAD_STD 32 // for quadratic prismatic elements [element type 3 / approximation degree 2] 00039 //#define APC_BASE_TETRA_STD 4 // for linear tetrahedral elements 00040 //#define APC_BASE_TETRA_QUAD_STD 42 // for quadratic tetrahedral elements [element type 4 / approximation degree 2] 00041 00043 //#define APC_LINEAR_APPROXIMATION_PDEG 1 00044 //#define APC_QUADRATIC_HIERACHICAL_APPROXIMATION_PDEG 15 00045 //#define APC_QUADRATIC_APPROXIMATION_PDEG 2 00046 //#define APC_MIXED_P1_P1_APPROXIMATION_PDEG 11 00047 //#define APC_MIXED_P2_P1_APPROXIMATION_PDEG 21 00048 //#define APC_MIXED_P2_P2_APPROXIMATION_PDEG 22 00049 00051 //#define APC_P1_NONE_G1 101 00052 //#define APC_P15_NONE_G1 151 00053 //#define APC_P2_NONE_G1 201 00054 //#define APC_P2_NONE_G2 202 00055 //#define APC_P1_P1_G1 111 00056 //#define APC_P2_P2_G2 222 00057 //#define APC_P2_P1_G2 212 00058 //#define APC_P2_P2_G2 222 00059 00061 //#define APC_TETRA_NUM_SHAP_LINEAR_APROXIMATION 4 00062 //#define APC_TETRA_NUM_SHAP_QUADRATIC_APROXIMATION 10 00063 //#define APC_PRISM_NUM_SHAP_LINEAR_APROXIMATION 6 00064 //#define APC_PRISM_NUM_SHAP_QUADRATIC_APROXIMATION 18 00065 //#define APC_PRISM_NUM_SHAP_HIERARCHICAL_QUADRATIC_APROXIMATION 15 00066 // +++++++++++++++ Move to aph_intf.h +++++++++++++++++++++ 00067 00068 #define APC_NONE_DOF_NUM_SHAP 0 00069 #define APC_TETRA_VERT_DOF_NUM_SHAP 4 00070 #define APC_TETRA_EDGE_DOF_NUM_SHAP 6 00071 #define APC_PRISM_VERT_DOF_NUM_SHAP 6 00072 #define APC_PRISM_EDGE_DOF_NUM_SHAP 9 00073 #define APC_PRISM_FACE_DOF_NUM_SHAP 3 00074 00075 #define APC_MAX_NUM_FIELD 10 /* maximal number of fields */ 00076 #define APC_MAXELVNO_ON_EDGES 9 /* maximal number of nodes on edges of an element with quadratic approximation */ 00077 #define APC_MAXEQ 12 /* maximal number of solution vector components */ 00078 #define APC_TRUE 1 00079 #define APC_FALSE 0 00080 00081 /*** DATA TYPES ***/ 00082 00083 /* solution degrees of freedom structure (parallel to vertex structure) */ 00084 typedef struct { 00085 double *vec_dof_1; /* pointer to array of field dofs */ 00086 double *vec_dof_2; /* pointer to array of field dofs */ 00087 double *vec_dof_3; /* pointer to array of field dofs */ 00088 00089 int *constr; /* costraints data */ 00090 int *constr_type; /* constraints type */ 00091 00092 int active; /* dof status flag (true/false) - true = has DOFs */ 00093 00094 /* for future extensions (for all mesh entities): 00095 int mesh_ent_type; 00096 int mesh_ent_id; 00097 int pdeg; 00098 */ 00099 } apt_dof_ent; 00100 00101 /* geometric data structure for edges and faces */ 00102 typedef struct { 00103 double *coord; /* coordinates X,Y,Z */ 00104 int active; /* active 0 - no ; 1 - yes */ 00105 } apt_geo_dof_ent; 00106 00107 /* discretization structure, including time and space discretization */ 00108 typedef struct { 00109 int mesh_id; /* identifier of the mesh associated with the field */ 00110 int nreq; /* number of components in solution vector */ 00111 /* APC_MAXELSD = APC_MAXELVD * APC_MAXEQ */ 00112 int mixed_nreq; /* number of components for both fields in mixed approximation */ 00113 int nr_sol; /* number of solution vectors for each element */ 00114 int uniform; /* indicator whether the approximation field is uniform */ 00115 00116 int constr; /* indicator whether there are constrained nodes (vertex or edge or face) */ 00117 int curvilinear; /* has curvilinear elements, very useful when used quadratic geometry; 0 - no ; 1 - yes */ 00118 00119 int nr_dof_vert; /* number of vector vertex DOFs */ 00120 int nr_dof_edge; /* number of vector edge DOFs */ 00121 int nr_dof_face; /* number of vector face DOFs */ 00122 int nr_dof_elem; /* number of vector element DOFs */ 00123 int nr_dof_ents; /* number of vector (vertex + active edge + element) DOFs */ 00124 00125 int nr_geo_dof_edge; /* number of vector geometry edge DOFs */ 00126 int nr_geo_dof_face; /* number of vector geometry face DOFs */ 00127 /* number_of_scalar_DOFs = number_of_vector_DOFs * nreq */ 00128 00129 int capacity_dof_vert; /* max number of vector vertex DOFs (size of DOF arrays) */ 00130 int capacity_dof_edge; /* max number of vector edge DOFs (size of DOF arrays) */ 00131 int capacity_dof_face; /* max number of vector face DOFs (size of DOF arrays) */ 00132 int capacity_dof_elem; /* max number of vector element DOFs (size DOF of arrays) */ 00133 00134 int nr_constr_dof_vert; /* number of constrained vertices */ 00135 int nr_constr_dof_edge; /* number of constrained edges */ 00136 int nr_constr_dof_face; /* number of constrained edges */ 00137 00138 apt_dof_ent *dof_vert; /* pointer to structure with the vectors of vertexes */ 00139 apt_dof_ent *dof_edge; /* pointer to structure with the vectors of edges */ 00140 apt_dof_ent *dof_face; /* pointer to structure with the vectors of faces */ 00141 apt_dof_ent *dof_elem; /* pointer to structure with the vectors of elements */ 00142 00143 int geo_pdeg; /* 1 - (multi-)linear elements, 2 - quadratic geometry */ 00144 /* -1 - second field in mixed approximation */ 00145 apt_geo_dof_ent *geo_dof_edge; /* geometric data for edges */ 00146 apt_geo_dof_ent *geo_dof_face; /* geometric data for faces */ 00147 00148 int pdeg; /* polynomial degree of freedom (for whole space)*/ 00149 int mixed_pdeg; /* pdeg of field which is part of mixed approximation */ 00150 00151 int mixed_second_field; /* information about second field in mixed approximation */ 00152 int mixed_primary_field; /* information about primary field in mixed approximation */ 00153 00154 /* Old entries from apd_std_lin */ 00155 // int constr; /* indicator whether there are constrained nodes */ 00156 // int nr_dof_ents; /* number of dof entities */ 00157 // int capacity_dof_ents; /* max numer of dof ent. with this malloc*/ 00158 // apt_dof_ent *dof_ents; /* pointer to structure with the vectors of dofs */ 00159 /* Old entries from apd_std_lin */ 00160 } apt_field; 00161 00162 /*** GLOBAL VARIABLES for the whole module ***/ 00163 00164 extern int apv_nr_fields; /* the number of fields in the problem */ 00165 extern int apv_cur_field_id; /* ID of the current field */ 00166 extern apt_field apv_fields[APC_MAX_NUM_FIELD]; /* array of fields */ 00167 00168 00169 /*** FUNCTIONS DECLARATIONS - headers for internal functions ***/ 00173 extern apt_field* apr_select_field( 00174 /* returns: pointer to the selected field */ 00175 int Field_id /* in: field ID */ 00176 ); 00177 00178 /* 00179 In file aps_gauss_util.c - utility routines for Gauss-Legendre integration 00180 */ 00181 00185 extern void apr_gauss_init(); 00186 00190 extern void apr_gauss_select( 00191 int Type, /* in: type of integration considered */ 00192 /* 1 - LINE INTEGRATION over [-1,1] */ 00193 /* 2 - GAUSSIAN QUADRATURES FOR TRIANGLES */ 00194 int Order, /* in: order of integration == degree */ 00195 /* of polynomials integrated exactly */ 00196 int *Ng, /* out: number of gauss points */ 00197 double **Xg_p, /* out: gauss points */ 00198 double **Wg_p /* out: gauss weights */ 00199 ); 00200 00201 00202 #ifdef __cplusplus 00203 } 00204 #endif 00205 00206 #endif