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-2022 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  //=============================================================================
64  //=============================================================================
65  template<unsigned DIM>
66  class BiharmonicProblem : public Problem
67  {
68  public:
73  typedef void (*DirichletBCFctPt)(const double& s, double& u);
74 
75 
77  typedef void (*BiharmonicSourceFctPt)(const Vector<double>& x, double& f);
78 
81  {
84  }
85 
88  {
89  delete Bulk_element_mesh_pt;
90  delete Face_element_mesh_pt;
91  };
92 
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 
110 
113  void doc_solution(
114  DocInfo& doc_info,
115  FiniteElement::SteadyExactSolutionFctPt exact_soln_pt = 0);
116 
119  {
120  return Bulk_element_mesh_pt;
121  }
122 
123  protected:
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 
151  {
153  if (Face_element_mesh_pt != 0)
154  {
156  }
159  }
160 
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 
183  void set_dirichlet_boundary_condition(const unsigned& b,
184  DirichletBCFctPt u_fn = 0,
185  DirichletBCFctPt dudn_fn = 0);
186 
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 
202 
205  };
206 
207 
208  //=============================================================================
222  //=============================================================================
223  template<unsigned DIM>
225  {
226  public:
232  typedef void (*FluidBCFctPt)(const double& s, Vector<double>& u);
233 
234 
237  {
238  // initialise the number of non bulk elements
239  Npoint_element = 0;
240  }
241 
242 
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 
261 
262 
265  void doc_solution(
266  DocInfo& doc_info,
267  FiniteElement::SteadyExactSolutionFctPt exact_soln_pt = 0);
268 
269 
270  protected:
275  void impose_solid_boundary_on_edge(const unsigned& b,
276  const double& psi = 0);
277 
285  void impose_traction_free_edge(const unsigned& b);
286 
287 
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  //=============================================================================
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 
325  void output(std::ostream& outfile) {}
326 
327 
329  void output(std::ostream& outfile, const unsigned& n_plot) {}
330 
331 
333  void output_fluid_velocity(std::ostream& outfile, const unsigned& n_plot) {}
334 
335 
337  void output(FILE* file_pt) {}
338 
339 
341  void output(FILE* file_pt, const unsigned& n_plot) {}
342 
343 
345  void compute_error(std::ostream& outfile,
347  double& error,
348  double& norm)
349  {
350  }
351 
352 
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 
369  DenseMatrix<double>& jacobian)
370  {
371  // call generic routine with flag set to 1
373  residuals, jacobian, 1);
374  }
375 
376 
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
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Scalar * b
Definition: benchVecAdd.cpp:17
biharmonic element class
Definition: biharmonic_elements.h:527
SourceFctPt & source_fct_pt()
Access functions for the source function pointer.
Definition: biharmonic_elements.h:484
Definition: biharmonic_problem.h:312
void output(std::ostream &outfile)
Output function – does nothing.
Definition: biharmonic_problem.h:325
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – does nothing.
Definition: biharmonic_problem.h:341
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Definition: biharmonic_problem.h:355
unsigned S_fixed_index
Definition: biharmonic_problem.h:385
void output(FILE *file_pt)
C-style output function – does nothing.
Definition: biharmonic_problem.h:337
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – does nothing.
Definition: biharmonic_problem.h:329
BiharmonicFluidBoundaryElement(Node *node_pt, const unsigned s_fixed_index)
Definition: biharmonic_problem.h:315
void output_fluid_velocity(std::ostream &outfile, const unsigned &n_plot)
Output function – does nothing.
Definition: biharmonic_problem.h:333
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: biharmonic_problem.h:368
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
compute_error – does nothing
Definition: biharmonic_problem.h:345
virtual void fill_in_generic_residual_contribution_biharmonic_boundary(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned JFLAG)
Definition: biharmonic_problem.cc:391
Definition: biharmonic_problem.h:225
void doc_solution(DocInfo &doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt=0)
Definition: biharmonic_problem.cc:810
void impose_fluid_flow_on_edge(const unsigned &b, FluidBCFctPt u_imposed_fn)
Definition: biharmonic_problem.cc:658
void impose_solid_boundary_on_edge(const unsigned &b, const double &psi=0)
Definition: biharmonic_problem.cc:457
void impose_traction_free_edge(const unsigned &b)
Definition: biharmonic_problem.cc:489
BiharmonicFluidProblem()
constructor
Definition: biharmonic_problem.h:236
unsigned Npoint_element
Definition: biharmonic_problem.h:297
void actions_after_newton_solve()
action after solve
Definition: biharmonic_problem.h:260
void(* FluidBCFctPt)(const double &s, Vector< double > &u)
Definition: biharmonic_problem.h:232
void actions_before_newton_solve()
actions before solve, performs self test
Definition: biharmonic_problem.h:244
Biharmonic Flux Element.
Definition: biharmonic_flux_elements.h:84
Definition: biharmonic_problem.h:67
void(* BiharmonicSourceFctPt)(const Vector< double > &x, double &f)
Definition of the Source Function.
Definition: biharmonic_problem.h:77
void actions_before_newton_solve()
actions before solve, performs self test
Definition: biharmonic_problem.h:94
void set_neumann_boundary_condition(const unsigned &b, BiharmonicFluxElement< 2 >::FluxFctPt flux0_fct_pt, BiharmonicFluxElement< 2 >::FluxFctPt flux1_fct_pt=0)
Definition: biharmonic_problem.cc:294
void build_bulk_mesh(const unsigned n_x, const unsigned n_y, TopologicallyRectangularDomain *domain_pt, HermiteQuadMesh< BiharmonicElement< 2 >>::MeshSpacingFnPtr spacing_fn=0)
Definition: biharmonic_problem.h:126
void build_global_mesh_and_assign_eqn_numbers()
Build global mesh and assign equation numbers.
Definition: biharmonic_problem.h:150
Mesh * bulk_element_mesh_pt()
Access function to the bulk element mesh pt.
Definition: biharmonic_problem.h:118
void doc_solution(DocInfo &doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt=0)
Definition: biharmonic_problem.cc:341
Mesh * Bulk_element_mesh_pt
Definition: biharmonic_problem.h:201
BiharmonicProblem()
Constructor.
Definition: biharmonic_problem.h:80
virtual ~BiharmonicProblem()
Destructor. Delete the meshes.
Definition: biharmonic_problem.h:87
void(* DirichletBCFctPt)(const double &s, double &u)
Definition: biharmonic_problem.h:73
void set_dirichlet_boundary_condition(const unsigned &b, DirichletBCFctPt u_fn=0, DirichletBCFctPt dudn_fn=0)
Definition: biharmonic_problem.cc:42
void actions_after_newton_solve()
action after solve
Definition: biharmonic_problem.h:109
void set_source_function(const BiharmonicSourceFctPt source_fct_pt)
Definition: biharmonic_problem.h:164
Mesh * Face_element_mesh_pt
mesh for face elements
Definition: biharmonic_problem.h:204
Definition: matrices.h:1271
Definition: oomph_utilities.h:499
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Definition: hermite_element_quad_mesh.template.h:53
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
Definition: nodes.h:906
Definition: elements.h:3439
Definition: problem.h:151
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
unsigned self_test()
Self-test: Check meshes and global data. Return 0 for OK.
Definition: problem.cc:13276
void build_global_mesh()
Definition: problem.cc:1493
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
Definition: topologically_rectangular_domain.h:51
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
int error
Definition: calibrate.py:297
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
GeneralisedAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: gen_axisym_advection_diffusion_elements.h:229
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
list x
Definition: plotDoE.py:28