axisym_fvk_elements.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 // Header file for axisymmetric FoepplvonKarman elements
27 #ifndef OOMPH_AXISYM_FOEPPLVONKARMAN_ELEMENTS_HEADER
28 #define OOMPH_AXISYM_FOEPPLVONKARMAN_ELEMENTS_HEADER
29 
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 #include "../generic/nodes.h"
37 #include "../generic/Qelements.h"
38 #include "../generic/oomph_utilities.h"
39 #include "../generic/element_with_external_element.h"
40 
41 namespace oomph
42 {
43  //=============================================================
49  //=============================================================
50  class AxisymFoepplvonKarmanEquations : public virtual FiniteElement
51  {
52  public:
55  typedef void (*AxisymFoepplvonKarmanPressureFctPt)(const double& r,
56  double& f);
57 
63  {
64  // Set all the physical constants to the default value (zero)
66  Linear_bending_model = false;
67  }
68 
71  const AxisymFoepplvonKarmanEquations& dummy) = delete;
72 
75 
77  const double& eta() const
78  {
79  return *Eta_pt;
80  }
81 
83  double*& eta_pt()
84  {
85  return Eta_pt;
86  }
87 
99  virtual inline unsigned nodal_index_fvk(const unsigned& i = 0) const
100  {
101  return i;
102  }
103 
105  void output(std::ostream& outfile)
106  {
107  const unsigned n_plot = 5;
108  output(outfile, n_plot);
109  }
110 
113  void output(std::ostream& outfile, const unsigned& n_plot);
114 
116  void output(FILE* file_pt)
117  {
118  const unsigned n_plot = 5;
119  output(file_pt, n_plot);
120  }
121 
124  void output(FILE* file_pt, const unsigned& n_plot);
125 
127  void output_fct(std::ostream& outfile,
128  const unsigned& n_plot,
130 
134  virtual void output_fct(
135  std::ostream& outfile,
136  const unsigned& n_plot,
137  const double& time,
139  {
140  throw OomphLibError(
141  "There is no time-dependent output_fct() for Foeppl von Karman"
142  "elements ",
145  }
146 
148  void compute_error(std::ostream& outfile,
150  double& error,
151  double& norm);
152 
153 
155  void compute_error(std::ostream& outfile,
157  const double& time,
158  double& error,
159  double& norm)
160  {
161  throw OomphLibError(
162  "There is no time-dependent compute_error() for Foeppl von Karman"
163  "elements",
166  }
167 
170  {
171  return Pressure_fct_pt;
172  }
173 
176  {
177  return Pressure_fct_pt;
178  }
179 
182  {
183  return Airy_forcing_fct_pt;
184  }
185 
188  {
189  return Airy_forcing_fct_pt;
190  }
191 
196  inline virtual void get_pressure_fvk(const unsigned& ipt,
197  const double& r,
198  double& pressure) const
199  {
200  // If no pressure function has been set, return zero
201  if (Pressure_fct_pt == 0)
202  {
203  pressure = 0.0;
204  }
205  else
206  {
207  // Get pressure strength
208  (*Pressure_fct_pt)(r, pressure);
209  }
210  }
211 
216  inline virtual void get_airy_forcing_fvk(const unsigned& ipt,
217  const double& r,
218  double& airy_forcing) const
219  {
220  // If no Airy forcing function has been set, return zero
221  if (Airy_forcing_fct_pt == 0)
222  {
223  airy_forcing = 0.0;
224  }
225  else
226  {
227  // Get Airy forcing strength
228  (*Airy_forcing_fct_pt)(r, airy_forcing);
229  }
230  }
231 
234  Vector<double>& gradient) const
235  {
236  // Find out how many nodes there are in the element
237  const unsigned n_node = nnode();
238 
239  // Get the index at which the unknown is stored
240  const unsigned w_nodal_index = nodal_index_fvk(0);
241 
242  // Set up memory for the shape and test functions
243  Shape psi(n_node);
244  DShape dpsidr(n_node, 1);
245 
246  // Call the derivatives of the shape and test functions
247  dshape_eulerian(s, psi, dpsidr);
248 
249  // Initialise to zero
250  gradient[0] = 0.0;
251 
252  // Loop over nodes
253  for (unsigned l = 0; l < n_node; l++)
254  {
255  gradient[0] += this->nodal_value(l, w_nodal_index) * dpsidr(l, 0);
256  }
257  }
258 
261 
262 
263  // hierher Jacobian not yet implemented
264  // void fill_in_contribution_to_jacobian(Vector<double> &residuals,
265  // DenseMatrix<double> &jacobian);
266 
269  inline double interpolated_w_fvk(const Vector<double>& s) const
270  {
271  // Find number of nodes
272  const unsigned n_node = nnode();
273 
274  // Get the index at which the unknown is stored
275  const unsigned w_nodal_index = nodal_index_fvk(0);
276 
277  // Local shape function
278  Shape psi(n_node);
279 
280  // Find values of shape function
281  shape(s, psi);
282 
283  // Initialise value of u
284  double interpolated_w = 0.0;
285 
286  // Loop over the local nodes and sum
287  for (unsigned l = 0; l < n_node; l++)
288  {
289  interpolated_w += this->nodal_value(l, w_nodal_index) * psi[l];
290  }
291 
292  return (interpolated_w);
293  }
294 
298  double& sigma_r_r,
299  double& sigma_phi_phi);
300 
301 
303  unsigned self_test();
304 
309  {
310  // Set the boolean flag
311  Linear_bending_model = true;
312 
313  // Get the index of the first FvK nodal value
314  unsigned first_fvk_nodal_index = nodal_index_fvk();
315 
316  // Get the total number of FvK nodal values (assuming they are stored
317  // contiguously) at node 0 (it's the same at all nodes anyway)
318  unsigned total_fvk_nodal_indices = 6;
319 
320  // Get the number of nodes in this element
321  unsigned n_node = nnode();
322 
323  // Loop over the appropriate nodal indices
324  for (unsigned index = first_fvk_nodal_index + 2;
325  index < first_fvk_nodal_index + total_fvk_nodal_indices;
326  index++)
327  {
328  // Loop over the nodes in the element
329  for (unsigned inod = 0; inod < n_node; inod++)
330  {
331  // Pin the nodal value at the current index
332  node_pt(inod)->pin(index);
333  }
334  }
335  }
336 
337 
338  protected:
342  const Vector<double>& s,
343  Shape& psi,
344  DShape& dpsidr,
345  Shape& test,
346  DShape& dtestdr) const = 0;
347 
348 
352  const unsigned& ipt,
353  Shape& psi,
354  DShape& dpsidr,
355  Shape& test,
356  DShape& dtestdr) const = 0;
357 
359  double* Eta_pt;
360 
363 
366 
368  static double Default_Physical_Constant_Value;
369 
373  };
374 
375 
379 
380 
381  //======================================================================
385  //======================================================================
386  template<unsigned NNODE_1D>
388  : public virtual QElement<1, NNODE_1D>,
389  public virtual AxisymFoepplvonKarmanEquations
390  {
391  private:
394  static const unsigned Initial_Nvalue;
395 
396  public:
400  : QElement<1, NNODE_1D>(), AxisymFoepplvonKarmanEquations()
401  {
402  }
403 
406  const AxisymFoepplvonKarmanElement<NNODE_1D>& dummy) = delete;
407 
410 
413  inline unsigned required_nvalue(const unsigned& n) const
414  {
415  return Initial_Nvalue;
416  }
417 
418 
421  void output(std::ostream& outfile)
422  {
424  }
425 
428  void output(std::ostream& outfile, const unsigned& n_plot)
429  {
431  }
432 
435  void output(FILE* file_pt)
436  {
438  }
439 
442  void output(FILE* file_pt, const unsigned& n_plot)
443  {
445  }
446 
449  void output_fct(std::ostream& outfile,
450  const unsigned& n_plot,
452  {
454  outfile, n_plot, exact_soln_pt);
455  }
456 
460  void output_fct(std::ostream& outfile,
461  const unsigned& n_plot,
462  const double& time,
464  {
466  outfile, n_plot, time, exact_soln_pt);
467  }
468 
469 
470  protected:
474  Shape& psi,
475  DShape& dpsidr,
476  Shape& test,
477  DShape& dtestdr) const;
478 
482  const unsigned& ipt,
483  Shape& psi,
484  DShape& dpsidr,
485  Shape& test,
486  DShape& dtestdr) const;
487  };
488 
489 
490  // Inline functions:
491 
492  //======================================================================
497  //======================================================================
498  template<unsigned NNODE_1D>
500  NNODE_1D>::dshape_and_dtest_eulerian_axisym_fvk(const Vector<double>& s,
501  Shape& psi,
502  DShape& dpsidr,
503  Shape& test,
504  DShape& dtestdr) const
505 
506  {
507  // Call the geometrical shape functions and derivatives
508  const double J = this->dshape_eulerian(s, psi, dpsidr);
509 
510  // Set the test functions equal to the shape functions
511  test = psi;
512  dtestdr = dpsidr;
513 
514  // Return the jacobian
515  return J;
516  }
517 
518 
519  //======================================================================
524  //======================================================================
525  template<unsigned NNODE_1D>
528  Shape& psi,
529  DShape& dpsidr,
530  Shape& test,
531  DShape& dtestdr) const
532  {
533  // Call the geometrical shape functions and derivatives
534  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidr);
535 
536  // Set the pointers of the test functions
537  test = psi;
538  dtestdr = dpsidr;
539 
540  // Return the jacobian
541  return J;
542  }
543 
544 
548 
549 
550  //======================================================================
555  //======================================================================
556  template<unsigned NNODE_1D, class FLUID_ELEMENT>
558  : public virtual AxisymFoepplvonKarmanElement<NNODE_1D>,
559  public virtual ElementWithExternalElement
560  {
561  public:
565  {
566  // Set source element storage: one interaction with an external
567  // element -- the fluid bulk element that provides the pressure
568  this->set_ninteraction(1);
569  }
570 
573 
576  const double& q() const
577  {
578  return *Q_pt;
579  }
580 
583  double*& q_pt()
584  {
585  return Q_pt;
586  }
587 
590  virtual unsigned ngeom_data() const
591  {
592  return this->nnode();
593  }
594 
597  virtual Data* geom_data_pt(const unsigned& j)
598  {
599  return this->node_pt(j);
600  }
601 
606  {
607  const unsigned t = 0;
608  this->position(t, zeta, r);
609  }
610 
614  void position(const unsigned& t,
615  const Vector<double>& zeta,
616  Vector<double>& r) const
617  {
618  // Find number of nodes
619  const unsigned n_node = this->nnode();
620 
621  // Get the index at which the poisson unknown is stored
622  const unsigned w_nodal_index = this->nodal_index_fvk(0);
623 
624  // Local shape function
625  Shape psi(n_node);
626 
627  // Find values of shape function
628  this->shape(zeta, psi);
629 
630  // Initialise
631  double interpolated_w = 0.0;
632  double interpolated_r = 0.0;
633 
634  // Loop over the local nodes and sum
635  for (unsigned l = 0; l < n_node; l++)
636  {
637  interpolated_w += this->nodal_value(t, l, w_nodal_index) * psi[l];
638  interpolated_r += this->node_pt(l)->x(t, 0) * psi[l];
639  }
640 
641  // Assign
642  r[0] = interpolated_r;
643  r[1] = interpolated_w;
644  }
645 
649  const unsigned& j,
650  Vector<double>& drdt)
651  {
652  // Find number of nodes
653  const unsigned n_node = this->nnode();
654 
655  // Get the index at which the poisson unknown is stored
656  const unsigned w_nodal_index = this->nodal_index_fvk(0);
657 
658  // Local shape function
659  Shape psi(n_node);
660 
661  // Find values of shape function
662  this->shape(zeta, psi);
663 
664  // Initialise
665  double interpolated_dwdt = 0.0;
666  double interpolated_drdt = 0.0;
667 
668  // Loop over the local nodes and sum
669  for (unsigned l = 0; l < n_node; l++)
670  {
671  // Get the timestepper
673 
674  // If we are doing an unsteady solve then calculate the derivative
675  if (!time_stepper_pt->is_steady())
676  {
677  // Get the number of values required to represent history
678  const unsigned n_time = time_stepper_pt->ntstorage();
679 
680  // Loop over history values
681  for (unsigned t = 0; t < n_time; t++)
682  {
683  // Add the contribution to the derivative
684  interpolated_dwdt += time_stepper_pt->weight(1, t) *
685  this->nodal_value(t, l, w_nodal_index) *
686  psi[l];
687  }
688  }
689  }
690 
691  // Assign
692  drdt[0] = interpolated_drdt;
693  drdt[1] = interpolated_dwdt;
694  }
695 
696 
701  inline virtual void get_pressure_fvk(const unsigned& ipt,
702  const double& r,
703  double& pressure) const
704  {
705  pressure = 0.0;
706 
707  // Get underlying version
709 
710  // Get FSI parameter
711  const double q_value = q();
712 
713  // Get fluid element
714  FLUID_ELEMENT* ext_el_pt =
715  dynamic_cast<FLUID_ELEMENT*>(external_element_pt(0, ipt));
717 
718  // Outer unit normal is vertically upward (in z direction)
719  // (within an essentiall flat) model for the wall)
721  normal[0] = 0.0;
722  normal[1] = 1.0;
723 
724  // Get traction
725  Vector<double> traction(3);
726  ext_el_pt->traction(s_ext, normal, traction);
727 
728  // Add z-component of traction
729  pressure -= q_value * traction[1];
730  }
731 
732 
734  void output_integration_points(std::ostream& outfile)
735  {
736  // How many nodes do we have?
737  unsigned nnod = this->nnode();
738  Shape psi(nnod);
739 
740  // Get the index at which the unknown is stored
741  const unsigned w_nodal_index = this->nodal_index_fvk(0);
742 
743  // Loop over the integration points
744  const unsigned n_intpt = this->integral_pt()->nweight();
745  outfile << "ZONE I=" << n_intpt << std::endl;
746  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
747  {
748  // Get shape fct
749  Vector<double> s(1);
750  s[0] = this->integral_pt()->knot(ipt, 0);
751  shape(s, psi);
752 
753  // Initialise
754  double interpolated_w = 0.0;
755  double interpolated_r = 0.0;
756 
757  // Loop over the local nodes and sum
758  for (unsigned l = 0; l < nnod; l++)
759  {
760  interpolated_w += this->nodal_value(l, w_nodal_index) * psi[l];
761  interpolated_r += this->node_pt(l)->x(0) * psi[l];
762  }
763 
764  // Get fluid element
765  FLUID_ELEMENT* ext_el_pt =
766  dynamic_cast<FLUID_ELEMENT*>(external_element_pt(0, ipt));
768 
769  // Get veloc
770  Vector<double> veloc(3);
771  ext_el_pt->interpolated_u_axi_nst(s_ext, veloc);
772  Vector<double> x(2);
773  ext_el_pt->interpolated_x(s_ext, x);
774 
775  outfile << interpolated_r << " " << interpolated_w << " " << veloc[0]
776  << " " << veloc[1] << " " << x[0] << " " << x[1] << " "
777  << std::endl;
778  }
779  }
780 
781 
783  void output_adjacent_fluid_elements(std::ostream& outfile,
784  const unsigned& nplot)
785  {
786  // Loop over the integration points
787  const unsigned n_intpt = this->integral_pt()->nweight();
788  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
789  {
790  // Get fluid element
791  FLUID_ELEMENT* ext_el_pt =
792  dynamic_cast<FLUID_ELEMENT*>(external_element_pt(0, ipt));
793 
794  // Dump it
795  ext_el_pt->output(outfile, nplot);
796  }
797  }
798 
802  {
804  }
805 
809  {
811  }
812 
813 
817  {
819  }
820 
821 
825  {
827  }
828 
829 
833  {
834  // Don't update elements repeatedly
835  std::map<FLUID_ELEMENT*, bool> done;
836 
837  // Number of integration points
838  unsigned n_intpt = integral_pt()->nweight();
839 
840  // Loop over all integration points in wall element
841  for (unsigned iint = 0; iint < n_intpt; iint++)
842  {
843  // Get fluid element that affects this integration point
844  FLUID_ELEMENT* el_f_pt =
845  dynamic_cast<FLUID_ELEMENT*>(external_element_pt(0, iint));
846 
847  // Is there an adjacent fluid element?
848  if (el_f_pt != 0)
849  {
850  // Have we updated its positions yet?
851  if (!done[el_f_pt])
852  {
853  // Update nodal positions
854  el_f_pt->node_update();
855  done[el_f_pt] = true;
856  }
857  }
858  }
859  }
860 
861 
864  void output(std::ostream& outfile, const unsigned& n_plot)
865  {
866  // Vector of local coordinates
867  Vector<double> s(1);
868 
869  // Tecplot header info
870  outfile << "ZONE\n";
871 
872  // Loop over plot points
873  unsigned num_plot_points = nplot_points(n_plot);
874  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
875  {
876  // Get local coordinates of plot point
877  get_s_plot(iplot, n_plot, s);
878 
879  // Get velocity
880  Vector<double> drdt(2);
881  dposition_dt(s, 1, drdt);
882 
883  // Get stress
884  double sigma_r_r = 0.0;
885  double sigma_phi_phi = 0.0;
886  bool success = this->interpolated_stress(s, sigma_r_r, sigma_phi_phi);
887  if (success)
888  {
889  outfile << this->interpolated_x(s, 0) << " "
890  << this->interpolated_w_fvk(s) << " " << drdt[0] << " "
891  << drdt[1] << " " << sigma_r_r << " " << sigma_phi_phi
892  << std::endl;
893  }
894  }
895  }
896 
897 
898  protected:
902  double* Q_pt;
903  };
904 
905 
906 } // namespace oomph
907 
908 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: axisym_displ_based_fvk_elements.h:404
void output(FILE *file_pt, const unsigned &n_plot)
Definition: axisym_fvk_elements.h:442
double dshape_and_dtest_eulerian_axisym_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const
void output(FILE *file_pt)
Definition: axisym_fvk_elements.h:435
double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: axisym_fvk_elements.h:460
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: axisym_fvk_elements.h:428
AxisymFoepplvonKarmanElement(const AxisymFoepplvonKarmanElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Definition: axisym_fvk_elements.h:421
static const unsigned Initial_Nvalue
Set the data for the number of Variables at each node - 3.
Definition: axisym_displ_based_fvk_elements.h:408
unsigned required_nvalue(const unsigned &n) const
Definition: axisym_fvk_elements.h:413
void operator=(const AxisymFoepplvonKarmanElement< NNODE_1D > &)=delete
Broken assignment operator.
AxisymFoepplvonKarmanElement()
Definition: axisym_fvk_elements.h:399
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: axisym_fvk_elements.h:449
Definition: axisym_displ_based_fvk_elements.h:51
AxisymFoepplvonKarmanPressureFctPt pressure_fct_pt() const
Access function: Pointer to pressure function. Const version.
Definition: axisym_fvk_elements.h:175
AxisymFoepplvonKarmanPressureFctPt Pressure_fct_pt
Pointer to pressure function:
Definition: axisym_displ_based_fvk_elements.h:375
static double Default_Physical_Constant_Value
Default value for physical constants.
Definition: axisym_displ_based_fvk_elements.h:386
void operator=(const AxisymFoepplvonKarmanEquations &)=delete
Broken assignment operator.
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: axisym_fvk_elements.h:134
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: axisym_fvk_elements.h:105
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals with this element's contribution.
bool interpolated_stress(const Vector< double > &s, double &sigma_r_r, double &sigma_phi_phi) const
Definition: axisym_displ_based_fvk_elements.cc:234
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: axisym_fvk_elements.h:116
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: axisym_fvk_elements.h:155
unsigned self_test()
Self-test: Return 0 for OK.
AxisymFoepplvonKarmanPressureFctPt Airy_forcing_fct_pt
Pointer to Airy forcing function.
Definition: axisym_fvk_elements.h:365
double * Eta_pt
Pointer to FvK parameter.
Definition: axisym_displ_based_fvk_elements.h:372
virtual double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
virtual unsigned nodal_index_fvk(const unsigned &i=0) const
Definition: axisym_fvk_elements.h:99
void(* AxisymFoepplvonKarmanPressureFctPt)(const double &r, double &f)
Definition: axisym_displ_based_fvk_elements.h:55
void output(FILE *file_pt, const unsigned &n_plot)
void get_gradient_of_deflection(const Vector< double > &s, Vector< double > &gradient) const
Get gradient of deflection: gradient[i] = dw/dr_i *‍/.
Definition: axisym_fvk_elements.h:233
AxisymFoepplvonKarmanPressureFctPt & airy_forcing_fct_pt()
Access function: Pointer to Airy forcing function.
Definition: axisym_fvk_elements.h:181
AxisymFoepplvonKarmanPressureFctPt & pressure_fct_pt()
Access function: Pointer to pressure function.
Definition: axisym_fvk_elements.h:169
void output(std::ostream &outfile, const unsigned &n_plot)
AxisymFoepplvonKarmanEquations()
Definition: axisym_fvk_elements.h:61
virtual double dshape_and_dtest_eulerian_axisym_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
bool Linear_bending_model
Definition: axisym_displ_based_fvk_elements.h:382
const double & eta() const
FvK parameter.
Definition: axisym_fvk_elements.h:77
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,w_exact at n_plot plot points.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
void use_linear_bending_model()
Definition: axisym_fvk_elements.h:308
double interpolated_w_fvk(const Vector< double > &s) const
Definition: axisym_fvk_elements.h:269
virtual void get_pressure_fvk(const unsigned &ipt, const double &r, double &pressure) const
Definition: axisym_fvk_elements.h:196
double *& eta_pt()
Pointer to FvK parameter.
Definition: axisym_fvk_elements.h:83
AxisymFoepplvonKarmanEquations(const AxisymFoepplvonKarmanEquations &dummy)=delete
Broken copy constructor.
virtual void get_airy_forcing_fvk(const unsigned &ipt, const double &r, double &airy_forcing) const
Definition: axisym_fvk_elements.h:216
AxisymFoepplvonKarmanPressureFctPt airy_forcing_fct_pt() const
Access function: Pointer to Airy forcing function. Const version.
Definition: axisym_fvk_elements.h:187
Definition: shape.h:278
Definition: nodes.h:86
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Definition: element_with_external_element.h:56
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Definition: element_with_external_element.h:136
void set_ninteraction(const unsigned &n_interaction)
Definition: element_with_external_element.h:178
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Definition: element_with_external_element.h:107
Definition: axisym_fvk_elements.h:560
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: axisym_fvk_elements.h:864
void reset_after_external_interaction_geometric_fd()
Definition: axisym_fvk_elements.h:808
void node_update_adjacent_fluid_elements()
Definition: axisym_fvk_elements.h:832
virtual ~FSIAxisymFoepplvonKarmanElement()
Empty virtual destructor.
Definition: axisym_fvk_elements.h:572
double *& q_pt()
Definition: axisym_fvk_elements.h:583
virtual Data * geom_data_pt(const unsigned &j)
Definition: axisym_fvk_elements.h:597
void update_before_external_interaction_geometric_fd()
Definition: axisym_fvk_elements.h:801
double * Q_pt
Definition: axisym_fvk_elements.h:902
void output_adjacent_fluid_elements(std::ostream &outfile, const unsigned &nplot)
Output adjacent fluid elements (for checking of fsi setup)
Definition: axisym_fvk_elements.h:783
void dposition_dt(const Vector< double > &zeta, const unsigned &j, Vector< double > &drdt)
Definition: axisym_fvk_elements.h:648
void output_integration_points(std::ostream &outfile)
Output integration points (for checking of fsi setup)
Definition: axisym_fvk_elements.h:734
virtual void get_pressure_fvk(const unsigned &ipt, const double &r, double &pressure) const
Definition: axisym_fvk_elements.h:701
void position(const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
Definition: axisym_fvk_elements.h:614
void position(const Vector< double > &zeta, Vector< double > &r) const
Definition: axisym_fvk_elements.h:605
const double & q() const
Definition: axisym_fvk_elements.h:576
void update_in_external_interaction_geometric_fd(const unsigned &i)
Definition: axisym_fvk_elements.h:816
virtual unsigned ngeom_data() const
Definition: axisym_fvk_elements.h:590
FSIAxisymFoepplvonKarmanElement()
Constructor.
Definition: axisym_fvk_elements.h:563
void reset_in_external_interaction_geometric_fd(const unsigned &i)
Definition: axisym_fvk_elements.h:824
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
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 void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
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
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
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
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.
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
Definition: oomph_definitions.h:222
Definition: Qelements.h:459
Definition: shape.h:76
Definition: timesteppers.h:231
unsigned ntstorage() const
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Definition: timesteppers.h:389
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
#define FLUID_ELEMENT
RealScalar s
Definition: level1_cplx_impl.h:130
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152
r
Definition: UniformPSDSelfTest.py:20
void normal(const Vector< double > &x, Vector< double > &normal)
Definition: free_surface_rotation.cc:65
int error
Definition: calibrate.py:297
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#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