29 #ifndef OOMPH_SHAPE_HEADER
30 #define OOMPH_SHAPE_HEADER
35 #include <oomph-lib-config.h>
100 std::ostringstream error_stream;
101 error_stream <<
"Range Error: ";
104 error_stream <<
i <<
" is not in the range (0," <<
Index1 - 1 <<
")"
109 error_stream <<
j <<
" is not in the range (0," <<
Index2 - 1 <<
")"
147 std::ostringstream error_stream;
148 error_stream <<
"Cannot assign Shape object:" << std::endl
149 <<
"Indices do not match "
151 <<
", RHS: " <<
shape.Index1 <<
" " <<
shape.Index2
168 std::ostringstream error_stream;
169 error_stream <<
"Cannot assign Shape object:" << std::endl
170 <<
"Indices do not match "
172 <<
", RHS: " << shape_pt->
Index1 <<
" " << shape_pt->
Index2
189 void resize(
const unsigned&
N,
const unsigned&
M = 1)
206 #ifdef RANGE_CHECKING
215 #ifdef RANGE_CHECKING
224 #ifdef RANGE_CHECKING
233 #ifdef RANGE_CHECKING
242 #ifdef RANGE_CHECKING
250 inline const double&
operator()(
const unsigned&
i,
const unsigned&
j)
const
252 #ifdef RANGE_CHECKING
302 const unsigned&
k)
const
307 std::ostringstream error_stream;
308 error_stream <<
"Range Error: ";
311 error_stream <<
i <<
" is not in the range (0," <<
Index1 - 1 <<
")"
316 error_stream <<
j <<
" is not in the range (0," <<
Index2 - 1 <<
")"
321 error_stream <<
k <<
" is not in the range (0," <<
Index3 - 1 <<
")"
340 DShape(
const unsigned&
N,
const unsigned&
M,
const unsigned&
P)
363 std::ostringstream error_stream;
364 error_stream <<
"Cannot assign DShape object:" << std::endl
365 <<
"Indices do not match "
368 <<
" " <<
dshape.Index3 << std::endl;
385 std::ostringstream error_stream;
386 error_stream <<
"Cannot assign Shape object:" << std::endl
387 <<
"Indices do not match "
389 <<
", RHS: " << dshape_pt->
Index1 <<
" "
410 void resize(
const unsigned&
N,
const unsigned&
P,
const unsigned&
M = 1)
428 #ifdef RANGE_CHECKING
435 inline const double&
operator()(
const unsigned&
i,
const unsigned&
k)
const
437 #ifdef RANGE_CHECKING
448 #ifdef RANGE_CHECKING
457 const unsigned&
k)
const
459 #ifdef RANGE_CHECKING
487 unsigned offset(
const unsigned long&
i,
const unsigned long&
j)
const
558 namespace OneDimLagrange
563 template<
unsigned NNODE_1D>
566 std::ostringstream error_stream;
567 error_stream <<
"One dimensional Lagrange shape functions "
568 <<
"have not been defined "
569 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
577 template<
unsigned NNODE_1D>
580 std::ostringstream error_stream;
581 error_stream <<
"One dimensional Lagrange shape function derivatives "
582 <<
"have not been defined "
583 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
592 template<
unsigned NNODE_1D>
595 std::ostringstream error_stream;
596 error_stream <<
"One dimensional Lagrange shape function "
597 <<
"second derivatives "
598 <<
"have not been defined "
599 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
610 Psi[0] = 0.5 * (1.0 -
s);
611 Psi[1] = 0.5 * (1.0 +
s);
637 Psi[0] = 0.5 *
s * (
s - 1.0);
638 Psi[1] = 1.0 -
s *
s;
639 Psi[2] = 0.5 *
s * (
s + 1.0);
670 double t3 = 0.5625 * t2;
671 double t4 = 0.5625 * t1;
672 double t5 = 0.625E-1 *
s;
673 double t7 = 0.16875E1 * t2;
674 double t8 = 0.16875E1 *
s;
675 Psi[0] = -t3 + t4 + t5 - 0.625E-1;
676 Psi[1] = t7 - t4 - t8 + 0.5625;
677 Psi[2] = -t7 - t4 + t8 + 0.5625;
678 Psi[3] = t3 + t4 - t5 - 0.625E-1;
688 double t2 = 0.16875E1 * t1;
689 double t3 = 0.1125E1 *
s;
690 double t5 = 0.50625E1 * t1;
691 DPsi[0] = -t2 + t3 + 0.625E-1;
692 DPsi[1] = t5 - t3 - 0.16875E1;
693 DPsi[2] = -t5 - t3 + 0.16875E1;
694 DPsi[3] = t2 + t3 - 0.625E-1;
704 double t2 = 0.16875E1 * t1;
705 double t5 = 0.50625E1 * t1;
706 DPsi[0] = -t2 + 0.1125E1;
707 DPsi[1] = t5 - 0.1125E1;
708 DPsi[2] = -t5 - 0.1125E1;
709 DPsi[3] = t2 + 0.1125E1;
718 namespace OneDimDiscontinuousGalerkin
723 template<
unsigned NNODE_1D>
727 std::ostringstream error_stream;
730 error_stream <<
"One dimensional Lagrange shape functions "
731 <<
"have not been defined "
732 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
742 template<
unsigned NNODE_1D>
746 std::ostringstream error_stream;
749 error_stream <<
"One dimensional Lagrange shape function derivatives "
750 <<
"have not been defined "
751 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
762 template<
unsigned NNODE_1D>
766 std::ostringstream error_stream;
769 error_stream <<
"One dimensional Lagrange shape function "
770 <<
"second derivatives "
771 <<
"have not been defined "
772 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
813 Psi[1] = 0.5 * (1.0 -
s);
814 Psi[2] = 0.5 * (1.0 +
s);
842 Psi[1] = 0.5 *
s * (
s - 1.0);
843 Psi[2] = 1.0 -
s *
s;
844 Psi[3] = 0.5 *
s * (
s + 1.0);
875 namespace OneDimDiscontinuousGalerkinMixedOrderBasis
880 template<
unsigned NNODE_1D>
884 std::ostringstream error_stream;
887 error_stream <<
"One dimensional Lagrange shape functions "
888 <<
"have not been defined "
889 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
899 template<
unsigned NNODE_1D>
903 std::ostringstream error_stream;
906 error_stream <<
"One dimensional Lagrange shape function derivatives "
907 <<
"have not been defined "
908 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
919 template<
unsigned NNODE_1D>
923 std::ostringstream error_stream;
926 error_stream <<
"One dimensional Lagrange shape function "
927 <<
"second derivatives "
928 <<
"have not been defined "
929 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
942 Psi[0] = 0.5 * (1.0 -
s);
943 Psi[1] = 0.5 * (1.0 +
s);
969 Psi[0] = 0.5 * (1.0 -
s);
971 Psi[2] = 0.5 * (1.0 +
s);
998 Psi[0] = 0.5 * (1.0 -
s);
1001 Psi[3] = 0.5 * (1.0 +
s);
1032 namespace OneDimDiscontinuousGalerkinMixedOrderTest
1037 template<
unsigned NNODE_1D>
1041 std::ostringstream error_stream;
1044 error_stream <<
"One dimensional Lagrange shape functions "
1045 <<
"have not been defined "
1046 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
1056 template<
unsigned NNODE_1D>
1060 std::ostringstream error_stream;
1063 error_stream <<
"One dimensional Lagrange shape function derivatives "
1064 <<
"have not been defined "
1065 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
1076 template<
unsigned NNODE_1D>
1080 std::ostringstream error_stream;
1083 error_stream <<
"One dimensional Lagrange shape function "
1084 <<
"second derivatives "
1085 <<
"have not been defined "
1086 <<
"for " << NNODE_1D <<
" nodes." << std::endl;
1187 namespace OneDimHermite
1194 inline void shape(
const double&
s,
double Psi[2][2])
1197 Psi[0][0] = 0.25 * (
s *
s *
s - 3.0 *
s + 2.0);
1198 Psi[0][1] = 0.25 * (
s *
s *
s -
s *
s -
s + 1.0);
1200 Psi[1][0] = 0.25 * (2.0 + 3.0 *
s -
s *
s *
s);
1201 Psi[1][1] = 0.25 * (
s *
s *
s +
s *
s -
s - 1.0);
1206 inline void dshape(
const double&
s,
double DPsi[2][2])
1209 DPsi[0][0] = 0.75 * (
s *
s - 1.0);
1210 DPsi[0][1] = 0.25 * (3.0 *
s *
s - 2.0 *
s - 1.0);
1212 DPsi[1][0] = 0.75 * (1.0 -
s *
s);
1213 DPsi[1][1] = 0.25 * (3.0 *
s *
s + 2.0 *
s - 1.0);
1217 inline void d2shape(
const double&
s,
double DPsi[2][2])
1220 DPsi[0][0] = 1.5 *
s;
1221 DPsi[0][1] = 0.5 * (3.0 *
s - 1.0);
1223 DPsi[1][0] = -1.5 *
s;
1224 DPsi[1][1] = 0.5 * (3.0 *
s + 1.0);
1232 template<
unsigned NNODE_1D>
1258 using namespace Orthpoly;
1260 unsigned p = NNODE_1D - 1;
1262 for (
unsigned i = 0;
i < NNODE_1D;
i++)
1279 template<
unsigned NNODE_1D>
1282 template<
unsigned NNODE_1D>
1286 template<
unsigned NNODE_1D>
1293 unsigned p = NNODE_1D - 1;
1299 for (
unsigned i = 0;
i < NNODE_1D;
i++)
1301 unsigned rootnum = 0;
1302 for (
unsigned j = 0;
j < NNODE_1D;
j++)
1313 if (
i == rootnum &&
i == 0)
1315 (*this)[
i] = -(1.0 +
p) *
p / 4.0;
1317 else if (
i == rootnum &&
i ==
p)
1319 (*this)[
i] = (1.0 +
p) *
p / 4.0;
1321 else if (
i == rootnum)
1356 (*this)[0] = 0.5 * (1.0 -
s);
1357 (*this)[1] = 0.5 * (1.0 +
s);
1358 for (
unsigned i = 2;
i < p_order;
i++)
1376 for (
unsigned i = 2;
i < p_order;
i++)
1378 (*this)[
i] = (0.5 * (1.0 -
s)) * (0.5 * (1.0 +
s)) *
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
float * p
Definition: Tutorial_Map_using.cpp:9
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:50
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
unsigned long nindex3() const
Return the range of index 3 of the derivatives of the shape functions.
Definition: shape.h:506
unsigned long nindex2() const
Return the range of index 2 of the derivatives of the shape functions.
Definition: shape.h:500
unsigned offset(const unsigned long &i, const unsigned long &j) const
Definition: shape.h:487
unsigned Index2
Size of the second index of the shape function.
Definition: shape.h:294
const double & raw_direct_access(const unsigned long &i) const
Definition: shape.h:478
~DShape()
Destructor, clean up the memory allocated by this object.
Definition: shape.h:401
void operator=(DShape *const &dshape_pt)
Definition: shape.h:378
void operator=(const DShape &dshape)
Definition: shape.h:356
const double & operator()(const unsigned &i, const unsigned &k) const
Overload the round bracket operator (const version)
Definition: shape.h:435
const double & operator()(const unsigned &i, const unsigned &j, const unsigned &k) const
Overload the round bracket operator (const version)
Definition: shape.h:455
unsigned long nindex1() const
Return the range of index 1 of the derivatives of the shape functions.
Definition: shape.h:494
double & raw_direct_access(const unsigned long &i)
Definition: shape.h:469
DShape()
Definition: shape.h:349
double * Allocated_storage
Definition: shape.h:288
void resize(const unsigned &N, const unsigned &P, const unsigned &M=1)
Definition: shape.h:410
unsigned Index1
Size of the first index of the shape function.
Definition: shape.h:291
unsigned Index3
Size of the third index of the shape function.
Definition: shape.h:297
DShape(const unsigned &N, const unsigned &M, const unsigned &P)
Constructor with three paramters: a two-index shape function.
Definition: shape.h:340
void range_check(const unsigned &i, const unsigned &j, const unsigned &k) const
Private function that checks whether the indices are in range.
Definition: shape.h:300
DShape(const unsigned &N, const unsigned &P)
Constructor with two parameters: a single-index shape function.
Definition: shape.h:332
double & operator()(const unsigned &i, const unsigned &k)
Overload the round bracket operator for access to the data.
Definition: shape.h:426
DShape(const DShape &dshape)=delete
Broken copy constructor.
double & operator()(const unsigned &i, const unsigned &j, const unsigned &k)
Overload the round bracket operator, with 3 indices.
Definition: shape.h:444
double * DPsi
Definition: shape.h:283
OneDimensionalLegendreDShape(const double &s)
Definition: shape.h:1291
Class that returns the shape functions associated with legendre.
Definition: shape.h:1234
OneDimensionalLegendreShape(const double &s)
Constructor.
Definition: shape.h:1256
static double nodal_position(const unsigned &n)
Definition: shape.h:1250
static Vector< double > z
Definition: shape.h:1238
static bool Nodes_calculated
Definition: shape.h:1235
static void calculate_nodal_positions()
Static function used to populate the stored positions.
Definition: shape.h:1241
OneDimensionalModalDShape(const unsigned p_order, const double &s)
Definition: shape.h:1370
OneDimensionalModalShape(const unsigned p_order, const double &s)
Constructor.
Definition: shape.h:1352
Definition: oomph_definitions.h:222
ShapeWithDeepCopy()
Default constructor.
Definition: shape.h:528
ShapeWithDeepCopy(const unsigned &N)
Constructor for a single-index set of shape functions.
Definition: shape.h:522
~ShapeWithDeepCopy()
Destructor, clear up the memory allocated by the object.
Definition: shape.h:544
ShapeWithDeepCopy(const unsigned &N, const unsigned &M)
Constructor for a two-index set of shape functions.
Definition: shape.h:525
void operator=(const ShapeWithDeepCopy &old_shape)=delete
Broken assignment operator.
ShapeWithDeepCopy(const ShapeWithDeepCopy &old_shape)
Deep copy constructor.
Definition: shape.h:531
const double & operator[](const unsigned &i) const
Overload the bracket operator (const version)
Definition: shape.h:213
const double & operator()(const unsigned &i) const
Overload the round bracket operator (const version)
Definition: shape.h:231
double * Allocated_storage
Definition: shape.h:86
unsigned nindex1() const
Return the range of index 1 of the shape function object.
Definition: shape.h:259
unsigned Index1
Size of the first index of the shape function.
Definition: shape.h:89
void resize(const unsigned &N, const unsigned &M=1)
Change the size of the storage.
Definition: shape.h:189
void operator=(const Shape &shape)
Definition: shape.h:141
unsigned Index2
Size of the second index of the shape function.
Definition: shape.h:92
double & operator()(const unsigned &i)
Overload the round bracket operator to provide access to values.
Definition: shape.h:222
double & operator()(const unsigned &i, const unsigned &j)
Overload the round bracket operator, allowing for two indices.
Definition: shape.h:240
const double & operator()(const unsigned &i, const unsigned &j) const
Definition: shape.h:250
Shape(const unsigned &N, const unsigned &M)
Constructor for a two-index set of shape functions.
Definition: shape.h:126
double & operator[](const unsigned &i)
Overload the bracket operator to provide access to values.
Definition: shape.h:204
double * Psi
Definition: shape.h:81
void operator=(Shape *const &shape_pt)
Definition: shape.h:162
unsigned nindex2() const
Return the range of index 2 of the shape function object.
Definition: shape.h:265
Shape(const Shape &shape)=delete
Broken copy constructor.
Shape(const unsigned &N)
Constructor for a single-index set of shape functions.
Definition: shape.h:119
void range_check(const unsigned &i, const unsigned &j) const
Private function that checks whether the index is in range.
Definition: shape.h:95
Shape()
Definition: shape.h:137
~Shape()
Destructor, clear up the memory allocated by the object.
Definition: shape.h:182
@ N
Definition: constructor.cpp:22
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:77
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
void dshape< 4 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:1007
void d2shape< 2 >(const double &s, double *DPsi)
Definition: shape.h:957
void dshape< 3 >(const double &s, double *DPsi)
Definition: shape.h:977
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:940
void d2shape< 3 >(const double &s, double *DPsi)
Definition: shape.h:987
void d2shape< 4 >(const double &s, double *DPsi)
Definition: shape.h:1018
void shape< 3 >(const double &s, double *Psi)
1D shape functions specialised to quadratic order (3 Nodes)
Definition: shape.h:967
void dshape(const double &s, double *DPsi)
Definition: shape.h:900
void shape< 4 >(const double &s, double *Psi)
1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:996
void d2shape(const double &s, double *DPsi)
Definition: shape.h:920
void dshape< 2 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:948
void shape(const double &s, double *Psi)
Definition: shape.h:881
void d2shape(const double &s, double *DPsi)
Definition: shape.h:1077
void d2shape< 3 >(const double &s, double *DPsi)
Definition: shape.h:1144
void dshape< 3 >(const double &s, double *DPsi)
Definition: shape.h:1134
void shape< 3 >(const double &s, double *Psi)
1D shape functions specialised to quadratic order (3 Nodes)
Definition: shape.h:1124
void dshape(const double &s, double *DPsi)
Definition: shape.h:1057
void dshape< 2 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:1105
void d2shape< 4 >(const double &s, double *DPsi)
Definition: shape.h:1175
void dshape< 4 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:1164
void shape(const double &s, double *Psi)
Definition: shape.h:1038
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:1097
void d2shape< 2 >(const double &s, double *DPsi)
Definition: shape.h:1114
void shape< 4 >(const double &s, double *Psi)
1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:1153
void dshape< 3 >(const double &s, double *DPsi)
Definition: shape.h:820
void dshape< 4 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:850
void shape< 3 >(const double &s, double *Psi)
1D shape functions specialised to quadratic order (3 Nodes)
Definition: shape.h:810
void dshape(const double &s, double *DPsi)
Definition: shape.h:743
void d2shape< 3 >(const double &s, double *DPsi)
Definition: shape.h:830
void dshape< 2 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:791
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:783
void d2shape(const double &s, double *DPsi)
Definition: shape.h:763
void shape< 4 >(const double &s, double *Psi)
1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:839
void d2shape< 4 >(const double &s, double *DPsi)
Definition: shape.h:861
void shape(const double &s, double *Psi)
Definition: shape.h:724
void d2shape< 2 >(const double &s, double *DPsi)
Definition: shape.h:800
void dshape(const double &s, double DPsi[2][2])
Derivatives of 1D Hermite shape functions.
Definition: shape.h:1206
void d2shape(const double &s, double DPsi[2][2])
Second derivatives of the Hermite shape functions.
Definition: shape.h:1217
void shape(const double &s, double Psi[2][2])
Constructor sets the values of the shape functions at the position s.
Definition: shape.h:1194
void shape(const double &s, double *Psi)
Definition: shape.h:564
void d2shape< 2 >(const double &s, double *DPsi)
Definition: shape.h:625
void dshape< 3 >(const double &s, double *DPsi)
Definition: shape.h:645
void d2shape(const double &s, double *DPsi)
Definition: shape.h:593
void dshape< 2 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:616
void d2shape< 3 >(const double &s, double *DPsi)
Definition: shape.h:656
void d2shape< 4 >(const double &s, double *DPsi)
Definition: shape.h:700
void dshape< 4 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:684
void shape< 3 >(const double &s, double *Psi)
1D shape functions specialised to quadratic order (3 Nodes)
Definition: shape.h:635
void shape< 4 >(const double &s, double *Psi)
1D shape functions specialised to cubic order (4 Nodes)
Definition: shape.h:665
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:608
void dshape(const double &s, double *DPsi)
Definition: shape.h:578
double dlegendre(const unsigned &p, const double &x)
Definition: orthpoly.h:121
const double eps
Definition: orthpoly.h:52
double ddlegendre(const unsigned &p, const double &x)
Definition: orthpoly.h:144
double legendre(const unsigned &p, const double &x)
Definition: orthpoly.h:57
void gll_nodes(const unsigned &Nnode, Vector< double > &x)
Calculates the Gauss Lobatto Legendre abscissas for degree p = NNode-1.
Definition: orthpoly.cc:33
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2