26 #ifndef OOMPH_BIHARMONIC_ELEMENTS_HEADER
27 #define OOMPH_BIHARMONIC_ELEMENTS_HEADER
31 #include <oomph-lib-config.h>
44 #include "../generic/matrices.h"
45 #include "../generic/elements.h"
46 #include "../generic/hermite_elements.h"
54 template<
unsigned DIM>
70 virtual double u(
const unsigned&
n,
const unsigned&
k)
const = 0;
113 void output(std::ostream& outfile,
const unsigned& nplot)
123 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
128 for (
unsigned i = 0;
i <
DIM;
i++)
154 void output(FILE* file_pt,
const unsigned& n_plot)
171 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
176 for (
unsigned i = 0;
i <
DIM;
i++)
184 outfile << dudx[1] <<
" " << -dudx[0] << std::endl;
196 unsigned n_node = this->
nnode();
202 Shape psi(n_node, n_value);
214 for (
unsigned n = 0;
n < n_node;
n++)
216 for (
unsigned k = 0;
k < n_value;
k++)
218 for (
unsigned d = 0; d <
DIM; d++)
229 const unsigned& nplot,
246 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
258 for (
unsigned i = 0;
i <
DIM;
i++)
260 outfile <<
x[
i] <<
" ";
274 const unsigned& nplot,
299 unsigned n_node = this->
nnode();
316 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
319 for (
unsigned i = 0;
i <
DIM;
i++)
343 for (
unsigned i = 0;
i <
DIM;
i++)
345 outfile <<
x[
i] <<
" ";
378 unsigned n_node = this->
nnode();
384 Shape psi(n_node, n_value);
390 for (
unsigned n = 0;
n < n_node;
n++)
392 for (
unsigned k = 0;
k < n_value;
k++)
394 uu +=
u(
n,
k) * psi(
n,
k);
448 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
451 int n_node = this->
nnode();
457 std::pair<unsigned long, unsigned> dof_lookup;
460 for (
int n = 0;
n < n_node;
n++)
463 for (
int k = 0;
k < n_value;
k++)
470 if (local_eqn_number >= 0)
473 dof_lookup.first = this->
eqn_number(local_eqn_number);
474 dof_lookup.second =
k;
475 dof_lookup_list.push_front(dof_lookup);
517 template<
unsigned DIM>
524 template<
unsigned DIM>
530 inline double u(
const unsigned&
n,
const unsigned&
k)
const
559 void output(std::ostream& outfile,
const unsigned& n_plot)
571 void output(FILE* file_pt,
const unsigned& n_plot)
579 const unsigned& nplot,
588 const unsigned& nplot,
603 outfile, exact_soln_pt,
error, norm);
614 outfile, exact_soln_pt, time,
error, norm);
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
biharmonic element class
Definition: biharmonic_elements.h:527
void output(std::ostream &outfile)
Output.
Definition: biharmonic_elements.h:553
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Call the equations-class overloaded unsteady error calculation.
Definition: biharmonic_elements.h:607
unsigned required_nvalue(const unsigned &n) const
Definition: biharmonic_elements.h:546
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
analytic solution wrapper
Definition: biharmonic_elements.h:578
void output(FILE *file_pt)
C-style output.
Definition: biharmonic_elements.h:565
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
computes error
Definition: biharmonic_elements.h:597
void output(FILE *file_pt, const unsigned &n_plot)
C_style output at n_plot points.
Definition: biharmonic_elements.h:571
void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Final override.
Definition: biharmonic_elements.h:587
double u(const unsigned &n, const unsigned &k) const
access function for value, kth dof of node n
Definition: biharmonic_elements.h:530
void output(std::ostream &outfile, const unsigned &n_plot)
output wrapper
Definition: biharmonic_elements.h:559
BiharmonicElement()
Definition: biharmonic_elements.h:538
~BiharmonicElement()
Definition: biharmonic_elements.h:541
Biharmonic Equation Class - contains the equations.
Definition: biharmonic_elements.h:56
void fill_in_contribution_to_residuals(Vector< double > &residuals)
wrapper function, adds contribution to residual
Definition: biharmonic_elements.h:93
void output_fluid_velocity(std::ostream &outfile, const unsigned &nplot)
output fluid velocity field
Definition: biharmonic_elements.h:161
SourceFctPt Source_fct_pt
Pointer to source function:
Definition: biharmonic_elements.h:504
void output(std::ostream &outfile)
Output at default number of plot points.
Definition: biharmonic_elements.h:142
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: biharmonic_elements.h:362
virtual void get_source(const unsigned &ipt, const Vector< double > &x, double &source) const
gets source strength
Definition: biharmonic_elements.h:74
BiharmonicEquations()
Constructor (must initialise the Source_fct_pt to null)
Definition: biharmonic_elements.h:63
static const unsigned d2_dof[]
Definition: biharmonic_elements.h:512
unsigned self_test()
self test wrapper
Definition: biharmonic_elements.h:404
void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: biharmonic_elements.h:273
virtual void fill_in_generic_residual_contribution_biharmonic(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned JFLAG)
Definition: biharmonic_elements.cc:41
void output(std::ostream &outfile, const unsigned &nplot)
output with nplot points
Definition: biharmonic_elements.h:113
~BiharmonicEquations()
Definition: biharmonic_elements.h:66
void interpolated_dudx(const Vector< double > &s, Vector< double > &dudx)
output with nplot points
Definition: biharmonic_elements.h:193
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
computes the error
Definition: biharmonic_elements.h:283
void output(FILE *file_pt)
C-style output.
Definition: biharmonic_elements.h:148
void(* SourceFctPt)(const Vector< double > &x, double &u)
source function type definition
Definition: biharmonic_elements.h:59
void fill_in_contribution_to_jacobian(Vector< double > &residual, DenseMatrix< double > &jacobian)
wrapper function, adds contribution to residual and generic
Definition: biharmonic_elements.h:104
unsigned ndof_types() const
Definition: biharmonic_elements.h:435
double interpolated_u_biharmonic(const Vector< double > &s)
calculates interpolated u at s
Definition: biharmonic_elements.h:372
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
output analytic solution
Definition: biharmonic_elements.h:228
SourceFctPt & source_fct_pt()
Access functions for the source function pointer.
Definition: biharmonic_elements.h:484
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: biharmonic_elements.h:447
unsigned get_d2_dof()
return number of second derivate degrees of freedom
Definition: biharmonic_elements.h:427
virtual double u(const unsigned &n, const unsigned &k) const =0
Vector< int > Local_eqn
Array to hold local eqn numbers: Local_eqn[n] (=-1 for BC)
Definition: biharmonic_elements.h:508
void output(FILE *file_pt, const unsigned &n_plot)
C_style output at n_plot points.
Definition: biharmonic_elements.h:154
SourceFctPt source_fct_pt() const
Access functions for the source function pointers (const version)
Definition: biharmonic_elements.h:490
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
Definition: matrices.h:1271
Definition: elements.h:1313
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output an exact solution over the element.
Definition: elements.h:3104
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
virtual unsigned required_nvalue(const unsigned &n) const
Definition: elements.h:2455
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual void compute_error(FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Calculate the norm of the error and that of the exact solution.
Definition: elements.h:3198
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
virtual unsigned self_test()
Definition: elements.cc:4440
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
unsigned long eqn_number(const unsigned &ieqn_local) const
Definition: elements.h:704
int local_eqn_number(const unsigned long &ieqn_global) const
Definition: elements.h:726
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
double value(const unsigned &i) const
Definition: nodes.cc:2408
Definition: hermite_elements.h:86
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28