biharmonic_problem.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 #ifndef OOMPH_BIHARMONIC_PROBLEM_HEADER
27 #define OOMPH_BIHARMONIC_PROBLEM_HEADER
28 
29 // Config header generated by autoconfig
30 #ifdef HAVE_CONFIG_H
31 #include <oomph-lib-config.h>
32 #endif
33 
34 #ifdef OOMPH_HAS_MPI
35 // mpi headers
36 #include "mpi.h"
37 #endif
38 
39 // Generic C++ headers
40 #include <iostream>
41 #include <math.h>
42 
43 // oomph-lib headers
44 #include "../generic/problem.h"
45 #include "../generic/hijacked_elements.h"
46 #include "../meshes/hermite_element_quad_mesh.template.h"
47 #include "../meshes/hermite_element_quad_mesh.template.cc"
48 #include "biharmonic_elements.h"
50 
51 
52 namespace oomph
53 {
54  //=============================================================================
55  /// Biharmonic Plate Problem Class - for problems where the load can be
56  /// assumed to be acting normal to the surface of the plate and the
57  /// deflections are small relative to the thickness of the plate. Developed
58  /// for the topologically rectangular Hermite Element Mesh. Contains functions
59  /// allowing the following boundary conditions to be applied (on a given
60  /// edge):
61  /// + clamped : u and du/dn imposed
62  /// + simply supported : u and laplacian(u) imposed
63  /// + free : laplacian(u) and dlaplacian(u)/dn imposed
64  //=============================================================================
65  template<unsigned DIM>
66  class BiharmonicProblem : public Problem
67  {
68  public:
69  /// Definition of a dirichlet boundary condition function pointer.
70  /// Takes the position along a boundary (s) in the macro element coordinate
71  /// scheme and returns the value of the boundary condition at that point
72  /// (u).
73  typedef void (*DirichletBCFctPt)(const double& s, double& u);
74 
75 
76  /// Definition of the Source Function.
77  typedef void (*BiharmonicSourceFctPt)(const Vector<double>& x, double& f);
78 
79  /// Constructor
81  {
84  }
85 
86  /// Destructor. Delete the meshes
88  {
89  delete Bulk_element_mesh_pt;
90  delete Face_element_mesh_pt;
91  };
92 
93  /// actions before solve, performs self test
95  {
96 #ifdef PARANOID
97  if (0 == self_test())
98  {
99  oomph_info << "self test passed" << std::endl;
100  }
101  else
102  {
103  oomph_info << "self test failed" << std::endl;
104  }
105 #endif
106  }
107 
108  /// action after solve
110 
111  /// documents the solution, and if an exact solution is provided,
112  /// then the error between the numerical and exact solution is presented
113  void doc_solution(
114  DocInfo& doc_info,
115  FiniteElement::SteadyExactSolutionFctPt exact_soln_pt = 0);
116 
117  /// Access function to the bulk element mesh pt
119  {
120  return Bulk_element_mesh_pt;
121  }
122 
123  protected:
124  /// builds the bulk mesh on a prescribed domain with a node spacing
125  /// defined by spacing fn with n_x by n_y elements
127  const unsigned n_x,
128  const unsigned n_y,
130  HermiteQuadMesh<BiharmonicElement<2>>::MeshSpacingFnPtr spacing_fn = 0)
131  {
132  if (spacing_fn == 0)
133  {
136  n_x, n_y, domain_pt);
137  }
138  else
139  {
142  n_x, n_y, domain_pt, spacing_fn);
143  }
144  // add_sub_mesh(Bulk_element_mesh_pt);
145  // build_global_mesh();
146  }
147 
148 
149  /// Build global mesh and assign equation numbers
151  {
153  if (Face_element_mesh_pt != 0)
154  {
156  }
159  }
160 
161  /// Impose a load to the surface of the plate.
162  /// note : MUST be called before neumann boundary conditions are imposed,
163  /// i.e. a free edge or a simply supported edge with laplacian(u) imposed
164  void set_source_function(const BiharmonicSourceFctPt source_fct_pt)
165  {
166  // number of elements in mesh
167  unsigned n_bulk_element = Bulk_element_mesh_pt->nelement();
168 
169  // loop over bulk elements
170  for (unsigned i = 0; i < n_bulk_element; i++)
171  {
172  // upcast from generalised element to specific element
173  BiharmonicElement<2>* element_pt = dynamic_cast<BiharmonicElement<2>*>(
175 
176  // set the source function pointer
177  element_pt->source_fct_pt() = source_fct_pt;
178  }
179  }
180 
181  /// Imposes the prescribed dirichlet BCs u (u_fn) and
182  /// du/dn (dudn_fn) dirichlet BCs by 'pinning'
183  void set_dirichlet_boundary_condition(const unsigned& b,
184  DirichletBCFctPt u_fn = 0,
185  DirichletBCFctPt dudn_fn = 0);
186 
187  /// Imposes the prescribed Neumann BCs laplacian(u) (flux0_fct_pt)
188  /// and dlaplacian(u)/dn (flux1_fct_pt) with flux edge elements
190  const unsigned& b,
192  BiharmonicFluxElement<2>::FluxFctPt flux1_fct_pt = 0);
193 
194  public:
195  // NOTE: these two private meshes are required for the block
196  // preconditioners.
197 
198  /// Mesh for BiharmonicElement<DIM> only - the block preconditioner
199  /// assemble the global equation number to block number mapping from
200  /// elements in this mesh only
202 
203  /// mesh for face elements
205  };
206 
207 
208  //=============================================================================
209  /// Biharmonic Fluid Problem Class - describes stokes flow in 2D.
210  /// Developed for the topologically rectangular Hermite Element Mesh. Contains
211  /// functions allowing the following boundary conditions to be applied (on a
212  /// given edge):
213  /// + wall : v_n = 0 and v_t = 0 (psi must also be prescribed)
214  /// + traction free : v_t = 0
215  /// + flow : v_n and v_t are prescribed
216  /// NOTE 1 : psi is the stream function
217  /// + fluid velocity normal to boundary v_n = +/- dpsi/dt
218  /// + fluid velocity tangential to boundary v_t = -/+ dpsi/dn
219  /// NOTE 2 : when a solid wall boundary condition is applied to ensure that
220  /// v_n = 0 the the streamfunction psi must also be prescribed (and
221  /// constant)
222  //=============================================================================
223  template<unsigned DIM>
225  {
226  public:
227  /// Definition of a dirichlet boundary condition function pointer.
228  /// Takes the position along a boundary (s) in the macro element coordinate
229  /// scheme and returns the fluid velocity normal (dpsi/dt) to the boundary
230  /// (u[0]) and the fluid velocity tangential (dpsidn) to the boundary
231  /// (u[1]).
232  typedef void (*FluidBCFctPt)(const double& s, Vector<double>& u);
233 
234 
235  /// constructor
237  {
238  // initialise the number of non bulk elements
239  Npoint_element = 0;
240  }
241 
242 
243  /// actions before solve, performs self test
245  {
246 #ifdef PARANOID
247  if (0 == self_test())
248  {
249  oomph_info << "self test passed" << std::endl;
250  }
251  else
252  {
253  oomph_info << "self test failed" << std::endl;
254  }
255 #endif
256  }
257 
258 
259  /// action after solve
261 
262 
263  /// documents the solution, and if an exact solution is provided,
264  /// then the error between the numerical and exact solution is presented
265  void doc_solution(
266  DocInfo& doc_info,
267  FiniteElement::SteadyExactSolutionFctPt exact_soln_pt = 0);
268 
269 
270  protected:
271  /// Imposes a solid boundary on boundary b - no flow into boundary
272  /// or along boundary v_n = 0 and v_t = 0. User must presribe the
273  /// streamfunction psi to ensure dpsi/dt = 0 is imposed at all points on the
274  /// boundary and not just at the nodes
275  void impose_solid_boundary_on_edge(const unsigned& b,
276  const double& psi = 0);
277 
278  /// Impose a traction free edge - i.e. v_t = 0 or dpsi/dn = 0. In
279  /// general dpsi/dn = 0 can only be imposed using equation elements to set
280  /// the DOFs dpsi/ds_n, however in the special case of dt/ds_n = 0, then
281  /// dpsi/ds_n = 0 and can be imposed using pinning - this is handled
282  /// automatically in this function. For a more detailed description of the
283  /// equations see the description of the class
284  /// BiharmonicFluidBoundaryElement
285  void impose_traction_free_edge(const unsigned& b);
286 
287 
288  /// Impose a prescribed fluid flow comprising the velocity normal to
289  /// the boundary (u_imposed_fn[0]) and the velocity tangential to the
290  /// boundary (u_imposed_fn[1])
291  void impose_fluid_flow_on_edge(const unsigned& b,
292  FluidBCFctPt u_imposed_fn);
293 
294 
295  private:
296  // number of non-bulk elements - i.e. biharmonic fluid boundary elements
297  unsigned Npoint_element;
298  };
299 
300 
301  //=============================================================================
302  /// Point equation element used to impose the traction free edge (i.e.
303  /// du/dn = 0) on the boundary when dt/ds_n != 0. The following equation is
304  /// implemented : du/ds_n = dt/ds_n * ds_t/dt * du/dt.
305  /// The bulk biharmonic elements on the boundary must be hijackable and the
306  /// du/ds_n and d2u/ds_nds_t boundary DOFs hijacked when these elements are
307  /// applied. At any node where dt/ds_n = 0 we can impose du/ds_n = 0 and
308  /// d2u/ds_nds_t = 0 using pinning - see
309  /// BiharmonicFluidProblem::impose_traction_free_edge()
310  //=============================================================================
312  {
313  public:
314  // constructor
315  BiharmonicFluidBoundaryElement(Node* node_pt, const unsigned s_fixed_index)
316  {
317  // set the node pt
318  this->node_pt(0) = node_pt;
319 
320  // store fixed index on the boundary
321  S_fixed_index = s_fixed_index;
322  }
323 
324  /// Output function -- does nothing
325  void output(std::ostream& outfile) {}
326 
327 
328  /// Output function -- does nothing
329  void output(std::ostream& outfile, const unsigned& n_plot) {}
330 
331 
332  /// Output function -- does nothing
333  void output_fluid_velocity(std::ostream& outfile, const unsigned& n_plot) {}
334 
335 
336  /// C-style output function -- does nothing
337  void output(FILE* file_pt) {}
338 
339 
340  /// C-style output function -- does nothing
341  void output(FILE* file_pt, const unsigned& n_plot) {}
342 
343 
344  /// compute_error -- does nothing
345  void compute_error(std::ostream& outfile,
347  double& error,
348  double& norm)
349  {
350  }
351 
352 
353  /// Compute the elemental residual vector - wrapper function called
354  /// by get_residuals in GeneralisedElement
356  {
357  // create a dummy matrix
358  DenseDoubleMatrix dummy(1);
359 
360  // call the generic residuals functions with flag set to zero
362  residuals, dummy, 0);
363  }
364 
365 
366  /// Compute the elemental residual vector and jacobian matrix -
367  /// wrapper function called by get_jacobian in GeneralisedElement
369  DenseMatrix<double>& jacobian)
370  {
371  // call generic routine with flag set to 1
373  residuals, jacobian, 1);
374  }
375 
376 
377  /// Computes the elemental residual vector and the elemental jacobian
378  /// matrix if JFLAG = 0
379  /// Imposes the equations : du/ds_n = dt/ds_n * ds_t/dt * du/dt
381  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned JFLAG);
382 
383  private:
384  // fixed local coordinate index on boundary
385  unsigned S_fixed_index;
386  };
387 
388 
389 } // namespace oomph
390 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
biharmonic element class
SourceFctPt & source_fct_pt()
Access functions for the source function pointer.
Point equation element used to impose the traction free edge (i.e. du/dn = 0) on the boundary when dt...
void output(std::ostream &outfile)
Output function – does nothing.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – does nothing.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the elemental residual vector - wrapper function called by get_residuals in GeneralisedElemen...
void output(FILE *file_pt)
C-style output function – does nothing.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – does nothing.
BiharmonicFluidBoundaryElement(Node *node_pt, const unsigned s_fixed_index)
void output_fluid_velocity(std::ostream &outfile, const unsigned &n_plot)
Output function – does nothing.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the elemental residual vector and jacobian matrix - wrapper function called by get_jacobian i...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
compute_error – does nothing
virtual void fill_in_generic_residual_contribution_biharmonic_boundary(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned JFLAG)
Computes the elemental residual vector and the elemental jacobian matrix if JFLAG = 0 Imposes the equ...
Biharmonic Fluid Problem Class - describes stokes flow in 2D. Developed for the topologically rectang...
void doc_solution(DocInfo &doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt=0)
documents the solution, and if an exact solution is provided, then the error between the numerical an...
void impose_fluid_flow_on_edge(const unsigned &b, FluidBCFctPt u_imposed_fn)
Impose a prescribed fluid flow comprising the velocity normal to the boundary (u_imposed_fn[0]) and t...
void impose_solid_boundary_on_edge(const unsigned &b, const double &psi=0)
Imposes a solid boundary on boundary b - no flow into boundary or along boundary v_n = 0 and v_t = 0....
void impose_traction_free_edge(const unsigned &b)
Impose a traction free edge - i.e. v_t = 0 or dpsi/dn = 0. In general dpsi/dn = 0 can only be imposed...
void actions_after_newton_solve()
action after solve
void(* FluidBCFctPt)(const double &s, Vector< double > &u)
Definition of a dirichlet boundary condition function pointer. Takes the position along a boundary (s...
void actions_before_newton_solve()
actions before solve, performs self test
Biharmonic Plate Problem Class - for problems where the load can be assumed to be acting normal to th...
void(* BiharmonicSourceFctPt)(const Vector< double > &x, double &f)
Definition of the Source Function.
void actions_before_newton_solve()
actions before solve, performs self test
void set_neumann_boundary_condition(const unsigned &b, BiharmonicFluxElement< 2 >::FluxFctPt flux0_fct_pt, BiharmonicFluxElement< 2 >::FluxFctPt flux1_fct_pt=0)
Imposes the prescribed Neumann BCs laplacian(u) (flux0_fct_pt) and dlaplacian(u)/dn (flux1_fct_pt) wi...
void build_bulk_mesh(const unsigned n_x, const unsigned n_y, TopologicallyRectangularDomain *domain_pt, HermiteQuadMesh< BiharmonicElement< 2 >>::MeshSpacingFnPtr spacing_fn=0)
builds the bulk mesh on a prescribed domain with a node spacing defined by spacing fn with n_x by n_y...
void build_global_mesh_and_assign_eqn_numbers()
Build global mesh and assign equation numbers.
Mesh * bulk_element_mesh_pt()
Access function to the bulk element mesh pt.
void doc_solution(DocInfo &doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt=0)
documents the solution, and if an exact solution is provided, then the error between the numerical an...
Mesh * Bulk_element_mesh_pt
Mesh for BiharmonicElement<DIM> only - the block preconditioner assemble the global equation number t...
virtual ~BiharmonicProblem()
Destructor. Delete the meshes.
void(* DirichletBCFctPt)(const double &s, double &u)
Definition of a dirichlet boundary condition function pointer. Takes the position along a boundary (s...
void set_dirichlet_boundary_condition(const unsigned &b, DirichletBCFctPt u_fn=0, DirichletBCFctPt dudn_fn=0)
Imposes the prescribed dirichlet BCs u (u_fn) and du/dn (dudn_fn) dirichlet BCs by 'pinning'.
void actions_after_newton_solve()
action after solve
void set_source_function(const BiharmonicSourceFctPt source_fct_pt)
Impose a load to the surface of the plate. note : MUST be called before neumann boundary conditions a...
Mesh * Face_element_mesh_pt
mesh for face elements
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
Definition: matrices.h:1271
Information for documentation of results: Directory and file number to enable output in the form RESL...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
A two dimensional Hermite bicubic element quadrilateral mesh for a topologically rectangular domain....
A general mesh class.
Definition: mesh.h:67
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:906
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3443
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: problem.h:154
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Add a submesh to the problem and return its number, i, by which it can be accessed via mesh_pt(i).
Definition: problem.h:1350
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Assign all equation numbers for problem: Deals with global data (= data that isn't attached to any el...
Definition: problem.cc:2076
void build_global_mesh()
Build the global mesh by combining the all the submeshes. Note: The nodes boundary information refers...
Definition: problem.cc:1580
unsigned self_test()
Self-test: Check meshes and global data. Return 0 for OK.
Definition: problem.cc:13339
Topologically Rectangular Domain - a domain dexcribing a topologically rectangular problem - primaril...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...