unsteady_heat_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 UnsteadyHeat elements
27 #ifndef OOMPH_UNSTEADY_HEAT_ELEMENTS_HEADER
28 #define OOMPH_UNSTEADY_HEAT_ELEMENTS_HEADER
29 
30 // Config header generated by autoconfig
31 #ifdef HAVE_CONFIG_H
32 #include <oomph-lib-config.h>
33 #endif
34 
35 
36 // OOMPH-LIB headers
37 #include "../generic/projection.h"
38 #include "../generic/nodes.h"
39 #include "../generic/Qelements.h"
40 #include "../generic/oomph_utilities.h"
41 
42 
43 namespace oomph
44 {
48  {
49  public:
52  typedef void (*UnsteadyHeatSourceFctPt)(const double& time,
53  const Vector<double>& x,
54  double& u);
55 
58  };
59 
60  //=============================================================
69  //=============================================================
70  template<unsigned DIM>
72  {
73  public:
76  typedef void (*UnsteadyHeatSourceFctPt)(const double& time,
77  const Vector<double>& x,
78  double& u);
79 
80 
86  {
87  // Set Alpha and Beta parameter to default (one for natural scaling of
88  // time)
91  }
92 
93 
96 
98  // Commented out broken assignment operator because this can lead to a
99  // conflict warning when used in the virtual inheritence hierarchy.
100  // Essentially the compiler doesn't realise that two separate
101  // implementations of the broken function are the same and so, quite
102  // rightly, it shouts.
103  /*void operator=(const UnsteadyHeatEquations&) = delete;*/
104 
112  virtual inline unsigned u_index_ust_heat() const
113  {
114  return 0;
115  }
116 
119  double du_dt_ust_heat(const unsigned& n) const
120  {
121  // Get the data's timestepper
123 
124  // Initialise dudt
125  double dudt = 0.0;
126 
127  // Loop over the timesteps, if there is a non Steady timestepper
128  if (!time_stepper_pt->is_steady())
129  {
130  // Find the index at which the variable is stored
131  const unsigned u_nodal_index = u_index_ust_heat();
132 
133  // Number of timsteps (past & present)
134  const unsigned n_time = time_stepper_pt->ntstorage();
135 
136  // Add the contributions to the time derivative
137  for (unsigned t = 0; t < n_time; t++)
138  {
139  dudt +=
140  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
141  }
142  }
143  return dudt;
144  }
145 
148  void disable_ALE()
149  {
150  ALE_is_disabled = true;
151  }
152 
153 
158  void enable_ALE()
159  {
160  ALE_is_disabled = false;
161  }
162 
164  void compute_norm(double& norm);
165 
167  void output(std::ostream& outfile)
168  {
169  unsigned nplot = 5;
170  output(outfile, nplot);
171  }
172 
173 
176  void output(std::ostream& outfile, const unsigned& nplot);
177 
179  void output(FILE* file_pt)
180  {
181  unsigned n_plot = 5;
182  output(file_pt, n_plot);
183  }
184 
185 
188  void output(FILE* file_pt, const unsigned& n_plot);
189 
190 
192  void output_fct(std::ostream& outfile,
193  const unsigned& nplot,
195 
196 
199  virtual void output_fct(
200  std::ostream& outfile,
201  const unsigned& nplot,
202  const double& time,
204 
205 
207  void compute_error(std::ostream& outfile,
209  double& error,
210  double& norm);
211 
212 
214  void compute_error(std::ostream& outfile,
216  const double& time,
217  double& error,
218  double& norm);
219 
220 
223  {
224  return Source_fct_pt;
225  }
226 
227 
230  {
231  return Source_fct_pt;
232  }
233 
234 
237  virtual inline void get_source_ust_heat(const double& t,
238  const unsigned& ipt,
239  const Vector<double>& x,
240  double& source) const
241  {
242  // If no source function has been set, return zero
243  if (Source_fct_pt == 0)
244  {
245  source = 0.0;
246  }
247  else
248  {
249  // Get source strength
250  (*Source_fct_pt)(t, x, source);
251  }
252  }
253 
255  const double& alpha() const
256  {
257  return *Alpha_pt;
258  }
259 
261  double*& alpha_pt()
262  {
263  return Alpha_pt;
264  }
265 
266 
268  const double& beta() const
269  {
270  return *Beta_pt;
271  }
272 
274  double*& beta_pt()
275  {
276  return Beta_pt;
277  }
278 
281  {
282  // Find out how many nodes there are in the element
283  unsigned n_node = nnode();
284 
285  // Find the index at which the variable is stored
286  unsigned u_nodal_index = u_index_ust_heat();
287 
288  // Set up memory for the shape and test functions
289  Shape psi(n_node);
290  DShape dpsidx(n_node, DIM);
291 
292  // Call the derivatives of the shape and test functions
293  dshape_eulerian(s, psi, dpsidx);
294 
295  // Initialise to zero
296  for (unsigned j = 0; j < DIM; j++)
297  {
298  flux[j] = 0.0;
299  }
300 
301  // Loop over nodes
302  for (unsigned l = 0; l < n_node; l++)
303  {
304  // Loop over derivative directions
305  for (unsigned j = 0; j < DIM; j++)
306  {
307  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
308  }
309  }
310  }
311 
312 
315  {
316  // Call the generic residuals function with flag set to 0
317  // using a dummy matrix argument
319  residuals, GeneralisedElement::Dummy_matrix, 0);
320  }
321 
322 
325  DenseMatrix<double>& jacobian)
326  {
327  // Call the generic routine with the flag set to 1
328  fill_in_generic_residual_contribution_ust_heat(residuals, jacobian, 1);
329  }
330 
331 
333  inline double interpolated_u_ust_heat(const Vector<double>& s) const
334  {
335  // Find number of nodes
336  unsigned n_node = nnode();
337 
338  // Find the index at which the variable is stored
339  unsigned u_nodal_index = u_index_ust_heat();
340 
341  // Local shape function
342  Shape psi(n_node);
343 
344  // Find values of shape function
345  shape(s, psi);
346 
347  // Initialise value of u
348  double interpolated_u = 0.0;
349 
350  // Loop over the local nodes and sum
351  for (unsigned l = 0; l < n_node; l++)
352  {
353  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
354  }
355 
356  return (interpolated_u);
357  }
358 
359 
362  inline double interpolated_u_ust_heat(const unsigned& t,
363  const Vector<double>& s) const
364  {
365  // Find number of nodes
366  unsigned n_node = nnode();
367 
368  // Find the index at which the variable is stored
369  unsigned u_nodal_index = u_index_ust_heat();
370 
371  // Local shape function
372  Shape psi(n_node);
373 
374  // Find values of shape function
375  shape(s, psi);
376 
377  // Initialise value of u
378  double interpolated_u = 0.0;
379 
380  // Loop over the local nodes and sum
381  for (unsigned l = 0; l < n_node; l++)
382  {
383  interpolated_u += nodal_value(t, l, u_nodal_index) * psi[l];
384  }
385 
386  return (interpolated_u);
387  }
388 
389 
392  inline double interpolated_du_dt_ust_heat(const Vector<double>& s) const
393  {
394  // Find number of nodes
395  unsigned n_node = nnode();
396 
397  // Local shape function
398  Shape psi(n_node);
399 
400  // Find values of shape function
401  shape(s, psi);
402 
403  // Initialise value of du/dt
404  double interpolated_dudt = 0.0;
405 
406  // Loop over the local nodes and sum
407  for (unsigned l = 0; l < n_node; l++)
408  {
409  interpolated_dudt += du_dt_ust_heat(l) * psi[l];
410  }
411 
412  return (interpolated_dudt);
413  }
414 
415 
417  unsigned self_test();
418 
419 
420  protected:
424  const Vector<double>& s,
425  Shape& psi,
426  DShape& dpsidx,
427  Shape& test,
428  DShape& dtestdx) const = 0;
429 
430 
434  const unsigned& ipt,
435  Shape& psi,
436  DShape& dpsidx,
437  Shape& test,
438  DShape& dtestdx) const = 0;
439 
443  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
444 
447 
452 
454  double* Alpha_pt;
455 
457  double* Beta_pt;
458 
459  private:
462  static double Default_alpha_parameter;
463 
466  static double Default_beta_parameter;
467  };
468 
469 
473 
474 
475  //======================================================================
478  //======================================================================
479  template<unsigned DIM, unsigned NNODE_1D>
480  class QUnsteadyHeatElement : public virtual QElement<DIM, NNODE_1D>,
481  public virtual UnsteadyHeatEquations<DIM>
482  {
483  private:
486  static const unsigned Initial_Nvalue;
487 
488  public:
492  : QElement<DIM, NNODE_1D>(), UnsteadyHeatEquations<DIM>()
493  {
494  }
495 
498  delete;
499 
501  /*void operator=(const QUnsteadyHeatElement<DIM,NNODE_1D>&) = delete;*/
502 
505  inline unsigned required_nvalue(const unsigned& n) const
506  {
507  return Initial_Nvalue;
508  }
509 
512  void output(std::ostream& outfile)
513  {
515  }
516 
517 
520  void output(std::ostream& outfile, const unsigned& n_plot)
521  {
522  UnsteadyHeatEquations<DIM>::output(outfile, n_plot);
523  }
524 
525 
528  void output(FILE* file_pt)
529  {
531  }
532 
533 
536  void output(FILE* file_pt, const unsigned& n_plot)
537  {
538  UnsteadyHeatEquations<DIM>::output(file_pt, n_plot);
539  }
540 
541 
544  void output_fct(std::ostream& outfile,
545  const unsigned& n_plot,
547  {
548  UnsteadyHeatEquations<DIM>::output_fct(outfile, n_plot, exact_soln_pt);
549  }
550 
551 
555  void output_fct(std::ostream& outfile,
556  const unsigned& n_plot,
557  const double& time,
559  {
561  outfile, n_plot, time, exact_soln_pt);
562  }
563 
564 
565  protected:
569  Shape& psi,
570  DShape& dpsidx,
571  Shape& test,
572  DShape& dtestdx) const;
573 
574 
578  const unsigned& ipt,
579  Shape& psi,
580  DShape& dpsidx,
581  Shape& test,
582  DShape& dtestdx) const;
583  };
584 
585 
586  // Inline functions:
587 
588 
589  //======================================================================
594  //======================================================================
595  template<unsigned DIM, unsigned NNODE_1D>
598  Shape& psi,
599  DShape& dpsidx,
600  Shape& test,
601  DShape& dtestdx) const
602  {
603  // Call the geometrical shape functions and derivatives
604  double J = this->dshape_eulerian(s, psi, dpsidx);
605 
606  // Loop over the test functions and derivatives and set them equal to the
607  // shape functions
608  for (unsigned i = 0; i < NNODE_1D; i++)
609  {
610  test[i] = psi[i];
611  for (unsigned j = 0; j < DIM; j++)
612  {
613  dtestdx(i, j) = dpsidx(i, j);
614  }
615  }
616 
617  // Return the jacobian
618  return J;
619  }
620 
621 
622  //======================================================================
627  //======================================================================
628  template<unsigned DIM, unsigned NNODE_1D>
631  Shape& psi,
632  DShape& dpsidx,
633  Shape& test,
634  DShape& dtestdx) const
635  {
636  // Call the geometrical shape functions and derivatives
637  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
638 
639  // Set the test functions equal to the shape functions
640  //(sets internal pointers)
641  test = psi;
642  dtestdx = dpsidx;
643 
644  // Return the jacobian
645  return J;
646  }
647 
648 
651 
652 
653  //=======================================================================
658  //=======================================================================
659  template<unsigned DIM, unsigned NNODE_1D>
661  : public virtual QElement<DIM - 1, NNODE_1D>
662  {
663  public:
666  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
667  };
668 
672 
673 
674  //=======================================================================
676  //=======================================================================
677  template<unsigned NNODE_1D>
678  class FaceGeometry<QUnsteadyHeatElement<1, NNODE_1D>>
679  : public virtual PointElement
680  {
681  public:
685  };
686 
687 
691 
692 
693  //==========================================================
695  //==========================================================
696  template<class UNSTEADY_HEAT_ELEMENT>
698  : public virtual ProjectableElement<UNSTEADY_HEAT_ELEMENT>
699  {
700  public:
704 
709  {
710 #ifdef PARANOID
711  if (fld != 0)
712  {
713  std::stringstream error_stream;
714  error_stream << "UnsteadyHeat elements only store a single field so "
715  "fld must be 0 rather"
716  << " than " << fld << std::endl;
717  throw OomphLibError(
718  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
719  }
720 #endif
721 
722  // Create the vector
723  unsigned nnod = this->nnode();
724  Vector<std::pair<Data*, unsigned>> data_values(nnod);
725 
726  // Loop over all nodes
727  for (unsigned j = 0; j < nnod; j++)
728  {
729  // Add the data value associated field: The node itself
730  data_values[j] = std::make_pair(this->node_pt(j), fld);
731  }
732 
733  // Return the vector
734  return data_values;
735  }
736 
739  {
740  return 1;
741  }
742 
745  unsigned nhistory_values_for_projection(const unsigned& fld)
746  {
747 #ifdef PARANOID
748  if (fld != 0)
749  {
750  std::stringstream error_stream;
751  error_stream << "UnsteadyHeat elements only store a single field so "
752  "fld must be 0 rather"
753  << " than " << fld << std::endl;
754  throw OomphLibError(
755  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
756  }
757 #endif
758  return this->node_pt(0)->ntstorage();
759  }
760 
764  {
765  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
766  }
767 
770  double jacobian_and_shape_of_field(const unsigned& fld,
771  const Vector<double>& s,
772  Shape& psi)
773  {
774 #ifdef PARANOID
775  if (fld != 0)
776  {
777  std::stringstream error_stream;
778  error_stream << "UnsteadyHeat elements only store a single field so "
779  "fld must be 0 rather"
780  << " than " << fld << std::endl;
781  throw OomphLibError(
782  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
783  }
784 #endif
785  unsigned n_dim = this->dim();
786  unsigned n_node = this->nnode();
787  Shape test(n_node);
788  DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
789  double J =
790  this->dshape_and_dtest_eulerian_ust_heat(s, psi, dpsidx, test, dtestdx);
791  return J;
792  }
793 
794 
797  double get_field(const unsigned& t,
798  const unsigned& fld,
799  const Vector<double>& s)
800  {
801 #ifdef PARANOID
802  if (fld != 0)
803  {
804  std::stringstream error_stream;
805  error_stream << "UnsteadyHeat elements only store a single field so "
806  "fld must be 0 rather"
807  << " than " << fld << std::endl;
808  throw OomphLibError(
809  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
810  }
811 #endif
812  // Find the index at which the variable is stored
813  unsigned u_nodal_index = this->u_index_ust_heat();
814 
815  // Local shape function
816  unsigned n_node = this->nnode();
817  Shape psi(n_node);
818 
819  // Find values of shape function
820  this->shape(s, psi);
821 
822  // Initialise value of u
823  double interpolated_u = 0.0;
824 
825  // Sum over the local nodes
826  for (unsigned l = 0; l < n_node; l++)
827  {
828  interpolated_u += this->nodal_value(t, l, u_nodal_index) * psi[l];
829  }
830  return interpolated_u;
831  }
832 
833 
835  unsigned nvalue_of_field(const unsigned& fld)
836  {
837 #ifdef PARANOID
838  if (fld != 0)
839  {
840  std::stringstream error_stream;
841  error_stream << "UnsteadyHeat elements only store a single field so "
842  "fld must be 0 rather"
843  << " than " << fld << std::endl;
844  throw OomphLibError(
845  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
846  }
847 #endif
848  return this->nnode();
849  }
850 
851 
853  int local_equation(const unsigned& fld, const unsigned& j)
854  {
855 #ifdef PARANOID
856  if (fld != 0)
857  {
858  std::stringstream error_stream;
859  error_stream << "UnsteadyHeat elements only store a single field so "
860  "fld must be 0 rather"
861  << " than " << fld << std::endl;
862  throw OomphLibError(
863  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
864  }
865 #endif
866  const unsigned u_nodal_index = this->u_index_ust_heat();
867  return this->nodal_local_eqn(j, u_nodal_index);
868  }
869 
870 
873  void output(std::ostream& outfile, const unsigned& nplot)
874  {
875  unsigned el_dim = this->dim();
876  // Vector of local coordinates
878 
879  // Tecplot header info
880  outfile << this->tecplot_zone_string(nplot);
881 
882  // Loop over plot points
883  unsigned num_plot_points = this->nplot_points(nplot);
884  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
885  {
886  // Get local coordinates of plot point
887  this->get_s_plot(iplot, nplot, s);
888 
889  for (unsigned i = 0; i < el_dim; i++)
890  {
891  outfile << this->interpolated_x(s, i) << " ";
892  }
893  outfile << this->interpolated_u_ust_heat(s) << " ";
894  outfile << this->interpolated_du_dt_ust_heat(s) << " ";
895 
896 
897  // History values of coordinates
898  unsigned n_prev =
900  for (unsigned t = 1; t < n_prev; t++)
901  {
902  for (unsigned i = 0; i < el_dim; i++)
903  {
904  outfile << this->interpolated_x(t, s, i) << " ";
905  }
906  }
907 
908  // History values of velocities
909  n_prev = this->node_pt(0)->time_stepper_pt()->ntstorage();
910  for (unsigned t = 1; t < n_prev; t++)
911  {
912  outfile << this->interpolated_u_ust_heat(t, s) << " ";
913  }
914  outfile << std::endl;
915  }
916 
917 
918  // Write tecplot footer (e.g. FE connectivity lists)
919  this->write_tecplot_zone_footer(outfile, nplot);
920  }
921  };
922 
923 
924  //=======================================================================
927  //=======================================================================
928  template<class ELEMENT>
930  : public virtual FaceGeometry<ELEMENT>
931  {
932  public:
933  FaceGeometry() : FaceGeometry<ELEMENT>() {}
934  };
935 
936 
937  //=======================================================================
940  //=======================================================================
941  template<class ELEMENT>
943  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
944  {
945  public:
947  };
948 
949 
950 } // namespace oomph
951 
952 #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: shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: unsteady_heat_elements.h:933
FaceGeometry()
Definition: unsteady_heat_elements.h:684
FaceGeometry()
Definition: unsteady_heat_elements.h:666
Definition: elements.h:4998
Definition: elements.h:1313
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 std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
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
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: projection.h:183
UnsteadyHeat upgraded to become projectable.
Definition: unsteady_heat_elements.h:699
ProjectableUnsteadyHeatElement()
Definition: unsteady_heat_elements.h:703
unsigned nhistory_values_for_coordinate_projection()
Definition: unsteady_heat_elements.h:763
void output(std::ostream &outfile, const unsigned &nplot)
Definition: unsteady_heat_elements.h:873
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: unsteady_heat_elements.h:770
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
Definition: unsteady_heat_elements.h:738
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: unsteady_heat_elements.h:835
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: unsteady_heat_elements.h:708
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: unsteady_heat_elements.h:745
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: unsteady_heat_elements.h:853
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: unsteady_heat_elements.h:797
Definition: Qelements.h:459
Definition: unsteady_heat_elements.h:482
QUnsteadyHeatElement()
Definition: unsteady_heat_elements.h:491
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: unsteady_heat_elements.h:544
void output(std::ostream &outfile)
Definition: unsteady_heat_elements.h:512
double dshape_and_dtest_eulerian_ust_heat(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: unsteady_heat_elements.h:597
QUnsteadyHeatElement(const QUnsteadyHeatElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(FILE *file_pt)
Definition: unsteady_heat_elements.h:528
double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: unsteady_heat_elements.h:630
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: unsteady_heat_elements.h:505
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: unsteady_heat_elements.h:555
static const unsigned Initial_Nvalue
Definition: unsteady_heat_elements.h:486
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: unsteady_heat_elements.h:520
void output(FILE *file_pt, const unsigned &n_plot)
Definition: unsteady_heat_elements.h:536
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
Definition: unsteady_heat_elements.h:48
virtual UnsteadyHeatSourceFctPt & source_fct_pt()=0
Access function: Pointer to source function.
void(* UnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Definition: unsteady_heat_elements.h:52
Definition: unsteady_heat_elements.h:72
void compute_norm(double &norm)
Compute norm of fe solution.
Definition: unsteady_heat_elements.cc:215
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: unsteady_heat_elements.cc:472
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: unsteady_heat_elements.h:280
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute element residual Vector and element Jacobian matrix (wrapper)
Definition: unsteady_heat_elements.h:324
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: unsteady_heat_elements.h:167
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: unsteady_heat_elements.h:179
void disable_ALE()
Definition: unsteady_heat_elements.h:148
double interpolated_u_ust_heat(const unsigned &t, const Vector< double > &s) const
Definition: unsteady_heat_elements.h:362
void enable_ALE()
Definition: unsteady_heat_elements.h:158
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute element residual Vector (wrapper)
Definition: unsteady_heat_elements.h:314
double interpolated_u_ust_heat(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: unsteady_heat_elements.h:333
unsigned self_test()
Self-test: Return 0 for OK.
Definition: unsteady_heat_elements.cc:264
virtual double dshape_and_dtest_eulerian_ust_heat(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
UnsteadyHeatEquations()
Definition: unsteady_heat_elements.h:85
UnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
Definition: unsteady_heat_elements.h:446
const double & alpha() const
Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:255
static double Default_beta_parameter
Default value for Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:466
virtual void get_source_ust_heat(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: unsteady_heat_elements.h:237
UnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: unsteady_heat_elements.h:222
void(* UnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Definition: unsteady_heat_elements.h:76
UnsteadyHeatSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: unsteady_heat_elements.h:229
UnsteadyHeatEquations(const UnsteadyHeatEquations &dummy)=delete
Broken copy constructor.
virtual double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
double *& alpha_pt()
Pointer to Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:261
double interpolated_du_dt_ust_heat(const Vector< double > &s) const
Definition: unsteady_heat_elements.h:392
const double & beta() const
Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:268
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points.
Definition: unsteady_heat_elements.cc:366
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: unsteady_heat_elements.cc:60
double *& beta_pt()
Pointer to Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:274
bool ALE_is_disabled
Definition: unsteady_heat_elements.h:451
double * Beta_pt
Pointer to Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:457
static double Default_alpha_parameter
2D UnsteadyHeat elements
Definition: unsteady_heat_elements.h:462
double du_dt_ust_heat(const unsigned &n) const
Definition: unsteady_heat_elements.h:119
virtual unsigned u_index_ust_heat() const
Broken assignment operator.
Definition: unsteady_heat_elements.h:112
double * Alpha_pt
Pointer to Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:454
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
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
unsigned el_dim
dimension
Definition: overloaded_cartesian_element_body.h:30
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2