biharmonic_problem.cc
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 // Config header generated by autoconfig
27 #ifdef HAVE_CONFIG_H
28 #include <oomph-lib-config.h>
29 #endif
30 
31 // oomph-lib includes
32 #include "biharmonic_problem.h"
33 
34 
35 namespace oomph
36 {
37  //=============================================================================
38  /// Impose a Clamped Edge. Imposes the prescribed dirichlet BCs u and
39  /// du/dn dirichlet BCs by pinning
40  //=============================================================================
41  template<unsigned DIM>
43  const unsigned& b, DirichletBCFctPt u_fn, DirichletBCFctPt dudn_fn)
44  {
45  // number of nodes on boundary b
46  unsigned n_node = Bulk_element_mesh_pt->nboundary_node(b);
47 
48  // fixed faced index for boundary
49  int face_index = Bulk_element_mesh_pt->face_index_at_boundary(b, 0);
50 
51  // Need to get the s_fixed_index
52  unsigned s_fixed_index = 0;
53  switch (face_index)
54  {
55  case -1:
56  case 1:
57  s_fixed_index = 0;
58  break;
59 
60  case -2:
61  case 2:
62  s_fixed_index = 1;
63  break;
64 
65  default:
66  throw OomphLibError("Face Index not +/-1 or +/-2: Need 2D QElements",
67  OOMPH_CURRENT_FUNCTION,
68  OOMPH_EXCEPTION_LOCATION);
69  }
70 
71  // node position along edge b in macro element boundary representation
72  // [-1,1]
73  Vector<double> s(1);
74 
75  // Set the edge sign
76  int edge_sign = 0;
77  switch (face_index)
78  {
79  case -1:
80  case 2:
81  edge_sign = -1;
82  break;
83 
84  case 1:
85  case -2:
86  edge_sign = 1;
87  break;
88  }
89 
90  // finite difference step
91  const double h = 10e-8;
92 
93  // node position along edge b in macro element boundary representation
94  // [-1,1]
95  Vector<double> m(2);
96 
97  // if u is prescribed then impose
98  if (u_fn != 0)
99  {
100  // loop over nodes on boundary b
101  for (unsigned n = 0; n < n_node; n++)
102  {
103  // find node position along edge [-1,1] in macro element representation
104  Bulk_element_mesh_pt->boundary_node_pt(b, n)
105  ->get_coordinates_on_boundary(b, m);
106 
107  // get u at node
108  double u;
109  (*u_fn)(m[0], u);
110 
111  // finite difference is used to compute dudm_t
112  // u0 and u1 store values of u left/below and right/above node n
113  double u_L, u_R;
114 
115  // if left/lower corner node
116  if (n == 0)
117  {
118  (*u_fn)(m[0], u_L);
119  (*u_fn)(m[0] + h, u_R);
120  }
121  // if right/upper corner node
122  else if (n == n_node - 1)
123  {
124  (*u_fn)(m[0] - h, u_L);
125  (*u_fn)(m[0], u_R);
126  }
127  // if other node
128  else
129  {
130  (*u_fn)(m[0] - 0.5 * h, u_L);
131  (*u_fn)(m[0] + 0.5 * h, u_R);
132  }
133 
134  // compute dudm_t
135  double dudm_t = (u_R - u_L) / h;
136 
137  // compute duds_t
138  double duds_t = m[1] * dudm_t;
139 
140  // pin and set u type dof
141  Bulk_element_mesh_pt->boundary_node_pt(b, n)->pin(0);
142  Bulk_element_mesh_pt->boundary_node_pt(b, n)->set_value(0, u);
143 
144  // pin and set duds_t type degree of freedom
145  Bulk_element_mesh_pt->boundary_node_pt(b, n)->pin(2 - s_fixed_index);
146  Bulk_element_mesh_pt->boundary_node_pt(b, n)->set_value(
147  2 - s_fixed_index, duds_t);
148  }
149  }
150 
151  // if dudn is prescribed then impose
152  if (dudn_fn != 0)
153  {
154  // loop over nodes on boundary b
155  for (unsigned n = 0; n < n_node; n++)
156  {
157  // vectors for dx_i/ds_n and dx_i/ds_t
158  Vector<double> dxds_n(2);
159  Vector<double> dxds_t(2);
160  dxds_n[0] = Bulk_element_mesh_pt->boundary_node_pt(b, n)->x_gen(
161  1 + s_fixed_index, 0);
162  dxds_n[1] = Bulk_element_mesh_pt->boundary_node_pt(b, n)->x_gen(
163  1 + s_fixed_index, 1);
164  dxds_t[0] = Bulk_element_mesh_pt->boundary_node_pt(b, n)->x_gen(
165  2 - s_fixed_index, 0);
166  dxds_t[1] = Bulk_element_mesh_pt->boundary_node_pt(b, n)->x_gen(
167  2 - s_fixed_index, 1);
168 
169  // vector for d2xi/ds_n ds_t
170  Vector<double> d2xds_nds_t(2);
171  d2xds_nds_t[0] =
172  Bulk_element_mesh_pt->boundary_node_pt(b, n)->x_gen(3, 0);
173  d2xds_nds_t[1] =
174  Bulk_element_mesh_pt->boundary_node_pt(b, n)->x_gen(3, 1);
175 
176  // compute dn/ds_n
177  double dnds_n =
178  ((dxds_n[0] * dxds_t[1]) - (dxds_n[1] * dxds_t[0])) /
179  (sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]) * edge_sign);
180 
181  // compute dt/ds_n
182  double dtds_n = ((dxds_n[0] * dxds_t[0]) + (dxds_n[1] * dxds_t[1])) /
183  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
184 
185  // compute dt/ds_t
186  double dtds_t = sqrt(pow(dxds_t[0], 2) + pow(dxds_t[1], 2));
187 
188  // compute ds_n/dn
189  double ds_ndn =
190  -1.0 * (edge_sign * sqrt(pow(dxds_t[0], 2) + pow(dxds_t[1], 2))) /
191  (dxds_t[0] * dxds_n[1] - dxds_n[0] * dxds_t[1]);
192 
193  // compute ds_t/dt
194  double ds_tdt = 1 / sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
195 
196  // compute d2t/ds_nds_t
197  double d2tds_nds_t =
198  (dxds_t[0] * d2xds_nds_t[0] + dxds_t[1] * d2xds_nds_t[1]) /
199  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
200 
201  // compute d2s_t/ds_ndt
202  double d2s_tds_ndt =
203  (dxds_t[0] * d2xds_nds_t[0] + dxds_t[1] * d2xds_nds_t[1]) /
204  pow(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1], 3.0 / 2.0);
205 
206  // get m_t and dm_t/ds_t for this node
207  Vector<double> m_N(2);
208  Bulk_element_mesh_pt->boundary_node_pt(b, n)
209  ->get_coordinates_on_boundary(b, m_N);
210 
211  // compute d2u/dm_t2 and d/dm_t(dudn) by finite difference
212  double d2udm_t2 = 0;
213  double ddm_tdudn = 0;
214  double u_2L, u_N, u_2R, dudn_L, dudn_R;
215  // evaluate u_fn and dudn_fn for current node
216  if (n == 0)
217  {
218  (*u_fn)(m_N[0], u_2L);
219  (*u_fn)(m_N[0] + h, u_N);
220  (*u_fn)(m_N[0] - 2 * h, u_2R);
221  (*dudn_fn)(m_N[0], dudn_L);
222  (*dudn_fn)(m_N[0] + h, dudn_R);
223  }
224  else if (n == n_node - 1)
225  {
226  (*u_fn)(m_N[0] - 2 * h, u_2L);
227  (*u_fn)(m_N[0] - h, u_N);
228  (*u_fn)(m_N[0], u_2R);
229  (*dudn_fn)(m_N[0] - h, dudn_L);
230  (*dudn_fn)(m_N[0], dudn_R);
231  }
232  else
233  {
234  (*u_fn)(m_N[0] - h, u_2L);
235  u_N = Bulk_element_mesh_pt->boundary_node_pt(b, n)->value(0);
236  (*u_fn)(m_N[0] + h, u_2R);
237  (*dudn_fn)(m_N[0] - 0.5 * h, dudn_L);
238  (*dudn_fn)(m_N[0] + 0.5 * h, dudn_R);
239  }
240  // compute
241  d2udm_t2 = (u_2L + u_2R - 2 * u_N) / (h * h);
242  ddm_tdudn = (dudn_R - dudn_L) / h;
243 
244  // get dudn at the node
245  double dudn;
246  (*dudn_fn)(m_N[0], dudn);
247 
248  // compute d2u/ds_t2
249  double d2uds_t2 = m_N[1] * m_N[1] * d2udm_t2;
250 
251  // get duds_t
252  double duds_t = Bulk_element_mesh_pt->boundary_node_pt(b, n)->value(
253  2 - s_fixed_index);
254 
255  // get du/dt
256  double dudt = ds_tdt * duds_t;
257 
258  // compute d2u/dndt
259  double d2udndt = ds_tdt = dtds_t * m_N[1] * ddm_tdudn;
260 
261  // compute d2udt2
262  double dtds_nd2udt2 =
263  edge_sign * (dxds_t[0] * dxds_n[1] - dxds_n[0] * dxds_t[1]) *
264  (ds_tdt *
265  (d2udndt - ds_ndn * (d2s_tds_ndt * dudt + ds_tdt * d2uds_t2)));
266 
267  // compute dds_n(dudt)
268  double dds_ndudt = dtds_nd2udt2 + dnds_n * d2udndt;
269 
270  // compute du/ds_n
271  double duds_n = dnds_n * dudn + dtds_n * ds_tdt * duds_t;
272 
273  // compute d2u/ds_nds_t
274  double d2uds_nds_t = d2tds_nds_t * dudt + dtds_t * dds_ndudt;
275 
276  // pin du/ds_n dof and set value
277  Bulk_element_mesh_pt->boundary_node_pt(b, n)->pin(1 + s_fixed_index);
278  Bulk_element_mesh_pt->boundary_node_pt(b, n)->set_value(
279  1 + s_fixed_index, duds_n);
280 
281  // pin d2u/ds_nds_t dof and set value
282  Bulk_element_mesh_pt->boundary_node_pt(b, n)->pin(3);
283  Bulk_element_mesh_pt->boundary_node_pt(b, n)->set_value(3, d2uds_nds_t);
284  }
285  }
286  }
287 
288 
289  //=============================================================================
290  /// Imposes a 'free' edge. Imposes the prescribed Neumann BCs
291  /// laplacian(u) and dlaplacian(u)/dn with flux edge elements
292  //=============================================================================
293  template<unsigned DIM>
295  const unsigned& b,
298  {
299  // if the face element mesh pt does not exist then build it
300  if (Face_element_mesh_pt == 0)
301  {
302  Face_element_mesh_pt = new Mesh();
303  }
304 
305  // How many bulk elements are adjacent to boundary b?
306  unsigned n_element = Bulk_element_mesh_pt->nboundary_element(b);
307 
308  // Loop over the bulk elements adjacent to boundary b?
309  for (unsigned e = 0; e < n_element; e++)
310  {
311  // Get pointer to the bulk element that is adjacent to boundary b
312  FiniteElement* bulk_elem_pt = dynamic_cast<FiniteElement*>(
313  Bulk_element_mesh_pt->boundary_element_pt(b, e));
314 
315  // What is the face index along the boundary
316  int face_index = Bulk_element_mesh_pt->face_index_at_boundary(b, e);
317 
318  // Build the corresponding prescribed-flux element
319  BiharmonicFluxElement<2>* flux_element_pt =
320  new BiharmonicFluxElement<2>(bulk_elem_pt, face_index, b);
321 
322 
323  // pass the flux BC pointers to the flux elements
324  flux_element_pt->flux0_fct_pt() = flux0_fct_pt;
325  if (flux1_fct_pt != 0)
326  {
327  flux_element_pt->flux1_fct_pt() = flux1_fct_pt;
328  }
329 
330  // Add the prescribed-flux element to the mesh
331  Face_element_mesh_pt->add_element_pt(flux_element_pt);
332  }
333  }
334 
335 
336  //=============================================================================
337  /// documents the solution, and if an exact solution is provided, then
338  /// the error between the numerical and exact solution is presented
339  //=============================================================================
340  template<unsigned DIM>
342  DocInfo& doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
343  {
344  std::ofstream some_file;
345  std::ostringstream filename;
346 
347  // Number of plot points: npts x npts
348  unsigned npts = 5;
349 
350  // Output solution
351  filename << doc_info.directory() << "/soln_" << doc_info.number() << ".dat";
352  some_file.open(filename.str().c_str());
353  Bulk_element_mesh_pt->output(some_file, npts);
354  some_file.close();
355 
356  // if the exact solution is not provided
357  if (exact_soln_pt != 0)
358  {
359  // Output exact solution
360  filename.str("");
361  filename << doc_info.directory() << "/exact_soln_" << doc_info.number()
362  << ".dat";
363  some_file.open(filename.str().c_str());
364  Bulk_element_mesh_pt->output_fct(some_file, npts, exact_soln_pt);
365  some_file.close();
366 
367  // Doc error and return of the square of the L2 error
368  double error, norm;
369  filename.str("");
370  filename << doc_info.directory() << "/error_" << doc_info.number()
371  << ".dat";
372  some_file.open(filename.str().c_str());
373  Bulk_element_mesh_pt->compute_error(
374  some_file, exact_soln_pt, error, norm);
375  some_file.close();
376 
377  // Doc L2 error and norm of solution
378  oomph_info << "\nNorm of error : " << sqrt(error) << std::endl;
379  oomph_info << "Norm of solution: " << sqrt(norm) << std::endl
380  << std::endl;
381  }
382  }
383 
384 
385  //=============================================================================
386  /// Computes the elemental residual vector and the elemental jacobian
387  /// matrix if JFLAG = 0
388  /// Imposes the equations : du/ds_n = dt/ds_n * ds_t/dt * du/dt
389  //=============================================================================
392  Vector<double>& residual, DenseMatrix<double>& jacobian, unsigned JFLAG)
393  {
394  // dof # corresponding to d/ds_n
395  unsigned k_normal = 1 + S_fixed_index;
396 
397  // dof # corresponding to d/ds_t
398  unsigned k_tangential = 2 - S_fixed_index;
399 
400  // vectors for dx_i/ds_n and dx_i/ds_t
401  Vector<double> dxds_n(2);
402  Vector<double> dxds_t(2);
403  dxds_n[0] = this->node_pt(0)->x_gen(1 + S_fixed_index, 0);
404  dxds_n[1] = this->node_pt(0)->x_gen(1 + S_fixed_index, 1);
405  dxds_t[0] = this->node_pt(0)->x_gen(2 - S_fixed_index, 0);
406  dxds_t[1] = this->node_pt(0)->x_gen(2 - S_fixed_index, 1);
407 
408  // vector for d2xi/ds_n ds_t
409  Vector<double> d2xds_nds_t(2);
410  d2xds_nds_t[0] = this->node_pt(0)->x_gen(3, 0);
411  d2xds_nds_t[1] = this->node_pt(0)->x_gen(3, 1);
412 
413  // compute dt/ds_n
414  double dtds_n = ((dxds_n[0] * dxds_t[0]) + (dxds_n[1] * dxds_t[1])) /
415  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
416 
417  // compute ds_t/dt
418  double ds_tdt = 1 / sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
419 
420  // determine the local equation number
421  int local_eqn_number = this->nodal_local_eqn(0, k_normal);
422 
423  // if local equation number equal to -1 then its a boundary node(pinned)
424  if (local_eqn_number >= 0)
425  {
426  // additions to residual for duds_n
427  residual[local_eqn_number] +=
428  (this->node_pt(0)->value(k_normal) -
429  dtds_n * ds_tdt * this->node_pt(0)->value(k_tangential));
430 
431  // if contributions to jacobian are required
432  if (JFLAG == 1)
433  {
434  // additions to jacobian for du/ds_n
435  int local_dof_number = this->nodal_local_eqn(0, k_normal);
436  if (local_dof_number >= 0)
437  {
438  jacobian(local_eqn_number, local_dof_number) += 1.0;
439  }
440  local_dof_number = this->nodal_local_eqn(0, k_tangential);
441  if (local_dof_number >= 0)
442  {
443  jacobian(local_eqn_number, local_dof_number) -= dtds_n * ds_tdt;
444  }
445  }
446  }
447  }
448 
449 
450  //=============================================================================
451  /// Imposes a solid boundary - no flow into boundary or along boundary
452  /// v_n = 0 and v_t = 0. User must presribe the streamfunction psi to ensure
453  /// dpsi/dt = 0 is imposed at all points on the boundary and not just at the
454  /// nodes
455  //=============================================================================
456  template<unsigned DIM>
458  const unsigned& b, const double& psi)
459  {
460  // number of nodes on boundary b
461  unsigned n_node = mesh_pt()->nboundary_node(b);
462 
463  // loop over nodes on boundary b
464  for (unsigned n = 0; n < n_node; n++)
465  {
466  // loop over DOFs to be pinned du/ds_n, du/ds_t and d2u/ds_nds_t
467  for (unsigned k = 1; k < 4; k++)
468  {
469  // pin and zero DOF k
470  mesh_pt()->boundary_node_pt(b, n)->pin(k);
471  mesh_pt()->boundary_node_pt(b, n)->set_value(k, 0.0);
472  }
473  mesh_pt()->boundary_node_pt(b, n)->pin(0);
474  mesh_pt()->boundary_node_pt(b, n)->set_value(0, psi);
475  }
476  }
477 
478 
479  //=============================================================================
480  /// Impose a traction free edge - i.e. v_t = 0 or dpsi/dn = 0. In
481  /// general dpsi/dn = 0 can only be imposed using equation elements to set
482  /// the DOFs dpsi/ds_n, however in the special case of dt/ds_n = 0, then
483  /// dpsi/ds_n = 0 and can be imposed using pinning - this is handled
484  /// automatically in this function. For a more detailed description of the
485  /// equations see the description of the class :
486  /// BiharmonicFluidBoundaryElement
487  //=============================================================================
488  template<unsigned DIM>
490  {
491  // fixed faced index for boundary
492  int face_index = mesh_pt()->face_index_at_boundary(b, 0);
493 
494  // Need to get the s_fixed_index
495  unsigned s_fixed_index = 0;
496  switch (face_index)
497  {
498  case -1:
499  case 1:
500  s_fixed_index = 0;
501  break;
502 
503  case -2:
504  case 2:
505  s_fixed_index = 1;
506  break;
507 
508  default:
509  throw OomphLibError("Face Index not +/-1 or +/-2: Need 2D QElements",
510  OOMPH_CURRENT_FUNCTION,
511  OOMPH_EXCEPTION_LOCATION);
512  }
513 
514  // create a point to a hijacked biharmonic element
515  Hijacked<BiharmonicElement<2>>* hijacked_element_pt;
516 
517  // vectors for dx_i/ds_n and dx_i/ds_t
518  Vector<double> dxds_n(2);
519  Vector<double> dxds_t(2);
520 
521  // number of nodes along edge
522  unsigned n_node = mesh_pt()->nboundary_node(b);
523 
524  // loop over boudnary nodes
525  for (unsigned n = 0; n < n_node; n++)
526  {
527  // get dx_i/ds_t and dx_i/ds_n at node n
528  dxds_n[0] =
529  mesh_pt()->boundary_node_pt(b, n)->x_gen(1 + s_fixed_index, 0);
530  dxds_n[1] =
531  mesh_pt()->boundary_node_pt(b, n)->x_gen(1 + s_fixed_index, 1);
532  dxds_t[0] =
533  mesh_pt()->boundary_node_pt(b, n)->x_gen(2 - s_fixed_index, 0);
534  dxds_t[1] =
535  mesh_pt()->boundary_node_pt(b, n)->x_gen(2 - s_fixed_index, 1);
536 
537  // compute dt/ds_n
538  double dtds_n = ((dxds_n[0] * dxds_t[0]) + (dxds_n[1] * dxds_t[1])) /
539  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
540 
541  // if dt/ds_n = 0 we can impose the traction free edge at this node by
542  // pinning dpsi/ds_n = 0, otherwise the equation elements
543  // BiharmonicFluidBoundaryElement are used
544  if (dtds_n == 0.0)
545  {
546  // pin dpsi/dn
547  mesh_pt()->boundary_node_pt(b, n)->pin(1 + s_fixed_index);
548  mesh_pt()->boundary_node_pt(b, n)->set_value(1 + s_fixed_index, 0.0);
549  }
550 
551  // otherwise impose equation elements
552  else
553  {
554  // hijack DOFs in element either side of node
555  // boundary 0
556  if (b == 0)
557  {
558  // hijack DOFs in left element
559  if (n > 0)
560  {
561  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
562  mesh_pt()->boundary_element_pt(b, n - 1));
563  delete hijacked_element_pt->hijack_nodal_value(1,
564  1 + s_fixed_index);
565  }
566  // hijack DOFs in right element
567  if (n < (n_node - 1))
568  {
569  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
570  mesh_pt()->boundary_element_pt(b, n));
571  delete hijacked_element_pt->hijack_nodal_value(0,
572  1 + s_fixed_index);
573  }
574  }
575  // boundary 1
576  else if (b == 1)
577  {
578  // hijack DOFs in left element
579  if (n > 0)
580  {
581  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
582  mesh_pt()->boundary_element_pt(b, n - 1));
583  delete hijacked_element_pt->hijack_nodal_value(3,
584  1 + s_fixed_index);
585  }
586  // hijack DOFs in right element
587  if (n < n_node - 1)
588  {
589  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
590  mesh_pt()->boundary_element_pt(b, n));
591  delete hijacked_element_pt->hijack_nodal_value(1,
592  1 + s_fixed_index);
593  }
594  }
595 
596  // boundary 2
597  else if (b == 2)
598  {
599  // hijack DOFs in left element
600  if (n > 0)
601  {
602  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
603  mesh_pt()->boundary_element_pt(b, n - 1));
604  delete hijacked_element_pt->hijack_nodal_value(3,
605  1 + s_fixed_index);
606  }
607  if (n < n_node - 1)
608  {
609  // hijack DOFs in right element
610  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
611  mesh_pt()->boundary_element_pt(b, n));
612  delete hijacked_element_pt->hijack_nodal_value(2,
613  1 + s_fixed_index);
614  }
615  }
616  // boundary 3
617  else if (b == 3)
618  {
619  // hijack DOFs in left element
620  if (n > 0)
621  {
622  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
623  mesh_pt()->boundary_element_pt(b, n - 1));
624  delete hijacked_element_pt->hijack_nodal_value(2,
625  1 + s_fixed_index);
626  }
627  if (n < n_node - 1)
628  {
629  // hijack DOFs in right element
630  hijacked_element_pt = dynamic_cast<Hijacked<BiharmonicElement<2>>*>(
631  mesh_pt()->boundary_element_pt(b, n));
632  delete hijacked_element_pt->hijack_nodal_value(0,
633  1 + s_fixed_index);
634  }
635  }
636 
637  // create the boundary point element
638  BiharmonicFluidBoundaryElement* boundary_point_element_pt =
639  new BiharmonicFluidBoundaryElement(mesh_pt()->boundary_node_pt(b, n),
640  s_fixed_index);
641 
642  // add element to mesh
643  mesh_pt()->add_element_pt(boundary_point_element_pt);
644 
645  // increment number of non bulk elements
646  Npoint_element++;
647  }
648  }
649  }
650 
651 
652  //=============================================================================
653  /// Impose a prescribed fluid flow comprising the velocity normal to
654  /// the boundary (u_imposed_fn[0]) and the velocity tangential to the
655  /// boundary (u_imposed_fn[1])
656  //=============================================================================
657  template<unsigned DIM>
659  const unsigned& b, FluidBCFctPt u_imposed_fn)
660  {
661  // number of nodes on boundary b
662  unsigned n_node = mesh_pt()->nboundary_node(b);
663 
664  // fixed faced index for boundary
665  int face_index = mesh_pt()->face_index_at_boundary(b, 0);
666 
667  // Need to get the s_fixed_index
668  unsigned s_fixed_index = 0;
669  // Also set the edge sign
670  int edge_sign = 0;
671 
672  switch (face_index)
673  {
674  case -1:
675  s_fixed_index = 0;
676  edge_sign = -1;
677  break;
678 
679  case 1:
680  s_fixed_index = 0;
681  edge_sign = 1;
682  break;
683 
684  case -2:
685  s_fixed_index = 1;
686  edge_sign = 1;
687  break;
688 
689  case 2:
690  s_fixed_index = 1;
691  edge_sign = -1;
692  break;
693 
694  default:
695  throw OomphLibError("Face Index not +/-1 or +/-2: Need 2D QElements",
696  OOMPH_CURRENT_FUNCTION,
697  OOMPH_EXCEPTION_LOCATION);
698  }
699 
700 
701  // node position along edge b in macro element boundary representation
702  // [-1,1]
703  Vector<double> s(1);
704 
705  // finite difference step
706  const double h = 10e-8;
707 
708  // loop over nodes on boundary b
709  for (unsigned n = 0; n < n_node; n++)
710  {
711  // get m_t and dm_t/ds_t for this node
712  Vector<double> m_N(2);
713  mesh_pt()->boundary_node_pt(b, n)->get_coordinates_on_boundary(b, m_N);
714 
715  // vectors for dx_i/ds_n and dx_i/ds_t
716  Vector<double> dxds_n(2);
717  Vector<double> dxds_t(2);
718  dxds_n[0] =
719  mesh_pt()->boundary_node_pt(b, n)->x_gen(1 + s_fixed_index, 0);
720  dxds_n[1] =
721  mesh_pt()->boundary_node_pt(b, n)->x_gen(1 + s_fixed_index, 1);
722  dxds_t[0] =
723  mesh_pt()->boundary_node_pt(b, n)->x_gen(2 - s_fixed_index, 0);
724  dxds_t[1] =
725  mesh_pt()->boundary_node_pt(b, n)->x_gen(2 - s_fixed_index, 1);
726 
727  // vector for d2xi/ds_n ds_t
728  Vector<double> d2xds_nds_t(2);
729  d2xds_nds_t[0] = mesh_pt()->boundary_node_pt(b, n)->x_gen(3, 0);
730  d2xds_nds_t[1] = mesh_pt()->boundary_node_pt(b, n)->x_gen(3, 1);
731 
732  // compute dn/ds_n
733  double dnds_n =
734  ((dxds_n[0] * dxds_t[1]) - (dxds_n[1] * dxds_t[0])) /
735  (sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]) * edge_sign);
736 
737  // compute dt/ds_n
738  double dtds_n = ((dxds_n[0] * dxds_t[0]) + (dxds_n[1] * dxds_t[1])) /
739  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
740 
741  // compute dt/ds_t
742  double dtds_t = sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
743 
744  // compute d2t/ds_nds_t
745  double d2tds_nds_t =
746  (dxds_t[0] * d2xds_nds_t[0] + dxds_t[1] * d2xds_nds_t[1]) /
747  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
748 
749  // get imposed velocities
750  Vector<double> u(2);
751  (*u_imposed_fn)(m_N[0], u);
752  u[0] *= edge_sign;
753  u[1] *= -edge_sign;
754 
755  // compute d/dm_t(dpsidn) by finite difference
756  double ddm_tdudn = 0;
757  double ddm_tdudt = 0;
758  Vector<double> u_L(2);
759  Vector<double> u_R(2);
760  // evaluate dudn_fn about current node
761  if (n == 0)
762  {
763  (*u_imposed_fn)(m_N[0], u_L);
764  (*u_imposed_fn)(m_N[0] + h, u_R);
765  }
766  else if (n == n_node - 1)
767  {
768  (*u_imposed_fn)(m_N[0] - h, u_L);
769  (*u_imposed_fn)(m_N[0], u_R);
770  }
771  else
772  {
773  (*u_imposed_fn)(m_N[0] - 0.5 * h, u_L);
774  (*u_imposed_fn)(m_N[0] + 0.5 * h, u_R);
775  }
776  // compute
777  ddm_tdudn = (u_R[1] - u_L[1]) / h;
778  ddm_tdudt = (u_R[0] - u_L[0]) / h;
779 
780  // compute du/ds_t
781  double duds_t = dtds_t * u[0];
782 
783  // compute du/ds_n
784  double duds_n = dnds_n * u[1] + dtds_n * u[0];
785 
786  // compute d2u/ds_n ds_t
787  double d2uds_nds_t = dnds_n * m_N[1] * ddm_tdudn + d2tds_nds_t * u[0] +
788  dtds_n * m_N[1] * ddm_tdudt;
789 
790  // pin du/ds_n dof and set value
791  mesh_pt()->boundary_node_pt(b, n)->pin(1 + s_fixed_index);
792  mesh_pt()->boundary_node_pt(b, n)->set_value(1 + s_fixed_index, duds_n);
793 
794  // pin du/ds_t dof and set value
795  mesh_pt()->boundary_node_pt(b, n)->pin(2 - s_fixed_index);
796  mesh_pt()->boundary_node_pt(b, n)->set_value(2 - s_fixed_index, duds_t);
797 
798  // pin du/ds_t dof and set value
799  mesh_pt()->boundary_node_pt(b, n)->pin(3);
800  mesh_pt()->boundary_node_pt(b, n)->set_value(3, d2uds_nds_t);
801  }
802  }
803 
804 
805  //=============================================================================
806  /// documents the solution, and if an exact solution is provided, then
807  /// the error between the numerical and exact solution is presented
808  //=============================================================================
809  template<unsigned DIM>
811  DocInfo& doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
812  {
813  // create an output stream
814  std::ofstream some_file;
815  std::ostringstream filename;
816 
817  // Number of plot points: npts x npts
818  unsigned npts = 5;
819 
820  // Output solution
821  filename << doc_info.directory() << "/soln_" << doc_info.label() << ".dat";
822  some_file.open(filename.str().c_str());
823  mesh_pt()->output(some_file, npts);
824  some_file.close();
825 
826  // Output fluid velocity solution
827  filename.str("");
828  filename << doc_info.directory() << "/soln_velocity_" << doc_info.label()
829  << ".dat";
830  some_file.open(filename.str().c_str());
831  unsigned n_element = mesh_pt()->nelement();
832  for (unsigned i = 0; i < n_element - Npoint_element; i++)
833  {
834  BiharmonicElement<2>* biharmonic_element_pt =
835  dynamic_cast<BiharmonicElement<2>*>(mesh_pt()->element_pt(i));
836  biharmonic_element_pt->output_fluid_velocity(some_file, npts);
837  }
838  some_file.close();
839 
840  // if the exact solution is not provided
841  if (exact_soln_pt != 0)
842  {
843  // Output exact solution
844  filename.str("");
845  filename << doc_info.directory() << "/exact_soln_" << doc_info.label()
846  << ".dat";
847  some_file.open(filename.str().c_str());
848  mesh_pt()->output_fct(some_file, npts, exact_soln_pt);
849  some_file.close();
850 
851  // Doc error and return of the square of the L2 error
852  double error, norm;
853  filename.str("");
854  filename << doc_info.directory() << "/error_" << doc_info.label()
855  << ".dat";
856  some_file.open(filename.str().c_str());
857  mesh_pt()->compute_error(some_file, exact_soln_pt, error, norm);
858  some_file.close();
859 
860  // Doc L2 error and norm of solution
861  oomph_info << "\nNorm of error : " << sqrt(error) << std::endl;
862  oomph_info << "Norm of solution: " << sqrt(norm) << std::endl
863  << std::endl;
864  }
865  }
866 
867  // ensure build
868  template class BiharmonicFluidProblem<2>;
869  template class BiharmonicProblem<2>;
870 
871 } // namespace oomph
e
Definition: cfortran.h:571
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
biharmonic element class
void output_fluid_velocity(std::ostream &outfile, const unsigned &nplot)
output fluid velocity field
Point equation element used to impose the traction free edge (i.e. du/dn = 0) on the boundary when dt...
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...
FluxFctPt & flux0_fct_pt()
Access function for the flux0 function pointer.
FluxFctPt & flux1_fct_pt()
Access function for the flux1 function pointer.
Biharmonic Plate Problem Class - for problems where the load can be assumed to be acting normal to th...
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 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 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'.
Information for documentation of results: Directory and file number to enable output in the form RESL...
std::string & label()
String used (e.g.) for labeling output files.
std::string directory() const
Output directory.
unsigned & number()
Number used (e.g.) for labeling output files.
A general Finite Element class.
Definition: elements.h:1317
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition: elements.h:1436
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
int local_eqn_number(const unsigned long &ieqn_global) const
Return the local equation number corresponding to the ieqn_global-th global equation number....
Definition: elements.h:730
Hijacked elements are elements in which one or more Data values that affect the element's residuals,...
Data * hijack_nodal_value(const unsigned &n, const unsigned &i, const bool &return_data=true)
Hijack the i-th value stored at node n. Optionally return a custom-made (copied) data object that con...
A general mesh class.
Definition: mesh.h:67
double & x_gen(const unsigned &k, const unsigned &i)
Reference to the generalised position x(k,i). ‘Type’: k; Coordinate direction: i.
Definition: nodes.h:1126
double value(const unsigned &i) const
Return i-th value (dofs or pinned) at this node either directly or via hanging node representation....
Definition: nodes.cc:2408
An OomphLibError object which should be thrown when an run-time error is encountered....
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...