Trilinos based (stochastic) FEM solvers
fepp.hpp
Go to the documentation of this file.
1 /*
2 Brian Staber (brian.staber@gmail.com)
3 */
4 
5 #ifndef FEPP_HPP
6 #define FEPP_HPP
7 
8 #include "meshpp.hpp"
9 
10 namespace tri3{
11  void shape_functions(Epetra_SerialDenseVector & N, double & xi, double & eta);
12  void d_shape_functions(Epetra_SerialDenseMatrix & D, double & xi, double & eta);
13  void dX_shape_functions(Epetra_SerialDenseMatrix & D, Epetra_SerialDenseMatrix & DX, double & jac, Epetra_SerialDenseMatrix & JacobianMatrix);
14 }
15 
16 namespace quad4{
17  void shape_functions(Epetra_SerialDenseVector & N, double & xi, double & eta);
18  void d_shape_functions(Epetra_SerialDenseMatrix & D, double & xi, double & eta);
19  void dX_shape_functions(Epetra_SerialDenseMatrix & D, Epetra_SerialDenseMatrix & DX, double & jac, Epetra_SerialDenseMatrix & JacobianMatrix);
20 }
21 
22 namespace tri6{
23  void shape_functions(Epetra_SerialDenseVector & N, double & xi, double & eta);
24  void d_shape_functions(Epetra_SerialDenseMatrix & D, double & xi, double & eta);
25  void dX_shape_functions(Epetra_SerialDenseMatrix & D, Epetra_SerialDenseMatrix & DX, double & jac, Epetra_SerialDenseMatrix & JacobianMatrix);
26 }
27 
28 namespace tetra4{
29  void shape_functions(Epetra_SerialDenseVector & N, double & xi, double & eta, double & zeta);
30  void d_shape_functions(Epetra_SerialDenseMatrix & D, double & xi, double & eta, double & zeta);
31 }
32 
33 namespace tetra10{
34  void shape_functions(Epetra_SerialDenseVector & N, double & xi, double & eta, double & zeta);
35  void d_shape_functions(Epetra_SerialDenseMatrix & D, double & xi, double & eta, double & zeta);
36 }
37 
38 namespace hexa8{
39  void shape_functions(Epetra_SerialDenseVector & N, double & xi, double & eta, double & zeta);
40  void d_shape_functions(Epetra_SerialDenseMatrix & D, double & xi, double & eta, double & zeta);
41 }
42 
43 void dX_shape_functions(Epetra_SerialDenseMatrix & D, Epetra_SerialDenseMatrix JacobianMatrix, double & jac, Epetra_SerialDenseMatrix & DX);
44 
45 void jacobian_matrix(Epetra_SerialDenseMatrix & X, Epetra_SerialDenseMatrix & D, Epetra_SerialDenseMatrix & JacobianMatrix);
46 void jacobian_det(Epetra_SerialDenseMatrix & JacobianMatrix, double & jac);
47 void jacobian_det_faces(Epetra_SerialDenseMatrix & JacobianMatrix, double & jac);
48 
49 void gauss_points_tri1(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta);
50 void gauss_points_tri3(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta);
51 void gauss_points_tri4(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta);
52 
53 void gauss_points_quad1(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta);
54 void gauss_points_quad4(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta);
55 void gauss_points_quad9(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta);
56 
57 void gauss_points_hexa4(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
58 void gauss_points_hexa8(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
59 void gauss_points_hexa27(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
60 
61 void gauss_points_tetra1(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
62 void gauss_points_tetra4(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
63 void gauss_points_tetra5(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
64 void gauss_points_tetra11(Epetra_SerialDenseVector & weight, Epetra_SerialDenseVector & xi, Epetra_SerialDenseVector & eta, Epetra_SerialDenseVector & zeta);
65 
66 int pnpoly(int & nvert, Epetra_SerialDenseVector & vertx, Epetra_SerialDenseVector & verty, double & testx, double & testy);
67 
68 
69 #endif
void gauss_points_tetra11(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:307
void gauss_points_tri4(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta)
Definition: fepp.cpp:195
Definition: fepp.hpp:10
void gauss_points_hexa8(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:234
void jacobian_det_faces(Epetra_SerialDenseMatrix &JacobianMatrix, double &jac)
Definition: fepp.cpp:179
void gauss_points_quad9(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta)
Definition: fepp.cpp:215
void gauss_points_tri3(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta)
Definition: fepp.cpp:189
void gauss_points_quad4(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta)
Definition: fepp.cpp:208
void jacobian_matrix(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &D, Epetra_SerialDenseMatrix &JacobianMatrix)
Definition: fepp.cpp:171
int pnpoly(int &nvert, Epetra_SerialDenseVector &vertx, Epetra_SerialDenseVector &verty, double &testx, double &testy)
Definition: fepp.cpp:315
void dX_shape_functions(Epetra_SerialDenseMatrix &D, Epetra_SerialDenseMatrix &DX, double &jac, Epetra_SerialDenseMatrix &JacobianMatrix)
Definition: fepp.cpp:17
Definition: fepp.hpp:16
Definition: fepp.hpp:28
void shape_functions(Epetra_SerialDenseVector &N, double &xi, double &eta)
Definition: fepp.cpp:7
void gauss_points_hexa4(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:223
void jacobian_det(Epetra_SerialDenseMatrix &JacobianMatrix, double &jac)
Definition: fepp.cpp:175
Definition: fepp.hpp:22
void gauss_points_tri1(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta)
Definition: fepp.cpp:183
void gauss_points_hexa27(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:248
void d_shape_functions(Epetra_SerialDenseMatrix &D, double &xi, double &eta)
Definition: fepp.cpp:12
Definition: fepp.hpp:38
void gauss_points_tetra5(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:299
void gauss_points_tetra1(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:281
Definition: fepp.hpp:33
void gauss_points_quad1(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta)
Definition: fepp.cpp:202
output eta
Definition: costFunction.m:33
void gauss_points_tetra4(Epetra_SerialDenseVector &weight, Epetra_SerialDenseVector &xi, Epetra_SerialDenseVector &eta, Epetra_SerialDenseVector &zeta)
Definition: fepp.cpp:289