#include "basefun.h"
Go to the source code of this file.
Functions | |
void | a_coeff (double *a, double *x, double *y) |
void | a_coeff_3d (double *a, double *x, double *y, double *z, double v) |
void | ac_bf_quad_3_2d (double *n, double *l) |
void | ac_dx_bf_quad_3_2d (double *n, double *b, double *l) |
void | ac_dy_bf_quad_3_2d (double *n, double *c, double *l) |
void | b_coeff (double *b, double *y) |
void | b_coeff_3d (double *b, double *y, double *z, double v) |
void | bf_cct (double *n, double *l, double *sx, double *sy) |
void | bf_lin_1d (double *n, double x) |
void | bf_lin_3_2d (double *n, double x, double y) |
void | bf_lin_4_2d (double *n, double x, double y) |
void | bf_lin_hex_3d (double *n, double x, double y, double z) |
void | bf_lin_tet (double *n, double x, double y, double z) |
void | bf_lin_wed_3d (double *n, double x, double y, double z) |
void | bf_quad_1d (double *n, double x) |
void | bf_quad_3_2d (double *n, double x, double y) |
void | bf_quad_4_2d (double *n, double x, double y) |
void | bf_quad_hex_3d (double *n, double x, double y, double z) |
void | bf_quad_hexrot_3d (double *n, double x, double y, double z) |
void | bf_quad_tet (double *n, double x, double y, double z) |
void | bf_quad_wed_3d (double *n, double x, double y, double z) |
void | bf_rot_3_2d (double *n, double *l, double *b, double *c) |
void | bf_rot_4_2d (double *n, double x, double y, double *nx, double *ny, double *l) |
void | c_coeff (double *c, double *x) |
void | c_coeff_3d (double *c, double *x, double *z, double v) |
void | d_coeff_3d (double *d, double *x, double *y, double v) |
void | defl2D_fun (double *n, double xi, double l, double kappa) |
void | der_defl2D_fun (double *n, double xi, double l, double kappa) |
void | der_dilat (double *n, double l) |
void | der_roty_fun (double *n, double xi, double l, double kappa) |
void | dilat (double *n, double xi) |
void | dx_bf_lin_1d (double *n) |
void | dx_bf_lin_3_2d (double *n) |
void | dx_bf_lin_4_2d (double *n, double y) |
void | dx_bf_lin_hex_3d (double *n, double y, double z) |
void | dx_bf_lin_tet (double *n) |
void | dx_bf_lin_wed_3d (double *n, double z) |
void | dx_bf_quad_1d (double *n, double x) |
void | dx_bf_quad_3_2d (double *n, double x, double y) |
void | dx_bf_quad_4_2d (double *n, double x, double y) |
void | dx_bf_quad_hex_3d (double *n, double x, double y, double z) |
void | dx_bf_quad_hexrot_3d (double *n, double x, double y, double z) |
void | dx_bf_quad_tet (double *n, double x, double y, double z) |
void | dx_bf_quad_wed_3d (double *n, double x, double y, double z) |
void | dx_bf_rot_3_2d (double *n, double *l, double *b, double *c) |
void | dx_bf_rot_4_2d (double *n, double x, double y, double *nx, double *ny, double *l) |
void | dx_cct (double *n, double *l, double *b, double *sx, double *sy) |
void | dy_bf_lin_3_2d (double *n) |
void | dy_bf_lin_4_2d (double *n, double x) |
void | dy_bf_lin_hex_3d (double *n, double x, double z) |
void | dy_bf_lin_tet (double *n) |
void | dy_bf_lin_wed_3d (double *n, double z) |
void | dy_bf_quad_3_2d (double *n, double x, double y) |
void | dy_bf_quad_4_2d (double *n, double x, double y) |
void | dy_bf_quad_hex_3d (double *n, double x, double y, double z) |
void | dy_bf_quad_hexrot_3d (double *n, double x, double y, double z) |
void | dy_bf_quad_tet (double *n, double x, double y, double z) |
void | dy_bf_quad_wed_3d (double *n, double x, double y, double z) |
void | dy_bf_rot_3_2d (double *n, double *l, double *b, double *c) |
void | dy_bf_rot_4_2d (double *n, double x, double y, double *nx, double *ny, double *l) |
void | dy_cct (double *n, double *l, double *c, double *sx, double *sy) |
void | dz_bf_lin_hex_3d (double *n, double x, double y) |
void | dz_bf_lin_tet (double *n) |
void | dz_bf_lin_wed_3d (double *n, double x, double y) |
void | dz_bf_quad_hex_3d (double *n, double x, double y, double z) |
void | dz_bf_quad_hexrot_3d (double *n, double x, double y, double z) |
void | dz_bf_quad_tet (double *n, double x, double y, double z) |
void | dz_bf_quad_wed_3d (double *n, double x, double y, double z) |
void | plsa (double *a, double *x, double *y, double det) |
void | plsb (double *b, double *y, double det) |
void | plsc (double *c, double *x, double det) |
void | roty_fun (double *n, double xi, double l, double kappa) |
void | vola_3d (double *a, double *x, double *y, double *z, double det) |
void | volb_3d (double *b, double *y, double *z, double det) |
void | volc_3d (double *c, double *x, double *z, double det) |
void | vold_3d (double *d, double *x, double *y, double det) |
void a_coeff | ( | double * | a, | |
double * | x, | |||
double * | y | |||
) |
function computes a_i coefficients of area-coordinates on triangles L(x,y)_i = (a_i + b_i.x + c_i.y)/area/2
a | - coefficients | |
x,y | - node coordinates |
Definition at line 162 of file basefun.cpp.
void a_coeff_3d | ( | double * | a, | |
double * | x, | |||
double * | y, | |||
double * | z, | |||
double | v | |||
) |
function computes coefficients a_i of volume coordinates L(x,y,z)_i = (a_i + b_i.x + c_i.y + d_i.z)/6/volume
a | - array of coefficients | |
x,y,z | - array of node coordinates | |
v | - volume of element |
Definition at line 774 of file basefun.cpp.
void ac_bf_quad_3_2d | ( | double * | n, | |
double * | l | |||
) |
function computes quadratic base functions on triangle
n | - array of base functions | |
l | - areacoordinates |
Definition at line 260 of file basefun.cpp.
void ac_dx_bf_quad_3_2d | ( | double * | n, | |
double * | b, | |||
double * | l | |||
) |
function computes derivatives of quadratic base functions on triangle with respect of x
n | - array of derivatives of base functions | |
b | - array of coefficients of area coordinates | |
l | - areacoordinates |
Definition at line 284 of file basefun.cpp.
void ac_dy_bf_quad_3_2d | ( | double * | n, | |
double * | c, | |||
double * | l | |||
) |
function computes derivatives of quadratic base functions on triangle with respect of y
n | - array of derivatives of base functions | |
c | - array of coefficients of area coordinates | |
l | - areacoordinates |
Definition at line 308 of file basefun.cpp.
void b_coeff | ( | double * | b, | |
double * | y | |||
) |
function computes b_i coefficients of area-coordinates on triangles L(x,y)_i = (a_i + b_i.x + c_i.y)/area/2
b | - coefficients | |
y | - node coordinates |
Definition at line 178 of file basefun.cpp.
Referenced by planeelemrotlt::load_matrix(), and planeelemrotlt::mass_matrix().
void b_coeff_3d | ( | double * | b, | |
double * | y, | |||
double * | z, | |||
double | v | |||
) |
function computes coefficients b_i of volume coordinates L(x,y,z)_i = (a_i + b_i.x + c_i.y + d_i.z)/6/volume
b | - array of coefficients | |
y,z | - array of node coordinates | |
v | - volume of element |
Definition at line 807 of file basefun.cpp.
void bf_cct | ( | double * | n, | |
double * | l, | |||
double * | sx, | |||
double * | sy | |||
) |
function computes base functions of constant curvature triangle (finite element for plate problems)
n | - array of base functions | |
l | - area coordinates | |
sx,sy | - arrays of components |
Definition at line 697 of file basefun.cpp.
Referenced by cctelem::bf_matrix(), cctelem::geom_matrix(), and cctelem::geom_matrix_block().
void bf_lin_1d | ( | double * | n, | |
double | x | |||
) |
function computes linear base functions on bar element natural coordinate is from segment <-1;1>
n | - array of base functions | |
x | - natural coodinate |
Definition at line 10 of file basefun.cpp.
Referenced by linbartax::approx(), linbart::approx(), plquadcontact::approx(), barel3d::approx(), barel2d::approx(), gtopology::approx_weights(), linbartax::bf_matrix(), and linbart::bf_matrix().
void bf_lin_3_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes linear approximation functions on triangle
n | - array containing approximation functions | |
x,y | - natural coordinates |
1.4.2002
Definition at line 332 of file basefun.cpp.
Referenced by gtopology::approx_weights(), planeelemlt::bf_matrix(), axisymlt::bf_matrix(), and give_valuesinpoints().
void bf_lin_4_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes bi-linear base functions on quadrilateral
n | - array of base functions | |
x,y | - natural coodinates |
Definition at line 429 of file basefun.cpp.
Referenced by quadlineart::approx(), quadlinaxisym::approx(), soilplateq::approx(), q4plate::approx(), planeelemrotlq::approx(), planeelemlq::approx(), axisymlq::approx(), gtopology::approx_weights(), quadlineart::bf_matrix(), quadlinaxisym::bf_matrix(), planeelemlq::bf_matrix(), axisymlq::bf_matrix(), quadlineart::compute_error(), planeelemlq::compute_error(), axisymlq::geom_matrix(), axisymlq::geom_matrix_block(), quadlinaxisym::give_approx_fun(), give_valuesinpoints(), patch_averaging::nsp_spcoord_maxcoord_assembling(), planeelemlq::ntdbr_vector(), planeelemlq::ntn_matrix(), and radius().
void bf_lin_hex_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes tri-linear base functions on hexahedron
n | - array of base functions | |
x,y,z | - natural coodinates |
Definition at line 1284 of file basefun.cpp.
Referenced by linhext::approx(), linhexrot::approx(), linhex::approx(), gtopology::approx_weights(), linhext::bf_matrix(), linhexrot::bf_matrix(), linhex::bf_matrix(), and linhexrot::geom_matrix_shear().
void bf_lin_tet | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes tri-linear base functions on tetrahedron
n | - array of base functions | |
x,y,z | - natural coodinates |
Definition at line 949 of file basefun.cpp.
Referenced by gtopology::approx_weights().
void bf_lin_wed_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes linear base functions on wedge
n | - array of basis functions | |
x,y,z | - natural coordinates |
JK, 9.5.2004
Definition at line 1095 of file basefun.cpp.
Referenced by linwedge::approx(), and linwedge::bf_matrix().
void bf_quad_1d | ( | double * | n, | |
double | x | |||
) |
function computes quadratic approximation functions on bar element
n | - array containing approximation functions | |
x | - natural coordinate |
Definition at line 33 of file basefun.cpp.
Referenced by quadbartax::approx(), quadbart::approx(), barelq3d::approx(), barelq2d::approx(), gtopology::approx_weights(), quadbartax::bf_matrix(), and quadbart::bf_matrix().
void bf_quad_3_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes quadratic approximation functions on triangle
n | - array containing approximation functions | |
x,y | - natural coordinates |
1.4.2002
Definition at line 375 of file basefun.cpp.
Referenced by planeelemsubqt::approx(), planeelemqt::approx(), gtopology::approx_weights(), planeelemsubqt::bf_matrix(), planeelemqt::bf_matrix(), planeelemqt::compute_error(), give_valuesinpoints(), least_square::nsp_spcoord_maxcoord_assembling(), patch_averaging::nsp_spcoord_maxcoord_assembling(), planeelemqt::ntdbr_vector(), and planeelemqt::ntn_matrix().
void bf_quad_4_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes bi-quadratic base functions on quadrilateral
n | - array of base functions | |
x,y | - natural coodinates |
Definition at line 477 of file basefun.cpp.
Referenced by quadquadrilattax::approx(), quadquadrilatt::approx(), planeelemqq::approx(), axisymqq::approx(), gtopology::approx_weights(), quadquadrilattax::bf_matrix(), quadquadrilatt::bf_matrix(), planeelemqq::bf_matrix(), axisymqq::bf_matrix(), planeelemqq::compute_error(), axisymqq::geom_matrix(), axisymqq::geom_matrix_block(), quadquadrilattax::give_approx_fun(), quadquadrilatt::give_approx_fun(), give_valuesinpoints(), least_square::nsp_spcoord_maxcoord_assembling(), patch_averaging::nsp_spcoord_maxcoord_assembling(), planeelemqq::ntdbr_vector(), and planeelemqq::ntn_matrix().
void bf_quad_hex_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes tri-quadratic base functions on hexahedron
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates |
Definition at line 1362 of file basefun.cpp.
Referenced by quadhext::approx(), quadhex::approx(), gtopology::approx_weights(), quadhext::bf_matrix(), and quadhex::bf_matrix().
void bf_quad_hexrot_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
Definition at line 1497 of file basefun.cpp.
void bf_quad_tet | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes tri-quadratic base functions on tetrahedron
n | - array of base functions | |
x,y,z | - natural coodinates TKo - 11.2008 |
Definition at line 996 of file basefun.cpp.
Referenced by quadtet::approx(), gtopology::approx_weights(), and quadtet::bf_matrix().
void bf_quad_wed_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes quadratic approximation functions on wedge
n | - array of basis functions | |
x,y,z | - natural coordinates |
JK, 19.9.2004
Definition at line 1172 of file basefun.cpp.
Referenced by quadwedge::approx(), and quadwedge::bf_matrix().
void bf_rot_3_2d | ( | double * | n, | |
double * | l, | |||
double * | b, | |||
double * | c | |||
) |
function computes modified base functions on triangle for element with rotational DOFs
n | - array of base functions | |
l | - area coordinates | |
b,c | - coefficients of area coordinates |
Definition at line 538 of file basefun.cpp.
Referenced by planeelemrotlt::addgeommat(), and planeelemrotlt::bf_matrix().
void bf_rot_4_2d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double * | nx, | |||
double * | ny, | |||
double * | l | |||
) |
function computes modified base functions on quadrilateral for element with rotational DOFs
n | - array of base functions | |
x,y | - natural coordinates | |
nx,ny | - array of x and y components of outer normal vectors | |
l | - array of lengths of element edges |
Definition at line 611 of file basefun.cpp.
Referenced by planeelemrotlq::addgeommat(), and planeelemrotlq::bf_matrix().
void c_coeff | ( | double * | c, | |
double * | x | |||
) |
function computes c_i coefficients of area-coordinates on triangles L(x,y)_i = (a_i + b_i.x + c_i.y)/area/2
c | - coefficients | |
x | - node coordinates |
Definition at line 194 of file basefun.cpp.
Referenced by planeelemrotlt::load_matrix(), and planeelemrotlt::mass_matrix().
void c_coeff_3d | ( | double * | c, | |
double * | x, | |||
double * | z, | |||
double | v | |||
) |
function computes coefficients c_i of volume coordinates L(x,y,z)_i = (a_i + b_i.x + c_i.y + d_i.z)/6/volume
c | - array of coefficients | |
x,z | - array of node coordinates | |
v | - volume of element |
Definition at line 832 of file basefun.cpp.
void d_coeff_3d | ( | double * | d, | |
double * | x, | |||
double * | y, | |||
double | v | |||
) |
function computes coefficients d_i of volume coordinates L(x,y,z)_i = (a_i + b_i.x + c_i.y + d_i.z)/6/volume
d | - array of coefficients | |
x,y | - array of node coordinates | |
v | - volume of element |
Definition at line 857 of file basefun.cpp.
void defl2D_fun | ( | double * | n, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function computes approximation function of deflection of 2D beam element
n | - array containing approximation functions | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - 6.0*E*I/k/G/A/l/l |
20.2.2002
Definition at line 92 of file basefun.cpp.
Referenced by beamel2d::bf_matrix().
void der_defl2D_fun | ( | double * | n, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function computes derivative of deflection function with respect of x
n | - array containing approximation functions | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - 6.0*E*I/k/G/A/l/l |
20.2.2002
Definition at line 110 of file basefun.cpp.
Referenced by beamel2d::dbf_matrix(), and beamel2d::geom_matrix().
void der_dilat | ( | double * | n, | |
double | l | |||
) |
function computes derivative of approximation function of dilatation of 2D beam element with respect of x
n | - array containing approximation functions | |
l | - length of the beam |
20.2.2002
Definition at line 76 of file basefun.cpp.
Referenced by beamel2d::geom_matrix().
void der_roty_fun | ( | double * | n, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function computes derivative of function of rotation of 2D beam element with respect of x
n | - array containing approximation functions | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - 6.0*E*I/k/G/A/l/l |
20.2.2002
Definition at line 146 of file basefun.cpp.
Referenced by beamel2d::geom_matrix().
void dilat | ( | double * | n, | |
double | xi | |||
) |
function computes approximation function of dilatation of 2D beam element
n | - array containing approximation functions | |
xi | - natural coordinate from segment <0;1> |
20.2.2002
Definition at line 62 of file basefun.cpp.
Referenced by beamel2d::bf_matrix().
void dx_bf_lin_1d | ( | double * | n | ) |
function computes derivatives of linear base functions on bar element with respect of natural coordinate x
n | - array of derivatives of base functions |
Definition at line 21 of file basefun.cpp.
Referenced by derivatives_1d(), barel2d::geom_matrix(), linbartax::grad_matrix(), linbart::grad_matrix(), jac1d2d(), jac1d3d(), and jac_1d().
void dx_bf_lin_3_2d | ( | double * | n | ) |
function computes derivatives of linear approximation functions with respect to natural coordinate xi
n | - array containing derivatives |
1.4.2002
Definition at line 346 of file basefun.cpp.
Referenced by derivatives_2d(), jac2d3d(), and jac_2d().
void dx_bf_lin_4_2d | ( | double * | n, | |
double | y | |||
) |
function computes derivatives of bi-linear base functions on quadrilateral with respect of natural coordinate x
n | - array of derivatives of base functions | |
y | - natural coodinate |
Definition at line 445 of file basefun.cpp.
Referenced by planeelemlq::bvectors(), derivatives_2d(), planeelemlq::geom_matrix(), axisymlq::geom_matrix(), planeelemlq::geom_matrix_block(), axisymlq::geom_matrix_block(), quadlineart::grad_matrix(), quadlinaxisym::grad_matrix(), jac2d3d(), and jac_2d().
void dx_bf_lin_hex_3d | ( | double * | n, | |
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-linear base functions on hexahedron with respect of natural coordinate x
n | - array of derivatives of base functions | |
y,z | - natural coodinates |
Definition at line 1302 of file basefun.cpp.
Referenced by linhexrot::bvectors(), linhex::bvectors(), derivatives_3d(), linhexrot::geom_matrix(), linhex::geom_matrix(), linhex::geom_matrix_block(), linhexrot::geom_matrix_shear(), linhext::grad_matrix(), and jac_3d().
void dx_bf_lin_tet | ( | double * | n | ) |
Definition at line 959 of file basefun.cpp.
Referenced by derivatives_3d(), lintet::geom_matrix(), and jac_3d().
void dx_bf_lin_wed_3d | ( | double * | n, | |
double | z | |||
) |
function computes derivatives of linear base functions on wedge with respect to x
n | - array of basis functions | |
z | - natural coordinate |
JK, 9.5.2004
Definition at line 1113 of file basefun.cpp.
Referenced by derivatives_3d(), linwedge::geom_matrix(), and jac_3d().
void dx_bf_quad_1d | ( | double * | n, | |
double | x | |||
) |
function computes derivatives of quadratic approximation functions with respect of x on bar element
n | - array containing approximation functions | |
x | - natural coordinate |
Definition at line 47 of file basefun.cpp.
Referenced by derivatives_1d(), barelq3d::geom_matrix(), barelq2d::geom_matrix(), quadbartax::grad_matrix(), quadbart::grad_matrix(), jac1d2d(), jac1d3d(), and jac_1d().
void dx_bf_quad_3_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes derivatives of quadratic approximation functions on triangle with respect to xi
n | - array containing approximation functions | |
x,y | - natural coordinates |
1.4.2002
Definition at line 393 of file basefun.cpp.
Referenced by derivatives_2d(), planeelemsubqt::geom_matrix(), planeelemqt::geom_matrix(), planeelemsubqt::geom_matrix_block(), planeelemqt::geom_matrix_block(), jac2d3d(), and jac_2d().
void dx_bf_quad_4_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes derivatives of bi-quadratic base functions on quadrilateral with respect of natural coordinate x
n | - array of derivatives of base functions | |
x,y | - natural coodinates |
Definition at line 497 of file basefun.cpp.
Referenced by derivatives_2d(), planeelemqq::geom_matrix(), axisymqq::geom_matrix(), axisymqq::geom_matrix_block(), quadquadrilattax::grad_matrix(), quadquadrilatt::grad_matrix(), jac2d3d(), and jac_2d().
void dx_bf_quad_hex_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-quadratic base functions on hexahedron with respect of natural coordinate x
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates |
Definition at line 1398 of file basefun.cpp.
Referenced by derivatives_3d(), quadhex::geom_matrix(), quadhext::grad_matrix(), and jac_3d().
void dx_bf_quad_hexrot_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
Definition at line 1516 of file basefun.cpp.
Referenced by linhexrot::geom_matrix(), and linhexrot::geom_matrix_shear().
void dx_bf_quad_tet | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-quadratic base functions on tetrahedron with respect of natural coordinate x
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates TKo - 11.2008 |
Definition at line 1019 of file basefun.cpp.
Referenced by derivatives_3d(), quadtet::geom_matrix(), and jac_3d().
void dx_bf_quad_wed_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of linear approximation functions on wedge with respect to x
n | - array of basis functions | |
x,y,z | - natural coordinates |
JK, 19.9.2004
Definition at line 1199 of file basefun.cpp.
Referenced by derivatives_3d(), quadwedge::geom_matrix(), and jac_3d().
void dx_bf_rot_3_2d | ( | double * | n, | |
double * | l, | |||
double * | b, | |||
double * | c | |||
) |
function computes derivatives of modified base functions on triangle for element with rotational DOFs with respect of x
n | - array of base functions | |
l | - area coordinates | |
b,c | - coefficients of area coordinates |
Definition at line 562 of file basefun.cpp.
Referenced by planeelemrotlt::addgeommat(), planeelemrotlt::geom_matrix(), and planeelemrotlt::geom_matrix_block().
void dx_bf_rot_4_2d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double * | nx, | |||
double * | ny, | |||
double * | l | |||
) |
function computes derivatives of modified base functions on quadrilateral for element with rotational DOFs with respect of natural coordinate x
n | - array of base functions | |
x,y | - natural coordinates | |
nx,ny | - array of x and y components of outer normal vectors | |
l | - array of lengths of element edges |
Definition at line 640 of file basefun.cpp.
Referenced by planeelemrotlq::addgeommat(), planeelemrotlq::geom_matrix(), and planeelemrotlq::geom_matrix_block().
void dx_cct | ( | double * | n, | |
double * | l, | |||
double * | b, | |||
double * | sx, | |||
double * | sy | |||
) |
function computes derivates of base functions of constant curvature triangle with respect of x (finite element for plate problems)
n | - array of base functions | |
l | - area coordinates | |
b | - array of coefficients of area coordinates | |
sx,sy | - arrays of components |
Definition at line 723 of file basefun.cpp.
Referenced by cctelem::geom_matrix(), and cctelem::geom_matrix_block().
void dy_bf_lin_3_2d | ( | double * | n | ) |
function computes derivatives of linear approximation functions with respect to natural coordinate eta
n | - array containing derivatives |
1.4.2002
Definition at line 360 of file basefun.cpp.
Referenced by derivatives_2d(), jac2d3d(), and jac_2d().
void dy_bf_lin_4_2d | ( | double * | n, | |
double | x | |||
) |
function computes derivatives of bi-linear base functions on quadrilateral with respect of natural coordinate y
n | - array of derivatives of base functions | |
x | - natural coodinate |
Definition at line 461 of file basefun.cpp.
Referenced by planeelemlq::bvectors(), derivatives_2d(), planeelemlq::geom_matrix(), axisymlq::geom_matrix(), planeelemlq::geom_matrix_block(), axisymlq::geom_matrix_block(), quadlineart::grad_matrix(), quadlinaxisym::grad_matrix(), jac2d3d(), and jac_2d().
void dy_bf_lin_hex_3d | ( | double * | n, | |
double | x, | |||
double | z | |||
) |
function computes derivatives of tri-linear base functions on hexahedron with respect of natural coordinate y
n | - array of derivatives of base functions | |
x,z | - natural coodinates |
Definition at line 1322 of file basefun.cpp.
Referenced by linhexrot::bvectors(), linhex::bvectors(), derivatives_3d(), linhexrot::geom_matrix(), linhex::geom_matrix(), linhex::geom_matrix_block(), linhexrot::geom_matrix_shear(), linhext::grad_matrix(), and jac_3d().
void dy_bf_lin_tet | ( | double * | n | ) |
Definition at line 969 of file basefun.cpp.
Referenced by derivatives_3d(), lintet::geom_matrix(), and jac_3d().
void dy_bf_lin_wed_3d | ( | double * | n, | |
double | z | |||
) |
function computes derivatives of linear base functions on wedge with respect to y
n | - array of basis functions | |
z | - natural coordinate |
JK, 9.5.2004
Definition at line 1131 of file basefun.cpp.
Referenced by derivatives_3d(), linwedge::geom_matrix(), and jac_3d().
void dy_bf_quad_3_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes derivatives of quadratic approximation functions on triangle with respect to eta
n | - array containing approximation functions | |
x,y | - natural coordinates |
1.4.2002
Definition at line 411 of file basefun.cpp.
Referenced by derivatives_2d(), planeelemsubqt::geom_matrix(), planeelemqt::geom_matrix(), planeelemsubqt::geom_matrix_block(), planeelemqt::geom_matrix_block(), jac2d3d(), and jac_2d().
void dy_bf_quad_4_2d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes derivatives of bi-quadratic base functions on quadrilateral with respect of natural coordinate y
n | - array of derivatives of base functions | |
x,y | - natural coodinates |
Definition at line 517 of file basefun.cpp.
Referenced by derivatives_2d(), planeelemqq::geom_matrix(), axisymqq::geom_matrix(), axisymqq::geom_matrix_block(), quadquadrilattax::grad_matrix(), quadquadrilatt::grad_matrix(), jac2d3d(), and jac_2d().
void dy_bf_quad_hex_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-quadratic base functions on hexahedron with respect of natural coordinate y
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates |
Definition at line 1434 of file basefun.cpp.
Referenced by derivatives_3d(), quadhex::geom_matrix(), quadhext::grad_matrix(), and jac_3d().
void dy_bf_quad_hexrot_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
Definition at line 1535 of file basefun.cpp.
Referenced by linhexrot::geom_matrix(), and linhexrot::geom_matrix_shear().
void dy_bf_quad_tet | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-quadratic base functions on tetrahedron with respect of natural coordinate y
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates TKo - 11.2008 |
Definition at line 1044 of file basefun.cpp.
Referenced by derivatives_3d(), quadtet::geom_matrix(), and jac_3d().
void dy_bf_quad_wed_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of quadratic approximation functions on wedge with respect to y
n | - array of basis functions | |
x,y,z | - natural coordinates |
JK, 19.9.2004
Definition at line 1226 of file basefun.cpp.
Referenced by derivatives_3d(), quadwedge::geom_matrix(), and jac_3d().
void dy_bf_rot_3_2d | ( | double * | n, | |
double * | l, | |||
double * | b, | |||
double * | c | |||
) |
function computes derivatives of modified base functions on triangle for element with rotational DOFs with respect of y
n | - array of base functions | |
l | - area coordinates | |
b,c | - coefficients of area coordinates |
Definition at line 586 of file basefun.cpp.
Referenced by planeelemrotlt::addgeommat(), planeelemrotlt::geom_matrix(), and planeelemrotlt::geom_matrix_block().
void dy_bf_rot_4_2d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double * | nx, | |||
double * | ny, | |||
double * | l | |||
) |
function computes derivatives of modified base functions on quadrilateral for element with rotational DOFs with respect of natural coordinate y
n | - array of base functions | |
x,y | - natural coordinates | |
nx,ny | - array of x and y components of outer normal vectors | |
l | - array of lengths of element edges |
Definition at line 669 of file basefun.cpp.
Referenced by planeelemrotlq::addgeommat(), planeelemrotlq::geom_matrix(), and planeelemrotlq::geom_matrix_block().
void dy_cct | ( | double * | n, | |
double * | l, | |||
double * | c, | |||
double * | sx, | |||
double * | sy | |||
) |
function computes derivates of base functions of constant curvature triangle with respect of y (finite element for plate problems)
n | - array of base functions | |
l | - area coordinates | |
b | - array of coefficients of area coordinates | |
sx,sy | - arrays of components |
Definition at line 749 of file basefun.cpp.
Referenced by cctelem::geom_matrix(), and cctelem::geom_matrix_block().
void dz_bf_lin_hex_3d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes derivatives of tri-linear base functions on hexahedron with respect of natural coordinate z
n | - array of derivatives of base functions | |
x,y | - natural coodinates |
Definition at line 1342 of file basefun.cpp.
Referenced by linhexrot::bvectors(), linhex::bvectors(), derivatives_3d(), linhexrot::geom_matrix(), linhex::geom_matrix(), linhex::geom_matrix_block(), linhexrot::geom_matrix_shear(), linhext::grad_matrix(), and jac_3d().
void dz_bf_lin_tet | ( | double * | n | ) |
Definition at line 979 of file basefun.cpp.
Referenced by derivatives_3d(), lintet::geom_matrix(), and jac_3d().
void dz_bf_lin_wed_3d | ( | double * | n, | |
double | x, | |||
double | y | |||
) |
function computes derivatives of linear base functions on wedge with respect to z
n | - array of basis functions | |
x,y | - natural coordinates |
JK, 9.5.2004
Definition at line 1149 of file basefun.cpp.
Referenced by derivatives_3d(), linwedge::geom_matrix(), and jac_3d().
void dz_bf_quad_hex_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-quadratic base functions on hexahedron with respect of natural coordinate z
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates |
Definition at line 1468 of file basefun.cpp.
Referenced by derivatives_3d(), quadhex::geom_matrix(), quadhext::grad_matrix(), and jac_3d().
void dz_bf_quad_hexrot_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
Definition at line 1555 of file basefun.cpp.
Referenced by linhexrot::geom_matrix(), and linhexrot::geom_matrix_shear().
void dz_bf_quad_tet | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of tri-quadratic base functions on tetrahedron with respect of natural coordinate z
n | - array of derivatives of base functions | |
x,y,z | - natural coodinates TKo - 11.2008 |
Definition at line 1069 of file basefun.cpp.
Referenced by derivatives_3d(), quadtet::geom_matrix(), and jac_3d().
void dz_bf_quad_wed_3d | ( | double * | n, | |
double | x, | |||
double | y, | |||
double | z | |||
) |
function computes derivatives of quadratic approximation functions on wedge with respect to z
n | - array of basis functions | |
x,y,z | - natural coordinates |
JK, 19.9.2004
Definition at line 1253 of file basefun.cpp.
Referenced by derivatives_3d(), quadwedge::geom_matrix(), and jac_3d().
void plsa | ( | double * | a, | |
double * | x, | |||
double * | y, | |||
double | det | |||
) |
function computes a_i coefficients of area-coordinates on triangles L(x,y)_i = a_i + b_i.x + c_i.y
a | - coefficients | |
x,y | - node coordinates | |
det | - determinant=2*area |
Definition at line 211 of file basefun.cpp.
void plsb | ( | double * | b, | |
double * | y, | |||
double | det | |||
) |
function computes b_i coefficients of area-coordinates on triangles L(x,y)_i = a_i + b_i.x + c_i.y
b | - coefficients | |
y | - node coordinates | |
det | - determinant=2*area |
Definition at line 228 of file basefun.cpp.
Referenced by planeelemrotlt::addgeommat(), trlineart::conductivity_matrix(), trlinaxisym::conductivity_matrix(), planeelemrotlt::geom_matrix(), planeelemlt::geom_matrix(), cctelem::geom_matrix(), axisymlt::geom_matrix(), planeelemrotlt::geom_matrix_block(), cctelem::geom_matrix_block(), axisymlt::geom_matrix_block(), trlineart::internal_fluxes(), trlinaxisym::internal_fluxes(), trlineart::intpointgrad(), trlinaxisym::intpointgrad(), trlineart::l_matrix(), and trlineart::l_t_matrix().
void plsc | ( | double * | c, | |
double * | x, | |||
double | det | |||
) |
function computes c_i coefficients of area-coordinates on triangles L(x,y)_i = a_i + b_i.x + c_i.y
c | - coefficients | |
x | - node coordinates | |
det | - determinant=2*area |
Definition at line 245 of file basefun.cpp.
Referenced by planeelemrotlt::addgeommat(), trlineart::conductivity_matrix(), trlinaxisym::conductivity_matrix(), planeelemrotlt::geom_matrix(), planeelemlt::geom_matrix(), cctelem::geom_matrix(), axisymlt::geom_matrix(), planeelemrotlt::geom_matrix_block(), cctelem::geom_matrix_block(), axisymlt::geom_matrix_block(), trlineart::internal_fluxes(), trlinaxisym::internal_fluxes(), trlineart::intpointgrad(), trlinaxisym::intpointgrad(), trlineart::l_matrix(), and trlineart::l_t_matrix().
void roty_fun | ( | double * | n, | |
double | xi, | |||
double | l, | |||
double | kappa | |||
) |
function computes approximation function of rotation of 2D beam element
n | - array containing approximation functions | |
xi | - natural coordinate from segment <0;1> | |
l | - length of the beam | |
kappa | - 6.0*E*I/k/G/A/l/l |
20.2.2002
Definition at line 128 of file basefun.cpp.
Referenced by beamel2d::bf_matrix(), and beamel2d::geom_matrix().
void vola_3d | ( | double * | a, | |
double * | x, | |||
double * | y, | |||
double * | z, | |||
double | det | |||
) |
function computes coefficients a_i of volume coordinates L(x,y,z)_i = a_i + b_i.x + c_i.y + d_i.z
a | - array of coefficients | |
x,y,z | - array of node coordinates | |
det | - determinat=6*volume |
Definition at line 881 of file basefun.cpp.
void volb_3d | ( | double * | b, | |
double * | y, | |||
double * | z, | |||
double | det | |||
) |
function computes coefficients b_i of volume coordinates L(x,y,z)_i = a_i + b_i.x + c_i.y + d_i.z
b | - array of coefficients | |
y,z | - array of node coordinates | |
det | - determinat=6*volume |
Definition at line 899 of file basefun.cpp.
Referenced by lintett::conductivity_matrix(), lintetrot::geom_matrix(), lintet::geom_matrix_old(), lintetrot::geom_matrix_shear(), lintett::internal_fluxes(), lintett::intpointgrad(), lintett::l_matrix(), and lintett::l_t_matrix().
void volc_3d | ( | double * | c, | |
double * | x, | |||
double * | z, | |||
double | det | |||
) |
function computes coefficients c_i of volume coordinates L(x,y,z)_i = a_i + b_i.x + c_i.y + d_i.z
c | - array of coefficients | |
x,z | - array of node coordinates | |
det | - determinat=6*volume |
Definition at line 917 of file basefun.cpp.
Referenced by lintett::conductivity_matrix(), lintetrot::geom_matrix(), lintet::geom_matrix_old(), lintetrot::geom_matrix_shear(), lintett::internal_fluxes(), lintett::intpointgrad(), lintett::l_matrix(), and lintett::l_t_matrix().
void vold_3d | ( | double * | d, | |
double * | x, | |||
double * | y, | |||
double | det | |||
) |
function computes coefficients d_i of volume coordinates L(x,y,z)_i = a_i + b_i.x + c_i.y + d_i.z
d | - array of coefficients | |
x,y | - array of node coordinates | |
det | - determinat=6*volume |
Definition at line 933 of file basefun.cpp.
Referenced by lintett::conductivity_matrix(), lintetrot::geom_matrix(), lintet::geom_matrix_old(), lintetrot::geom_matrix_shear(), lintett::internal_fluxes(), lintett::intpointgrad(), lintett::l_matrix(), and lintett::l_t_matrix().