AdH Kraken 0.0.0
Next generation Adaptive Hydraulics
Loading...
Searching...
No Matches
selem_2d.h
1#ifndef H_SELEM_2D_
2#define H_SELEM_2D_
3
4// dependencies :: SVECT, SVECT2D
5
6/***********************************************************/
7/***********************************************************/
8/***********************************************************/
9
10typedef struct {
11
12 int id; /* current 2d element number */
13 int gid; /*global element id */
14 int id_orig; /* original 2d element number (every 2d element should have a value) */
15 int id_3d; /* the 3d element that has a face coincident with this 2d element */
16 int nnodes; /* total # of nodes on element */
17 int nnodes_quad; /* the total # of quadratic nodes on element */
18 double djac; /* jacobian given only 2D coordinates */
19 double djac3d; /* the 2D jacobian surface in 3D */
20 double djac3d_fixed; /* the area of the original element in 3D */
21 SVECT nrml; /* the normal to the face */
22 int string; /* string # associated with this 2d element*/
23 int mat; /* 2d element physics material type */
24 int nvars; /* the total number of independent variables on this element */
25 int *vars; /* an array of independent variables on this element */
26 int bflag; /* element boundary flag :: -1 = body, 0 = surface, 1 = bottom, 2 = sidewall (3d) */
27 SVECT2D *grad_shp; /* the gradients of the shape functions */
28 int *nodes; /* the nodes in the element */
29 int *levels; /* the node levels in the element */
30 int nedges; /* the number of edges on the element */
31 double area;
32 int resident_pe; /* arbitrary since AdH decomp is node-based. Defines as lowest node PE */
33
34 /*interface flag to prevent adaption */
35 int interface;
36
37 FLUX_ELEM *flux_elem; /* supermodel flux coupling information */
38 int flux_elem_tot;
39
40} SELEM_2D;
41
42/*********************************************************/
43/* struct methods -------------------------------------- */
44
45void selem2d_get_triangle_local_shape(double xhat, double yhat, double zhat, double *lshape);
46void selem2d_get_quadrilateral_local_shape(double xhat, double yhat, double zhat, double *lshape);
47void selem2d_alloc(SELEM_2D *elem2d, int nnodes_on_elem);
48void selem2d_load(SELEM_2D *elem2d, int gid, int lid, int elem_nnodes, int *local_node_ids, int bflag, SVECT *nds, int mat);
49void selem2d_free(SELEM_2D *elem2d);
50void selem2d_alloc_array(SELEM_2D **elem2d, int nelems2d);
51void selem2d_free_array(SELEM_2D *elem2d, int nelems2d);
52void selem2d_init(SELEM_2D *elem2d);
53void selem2d_init_array(SELEM_2D *elem2d, int nelems2d);
54void selem2d_init_alloc_array(SELEM_2D **elem2d, int nelems2d);
55void selem2d_copy(SELEM_2D *to, SELEM_2D from);
56void selem2d_printScreen(SELEM_2D *elem2d);
57SVECT selem2d_get_triangle_centroid(SNODE nd1, SNODE nd2, SNODE nd3);
58double selem2d_get_triangle_orientation(SNODE nd1, SNODE nd2, SNODE nd3, SVECT ref_vec);
59void selem2d_get_triangle_linear_djac_nrml_gradPhi(SELEM_2D *elem2d, SNODE *nd_SNODE, SVECT *nd_SVECT);
60SVECT selem2d_get_elem2d_normals(SVECT *nd);
61double selem2d_get_quadrilateral_linear_djac2d(double xhat, double yhat, SVECT *nd);
62double selem2d_get_quadrilateral_linear_djac_gradPhi(double xhat, double yhat, SVECT *nd, SVECT *grad_shp);
63void selem2d_get_triangle_local_shape_quad(double xhat, double yhat, double zhat, double *lshape);
64void selem2d_get_quadrilateral_local_shape_quad(double xhat, double yhat, double zhat, double *lshape);
65void selem2d_get_triangle_local_shape_gradients(SVECT *lgrad_shp);
66void selem2d_get_quadrilateral_local_shape_gradients(double xhat, double yhat, double zhat, SVECT *lgrad_shp);
67double selem2d_get_elem2d_area2d(SELEM_2D *elem2d, SVECT *nds);
68
69/***********************************************************/
70/***********************************************************/
71/***********************************************************/
72
73#endif
Definition: selem_1d.h:11
Definition: selem_2d.h:10
Definition: snode.h:8
Definition: svect2d.h:4
Definition: svect.h:4