space_time_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-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 // Header file for SpaceTimeUnsteadyHeat elements
27 #ifndef OOMPH_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
28 #define OOMPH_SPACE_TIME_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 // Oomph-lib headers
36 #include "generic/Qelements.h"
37 #include "generic/shape.h"
38 #include "generic/projection.h"
39 
40 /// //////////////////////////////////////////////////////////////////////
41 /// //////////////////////////////////////////////////////////////////////
42 /// //////////////////////////////////////////////////////////////////////
43 
44 /// DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
45 namespace oomph
46 {
47  //============================================================================
48  /// Base class so that we don't need to know the dimension just to
49  /// set the source function!
50  //============================================================================
52  {
53  public:
54  /// Function pointer to source function fct(t,x,f(x,t)) -- x
55  /// is a Vector!
56  typedef void (*SpaceTimeUnsteadyHeatSourceFctPt)(const double& time,
57  const Vector<double>& x,
58  double& u);
59 
60  /// Access function: Pointer to source function
62  };
63 
64  //============================================================================
65  /// A class for all isoparametric elements that solve the
66  /// SpaceTimeUnsteadyHeat equations.
67  /// \f[ \frac{\partial^2 u}{\partial x_i^2}=\frac{\partial u}{\partial t}+f(t,x_j) \f]
68  /// This contains the generic maths. Shape functions, geometric
69  /// mapping etc. must get implemented in derived class.
70  /// Note that this class assumes an isoparametric formulation, i.e. that
71  /// the scalar unknown is interpolated using the same shape funcitons
72  /// as the position.
73  //============================================================================
74  template<unsigned SPATIAL_DIM>
76  : public virtual SpaceTimeUnsteadyHeatEquationsBase
77  {
78  public:
79  /// Function pointer to source function fct(t,x,f(x,t)) -- x
80  /// is a Vector!
81  /// DRAIG: Why is this here? There is already one in the base class!
82  typedef void (*SpaceTimeUnsteadyHeatSourceFctPt)(const double& time,
83  const Vector<double>& x,
84  double& u);
85 
86 
87  /// Constructor: Initialises the Source_fct_pt to null and sets
88  /// flag to use ALE formulation of the equations. Also, set Alpha (thermal
89  /// inertia) and Beta (thermal conductivity) parameters to defaults (both
90  /// one for natural scaling).
92  {
93  // Set Alpha parameter to default (one for natural scaling)
95 
96  // Set Beta parameter to default (one for natural scaling)
98  } // End of SpaceTimeUnsteadyHeatEquations
99 
100 
101  /// Broken copy constructor
103  const SpaceTimeUnsteadyHeatEquations& dummy) = delete;
104 
105  /// Disable ALE, i.e. assert the mesh is not moving -- you do this
106  /// at your own risk!
107  void disable_ALE()
108  {
109  // Set the flag to true
110  ALE_is_disabled = true;
111  } // End of disable_ALE
112 
113 
114  /// (Re-)enable ALE, i.e. take possible mesh motion into account
115  /// when evaluating the time-derivative. Note: By default, ALE is
116  /// enabled, at the expense of possibly creating unnecessary work
117  /// in problems where the mesh is, in fact, stationary.
118  void enable_ALE()
119  {
120  // Set the flag to false
121  ALE_is_disabled = false;
122  } // End of enable_ALE
123 
124 
125  /// Compute norm of FE solution
126  void compute_norm(double& norm);
127 
128 
129  /// Output with default number of plot points
130  void output(std::ostream& outfile)
131  {
132  // Number of plot points
133  unsigned nplot = 5;
134 
135  // Output the solution
136  output(outfile, nplot);
137  } // End of output
138 
139 
140  /// Output FE representation of soln: x,y,u or x,y,z,u at
141  /// nplot^SPATIAL_DIM plot points
142  void output(std::ostream& outfile, const unsigned& nplot);
143 
144 
145  /// C_style output with default number of plot points
146  void output(FILE* file_pt)
147  {
148  // Number of plot points
149  unsigned nplot = 5;
150 
151  // Output the solution
152  output(file_pt, nplot);
153  } // End of output
154 
155 
156  /// C-style output FE representation of soln: x,y,u or x,y,z,u at
157  /// nplot^SPATIAL_DIM plot points
158  void output(FILE* file_pt, const unsigned& nplot);
159 
160 
161  /// Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^SPATIAL_DIM
162  /// plot points
163  void output_fct(std::ostream& outfile,
164  const unsigned& nplot,
166 
167 
168  /// Output exact soln: x,y,u_exact or x,y,z,u_exact at
169  /// nplot^SPATIAL_DIM plot points (time-dependent version)
170  virtual void output_fct(
171  std::ostream& outfile,
172  const unsigned& nplot,
173  const double& time,
175 
176 
177  /// Get error and norm against exact solution
178  void compute_error(std::ostream& outfile,
180  double& error,
181  double& norm);
182 
183 
184  /// Get error and norm against exact solution
185  void compute_error(std::ostream& outfile,
187  const double& time,
188  double& error,
189  double& norm);
190 
191 
192  /// C-style output FE representation of soln: x,y,u or x,y,z,u at
193  /// nplot^SPATIAL_DIM plot points
194  void output_element_paraview(std::ofstream& outfile, const unsigned& nplot);
195 
196 
197  /// Number of scalars/fields output by this element. Reimplements
198  /// broken virtual function in base class.
199  unsigned nscalar_paraview() const
200  {
201  // Only one field to output
202  return 1;
203  } // End of nscalar_paraview
204 
205 
206  /// Write values of the i-th scalar field at the plot points. Needs
207  /// to be implemented for each new specific element type.
208  void scalar_value_paraview(std::ofstream& file_out,
209  const unsigned& i,
210  const unsigned& nplot) const
211  {
212 #ifdef PARANOID
213  if (i != 0)
214  {
215  std::stringstream error_stream;
216  error_stream << "Space-time unsteady heat elements only store a single "
217  << "field so i must be 0 rather than " << i << std::endl;
218  throw OomphLibError(
219  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
220  }
221 #endif
222 
223  // Get the number of plot points
224  unsigned local_loop = this->nplot_points_paraview(nplot);
225 
226  // Loop over the plot points
227  for (unsigned j = 0; j < local_loop; j++)
228  {
229  // Storage for the local coordinates
230  Vector<double> s(SPATIAL_DIM + 1);
231 
232  // Get the local coordinate of the required plot point
233  this->get_s_plot(j, nplot, s);
234 
235  // Output the interpolated solution value
236  file_out << this->interpolated_u_ust_heat(s) << std::endl;
237  }
238  } // End of scalar_value_paraview
239 
240 
241  /// Write values of the i-th scalar field at the plot points. Needs
242  /// to be implemented for each new specific element type.
244  std::ofstream& file_out,
245  const unsigned& i,
246  const unsigned& nplot,
247  FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
248  {
249 #ifdef PARANOID
250  if (i != 0)
251  {
252  std::stringstream error_stream;
253  error_stream << "Space-time unsteady heat elements only store a single "
254  << "field so i must be 0 rather than " << i << std::endl;
255  throw OomphLibError(
256  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
257  }
258 #endif
259 
260  // Get the number of plot points
261  unsigned local_loop = this->nplot_points_paraview(nplot);
262 
263  // Loop over the plot points
264  for (unsigned j = 0; j < local_loop; j++)
265  {
266  // Storage for the local coordinates
267  Vector<double> s(SPATIAL_DIM + 1);
268 
269  // Storage for the global coordinates
270  Vector<double> spatial_coordinates(SPATIAL_DIM);
271 
272  // Get the local coordinate of the required plot point
273  this->get_s_plot(j, nplot, s);
274 
275  // Loop over the spatial coordinates
276  for (unsigned i = 0; i < SPATIAL_DIM; i++)
277  {
278  // Assign the i-th spatial coordinate
279  spatial_coordinates[i] = interpolated_x(s, i);
280  }
281 
282  // Exact solution vector (here it's simply a scalar)
283  Vector<double> exact_soln(1, 0.0);
284 
285  // Get the exact solution at this point
286  (*exact_soln_pt)(spatial_coordinates, exact_soln);
287 
288  // Output the interpolated solution value
289  file_out << exact_soln[0] << std::endl;
290  } // for (unsigned j=0;j<local_loop;j++)
291  } // End of scalar_value_fct_paraview
292 
293 
294  /// Write values of the i-th scalar field at the plot points. Needs
295  /// to be implemented for each new specific element type.
297  std::ofstream& file_out,
298  const unsigned& i,
299  const unsigned& nplot,
300  const double& time,
301  FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
302  {
303 #ifdef PARANOID
304  if (i != 0)
305  {
306  std::stringstream error_stream;
307  error_stream << "Space-time unsteady heat elements only store a single "
308  << "field so i must be 0 rather than " << i << std::endl;
309  throw OomphLibError(
310  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
311  }
312 #endif
313 
314  // Get the number of plot points
315  unsigned local_loop = this->nplot_points_paraview(nplot);
316 
317  // Loop over the plot points
318  for (unsigned j = 0; j < local_loop; j++)
319  {
320  // Storage for the local coordinates
321  Vector<double> s(SPATIAL_DIM + 1);
322 
323  // Storage for the time value
324  double interpolated_t = 0.0;
325 
326  // Storage for the global coordinates
327  Vector<double> spatial_coordinates(SPATIAL_DIM);
328 
329  // Get the local coordinate of the required plot point
330  this->get_s_plot(j, nplot, s);
331 
332  // Loop over the spatial coordinates
333  for (unsigned i = 0; i < SPATIAL_DIM; i++)
334  {
335  // Assign the i-th spatial coordinate
336  spatial_coordinates[i] = interpolated_x(s, i);
337  }
338 
339  // Get the time value
340  interpolated_t = interpolated_x(s, SPATIAL_DIM);
341 
342  // Exact solution vector (here it's simply a scalar)
343  Vector<double> exact_soln(1, 0.0);
344 
345  // Get the exact solution at this point
346  (*exact_soln_pt)(interpolated_t, spatial_coordinates, exact_soln);
347 
348  // Output the interpolated solution value
349  file_out << exact_soln[0] << std::endl;
350  } // for (unsigned j=0;j<local_loop;j++)
351  } // End of scalar_value_fct_paraview
352 
353 
354  /// Name of the i-th scalar field. Default implementation
355  /// returns V1 for the first one, V2 for the second etc.
356  std::string scalar_name_paraview(const unsigned& i) const
357  {
358  // If we're outputting the solution
359  if (i == 0)
360  {
361  // There's only one field to output
362  return "U";
363  }
364  // Never get here
365  else
366  {
367  std::stringstream error_stream;
368  error_stream << "These unsteady heat elements only store 1 field, \n"
369  << "but i is currently " << i << std::endl;
370  throw OomphLibError(
371  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
372 
373  // Dummy return
374  return " ";
375  }
376  } // End of scalar_name_paraview
377 
378 
379  /// Access function: Pointer to source function
381  {
382  // Return the source function pointer
383  return Source_fct_pt;
384  } // End of source_fct_pt
385 
386 
387  /// Access function: Pointer to source function. Const version
389  {
390  // Return the source function pointer
391  return Source_fct_pt;
392  }
393 
394 
395  /// Get source term at continous time t and (Eulerian) position x.
396  /// Virtual so it can be overloaded in derived multi-physics elements.
397  virtual inline void get_source_ust_heat(const double& t,
398  const unsigned& ipt,
399  const Vector<double>& x,
400  double& source) const
401  {
402  // If no source function has been set, return zero
403  if (Source_fct_pt == 0)
404  {
405  // Set the source term value to zero
406  source = 0.0;
407  }
408  // Otherwise return the appropriate value
409  else
410  {
411  // Get source strength
412  (*Source_fct_pt)(t, x, source);
413  }
414  } // End of get_source_ust_heat
415 
416 
417  /// Alpha parameter (thermal inertia)
418  const double& alpha() const
419  {
420  // Return the value of Alpha
421  return *Alpha_pt;
422  } // End of alpha
423 
424 
425  /// Pointer to Alpha parameter (thermal inertia)
426  double*& alpha_pt()
427  {
428  // Return the pointer to Alpha
429  return Alpha_pt;
430  } // End of alpha_pt
431 
432 
433  /// Beta parameter (thermal conductivity)
434  const double& beta() const
435  {
436  // Return the pointer to Beta
437  return *Beta_pt;
438  } // End of beta
439 
440 
441  /// Pointer to Beta parameter (thermal conductivity)
442  double*& beta_pt()
443  {
444  // Return the pointer to Beta
445  return Beta_pt;
446  } // End of beta_pt
447 
448 
449  /// Get flux: flux[i]=du/dx_i
450  void get_flux(const Vector<double>& s, Vector<double>& flux) const
451  {
452  // Find out how many nodes there are in the element
453  unsigned n_node = nnode();
454 
455  // Find the index at which the variable is stored
456  unsigned u_nodal_index = u_index_ust_heat();
457 
458  // Set up memory for the shape and test functions
459  Shape psi(n_node);
460 
461  // Set up memory for the derivatives of the shape and test functions
462  DShape dpsidx(n_node, SPATIAL_DIM + 1);
463 
464  // Call the derivatives of the shape and test functions
465  dshape_eulerian(s, psi, dpsidx);
466 
467  // Loop over the entries of the flux vector
468  for (unsigned j = 0; j < SPATIAL_DIM; j++)
469  {
470  // Initialise j-th flux entry to zero
471  flux[j] = 0.0;
472  }
473 
474  // Loop over nodes
475  for (unsigned l = 0; l < n_node; l++)
476  {
477  // Loop over derivative directions
478  for (unsigned j = 0; j < SPATIAL_DIM; j++)
479  {
480  // Update the flux value
481  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
482  }
483  } // for (unsigned l=0;l<n_node;l++)
484  } // End of get_flux
485 
486 
487  /// Compute element residual Vector (wrapper)
489  {
490  // Call the generic residuals function with flag set to 0
491  // using a dummy matrix argument
493  residuals, GeneralisedElement::Dummy_matrix, 0);
494  } // End of fill_in_contribution_to_residuals
495 
496 
497  /// Compute element residual Vector and element Jacobian matrix (wrapper)
499  DenseMatrix<double>& jacobian)
500  {
501  // Call the generic routine with the flag set to 1
502  fill_in_generic_residual_contribution_ust_heat(residuals, jacobian, 1);
503  } // End of fill_in_contribution_to_jacobian
504 
505 
506  /// Return FE representation of function value u(s) at local coordinate s
507  inline double interpolated_u_ust_heat(const Vector<double>& s) const
508  {
509  // Find number of nodes
510  unsigned n_node = nnode();
511 
512  // Find the index at which the variable is stored
513  unsigned u_nodal_index = u_index_ust_heat();
514 
515  // Local shape function
516  Shape psi(n_node);
517 
518  // Find values of the shape functions at local coordinate s
519  shape(s, psi);
520 
521  // Initialise value of u
522  double interpolated_u = 0.0;
523 
524  // Loop over the local nodes and sum
525  for (unsigned l = 0; l < n_node; l++)
526  {
527  // Update the interpolated u value
528  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
529  }
530 
531  // Return the interpolated u value
532  return interpolated_u;
533  } // End of interpolated_u_ust_heat
534 
535 
536  /// Return the index at which the unknown value
537  /// is stored. The default value, 0, is appropriate for single-physics
538  /// problems, when there is only one variable, the value that satisfies the
539  /// unsteady heat equation.
540  /// In derived multi-physics elements, this function should be overloaded
541  /// to reflect the chosen storage scheme. Note that these equations require
542  /// that the unknown is always stored at the same index at each node.
543  virtual inline unsigned u_index_ust_heat() const
544  {
545  // Return the default value
546  return 0;
547  } // End of u_index_ust_heat
548 
549 
550  /// Return FE representation of function value u(s) at local
551  /// coordinate s at previous time t (t=0: present)
552  /// DRAIG: This needs to be broken; doesn't make sense in space-time
553  /// elements!
554  inline double interpolated_u_ust_heat(const unsigned& t,
555  const Vector<double>& s) const
556  {
557  // Find number of nodes
558  unsigned n_node = nnode();
559 
560  // Find the index at which the variable is stored
561  unsigned u_nodal_index = u_index_ust_heat();
562 
563  // Local shape function
564  Shape psi(n_node);
565 
566  // Find values of shape function
567  shape(s, psi);
568 
569  // Initialise value of u
570  double interpolated_u = 0.0;
571 
572  // Loop over the local nodes and sum
573  for (unsigned l = 0; l < n_node; l++)
574  {
575  // Update the interpolated u value
576  interpolated_u += nodal_value(t, l, u_nodal_index) * psi[l];
577  }
578 
579  // Return the interpolated u value
580  return interpolated_u;
581  } // End of interpolated_u_ust_heat
582 
583 
584  /// Calculate du/dt at the n-th local node. Uses suitably
585  /// interpolated value for hanging nodes.
586  double du_dt_ust_heat(const unsigned& n) const
587  {
588  // Storage for the local coordinates
589  Vector<double> s(SPATIAL_DIM + 1, 0.0);
590 
591  // Get the local coordinate at the n-th node
593 
594  // Return the interpolated du/dt value
596  } // End of du_dt_ust_heat
597 
598 
599  /// Return FE representation of function value du/dt(s) at local coordinate
600  /// s
601  inline double interpolated_du_dt_ust_heat(const Vector<double>& s) const
602  {
603  // Find number of nodes
604  unsigned n_node = nnode();
605 
606  // Find the index at which the variable is stored
607  unsigned u_nodal_index = u_index_ust_heat();
608 
609  // Local shape function
610  Shape psi(n_node);
611 
612  // Allocate space for the derivatives of the shape functions
613  DShape dpsidx(n_node, SPATIAL_DIM + 1);
614 
615  // Compute the geometric shape functions and also first derivatives
616  // w.r.t. global coordinates at local coordinate s
617  dshape_eulerian(s, psi, dpsidx);
618 
619  // Initialise value of du/dt
620  double interpolated_dudt = 0.0;
621 
622  // Loop over the local nodes and sum
623  for (unsigned l = 0; l < n_node; l++)
624  {
625  // Update the interpolated du/dt value
626  interpolated_dudt +=
627  nodal_value(l, u_nodal_index) * dpsidx(l, SPATIAL_DIM);
628  }
629 
630  // Return the interpolated du/dt value
631  return interpolated_dudt;
632  } // End of interpolated_du_dt_ust_heat
633 
634 
635  /// Self-test: Return 0 for OK
636  unsigned self_test();
637 
638  protected:
639  /// Shape/test functions and derivs w.r.t. to global coords at
640  /// local coordinate s; return Jacobian of mapping
642  const Vector<double>& s,
643  Shape& psi,
644  DShape& dpsidx,
645  Shape& test,
646  DShape& dtestdx) const = 0;
647 
648 
649  /// Shape/test functions and derivs w.r.t. to global coords at
650  /// integration point ipt; return Jacobian of mapping
652  const unsigned& ipt,
653  Shape& psi,
654  DShape& dpsidx,
655  Shape& test,
656  DShape& dtestdx) const = 0;
657 
658 
659  /// Compute element residual Vector only (if flag=and/or element
660  /// Jacobian matrix
662  Vector<double>& residuals,
663  DenseMatrix<double>& jacobian,
664  const unsigned& flag);
665 
666  /// Pointer to source function:
668 
669  /// Boolean flag to indicate if ALE formulation is disabled when
670  /// time-derivatives are computed. Only set to true if you're sure that
671  /// the mesh is stationary.
673 
674  /// Pointer to Alpha parameter (thermal inertia)
675  double* Alpha_pt;
676 
677  /// Pointer to Beta parameter (thermal conductivity)
678  double* Beta_pt;
679 
680  private:
681  /// Static default value for the Alpha parameter (thermal inertia):
682  /// One for natural scaling
683  static double Default_alpha_parameter;
684 
685  /// Static default value for the Beta parameter (thermal
686  /// conductivity): One for natural scaling
687  static double Default_beta_parameter;
688  };
689 
690 
691  /// ////////////////////////////////////////////////////////////////////////
692  /// ////////////////////////////////////////////////////////////////////////
693  /// ////////////////////////////////////////////////////////////////////////
694 
695 
696  //=========================================================================
697  /// QUnsteadyHeatSpaceTimeElement elements are quadrilateral/brick-
698  /// shaped UnsteadyHeat elements with isoparametric interpolation for
699  /// the function.
700  //=========================================================================
701  template<unsigned SPATIAL_DIM, unsigned NNODE_1D>
703  : public virtual QElement<SPATIAL_DIM + 1, NNODE_1D>,
704  public virtual SpaceTimeUnsteadyHeatEquations<SPATIAL_DIM>
705  {
706  public:
707  /// Constructor: Call constructors for QElement and
708  /// SpaceTimeUnsteadyHeat equations
710  : QElement<SPATIAL_DIM + 1, NNODE_1D>(),
711  SpaceTimeUnsteadyHeatEquations<SPATIAL_DIM>()
712  {
713  }
714 
715  /// Broken copy constructor
718  delete;
719 
720  /// Required number of 'values' (pinned or dofs) at node n
721  inline unsigned required_nvalue(const unsigned& n) const
722  {
723  // Return the appropriate value
724  return Initial_Nvalue;
725  } // End of required_nvalue
726 
727 
728  /// Output function:
729  /// x,t,u or x,y,t,u
730  void output(std::ostream& outfile)
731  {
732  // Call the function in the base class
734  } // End of output
735 
736 
737  /// Output function:
738  /// x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points
739  void output(std::ostream& outfile, const unsigned& n_plot)
740  {
741  // Call the function in the base class
743  } // End of output
744 
745 
746  /// C-style output function:
747  /// x,t,u or x,y,t,u
748  void output(FILE* file_pt)
749  {
750  // Call the function in the base class
752  } // End of output
753 
754 
755  /// C-style output function:
756  /// x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points
757  void output(FILE* file_pt, const unsigned& n_plot)
758  {
759  // Call the function in the base class
761  } // End of output
762 
763 
764  /// Output function for an exact solution:
765  /// x,t,u_exact or x,y,t,u_exact at n_plot^(SPATIAL_DIM+1) plot points
766  void output_fct(std::ostream& outfile,
767  const unsigned& n_plot,
769  {
770  // Call the function in the base class
772  outfile, n_plot, exact_soln_pt);
773  } // End of output_fct
774 
775 
776  /// Output function for a time-dependent exact solution.
777  /// x,t,u_exact or x,y,t,u_exact at n_plot^(SPATIAL_DIM+1) plot points
778  /// (Calls the unsteady version)
779  void output_fct(std::ostream& outfile,
780  const unsigned& n_plot,
781  const double& time,
783  {
784  // Call the function in the base class
786  outfile, n_plot, time, exact_soln_pt);
787  } // End of output_fct
788 
789  protected:
790  /// Shape/test functions & derivs. w.r.t. to global coords. Return Jacobian.
792  Shape& psi,
793  DShape& dpsidx,
794  Shape& test,
795  DShape& dtestdx) const;
796 
797 
798  /// Shape/test functions and derivs w.r.t. to global coords at
799  /// integration point ipt; return Jacobian of mapping
801  const unsigned& ipt,
802  Shape& psi,
803  DShape& dpsidx,
804  Shape& test,
805  DShape& dtestdx) const;
806 
807  private:
808  /// Static array of ints to hold number of variables at nodes:
809  /// Initial_Nvalue[n]
810  static const unsigned Initial_Nvalue;
811  };
812 
813 
814  //======================================================================
815  /// Define the shape functions and test functions and derivatives
816  /// w.r.t. global coordinates and return Jacobian of mapping.
817  ///
818  /// Galerkin: Test functions=shape functions
819  //======================================================================
820  template<unsigned SPATIAL_DIM, unsigned NNODE_1D>
823  Shape& psi,
824  DShape& dpsidx,
825  Shape& test,
826  DShape& dtestdx) const
827  {
828  // Call the geometrical shape functions and derivatives
829  double det = this->dshape_eulerian(s, psi, dpsidx);
830 
831  // The test functions are equal to the shape functions
832  test = psi;
833 
834  // The test function derivatives are equal to those of the shape functions
835  dtestdx = dpsidx;
836 
837  // Return the Jacobian of the mapping
838  return det;
839  } // End of dshape_and_dtest_eulerian_ust_heat
840 
841 
842  //======================================================================
843  /// Define the shape functions and test functions and derivatives
844  /// w.r.t. global coordinates and return Jacobian of mapping.
845  ///
846  /// Galerkin: Test functions=shape functions
847  //======================================================================
848  template<unsigned SPATIAL_DIM, unsigned NNODE_1D>
851  Shape& psi,
852  DShape& dpsidx,
853  Shape& test,
854  DShape& dtestdx) const
855  {
856  // Find the element dimension
857  const unsigned el_dim = SPATIAL_DIM + 1;
858 
859  // Storage for the local coordinates of the integration point
860  Vector<double> s(el_dim, 0.0);
861 
862  // Set the local coordinate
863  for (unsigned i = 0; i < el_dim; i++)
864  {
865  // Calculate the i-th local coordinate at the ipt-th knot point
866  s[i] = this->integral_pt()->knot(ipt, i);
867  }
868 
869  // Return the Jacobian of the geometrical shape functions and derivatives
870  return dshape_and_dtest_eulerian_ust_heat(s, psi, dpsidx, test, dtestdx);
871  } // End of dshape_and_dtest_eulerian_at_knot_ust_heat
872 
873 
874  /// /////////////////////////////////////////////////////////////////////
875  /// /////////////////////////////////////////////////////////////////////
876  /// /////////////////////////////////////////////////////////////////////
877 
878 
879  //=======================================================================
880  /// Face geometry for the QUnsteadyHeatSpaceTimeElement elements: The
881  /// spatial dimension of the face elements is one lower than that of
882  /// the bulk element but they have the same number of points along their
883  /// 1D edges.
884  //=======================================================================
885  template<unsigned SPATIAL_DIM, unsigned NNODE_1D>
886  class FaceGeometry<QUnsteadyHeatSpaceTimeElement<SPATIAL_DIM, NNODE_1D>>
887  : public virtual QElement<SPATIAL_DIM, NNODE_1D>
888  {
889  public:
890  /// Constructor: Call the constructor for the appropriate
891  /// lower-dimensional QElement
892  FaceGeometry() : QElement<SPATIAL_DIM, NNODE_1D>() {}
893  };
894 
895 
896  /// /////////////////////////////////////////////////////////////////////
897  /// /////////////////////////////////////////////////////////////////////
898  /// /////////////////////////////////////////////////////////////////////
899 
900 
901  //=======================================================================
902  /// Face geometry for the 1D QUnsteadyHeatSpaceTimeElement elements:
903  /// Point elements
904  //=======================================================================
905  template<unsigned NNODE_1D>
907  : public virtual PointElement
908  {
909  public:
910  /// Constructor: Call the constructor for the appropriate
911  /// lower-dimensional QElement
913  };
914 
915 
916  /// /////////////////////////////////////////////////////////////////////
917  /// /////////////////////////////////////////////////////////////////////
918  /// /////////////////////////////////////////////////////////////////////
919 
920 
921  //==========================================================
922  /// SpaceTimeUnsteadyHeat upgraded to become projectable
923  //==========================================================
924  template<class UNSTEADY_HEAT_ELEMENT>
926  : public virtual ProjectableElement<UNSTEADY_HEAT_ELEMENT>
927  {
928  public:
929  /// Constructor [this was only required explicitly
930  /// from gcc 4.5.2 onwards...]
932 
933 
934  /// Specify the values associated with field fld. The information
935  /// is returned in a vector of pairs which comprise the Data object and
936  /// the value within it, that correspond to field fld.
938  {
939 #ifdef PARANOID
940  // If we're not dealing with the first field
941  if (fld != 0)
942  {
943  // Create a stringstream object to create an error message
944  std::stringstream error_stream;
945 
946  // Create the error string
947  error_stream << "SpaceTimeUnsteadyHeat elements only store a single "
948  << "field so fld must be 0 rather than " << fld
949  << std::endl;
950 
951  // Throw an error
952  throw OomphLibError(
953  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
954  }
955 #endif
956 
957  // The number of nodes in this element
958  unsigned nnod = this->nnode();
959 
960  // Storage for the pairs
961  Vector<std::pair<Data*, unsigned>> data_values(nnod);
962 
963  // Loop over all nodes
964  for (unsigned j = 0; j < nnod; j++)
965  {
966  // Add the data value and associated field: The node itself
967  data_values[j] = std::make_pair(this->node_pt(j), fld);
968  }
969 
970  // Return the vector
971  return data_values;
972  } // End of data_values_of_field
973 
974 
975  /// Number of fields to be projected: Just one
977  {
978  // Return the appropriate value
979  return 1;
980  } // End of nfields_for_projection
981 
982 
983  /// Number of history values to be stored for fld-th field.
984  /// (Note: count includes current value!)
985  unsigned nhistory_values_for_projection(const unsigned& fld)
986  {
987 #ifdef PARANOID
988  // If we're not dealing with the first field
989  if (fld != 0)
990  {
991  // Create a stringstream object to create an error message
992  std::stringstream error_stream;
993 
994  // Create the error string
995  error_stream << "SpaceTimeUnsteadyHeat elements only store a single "
996  << "field so fld must be 0 rather than " << fld
997  << std::endl;
998 
999  // Throw an error
1000  throw OomphLibError(
1001  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1002  }
1003 #endif
1004 
1005  // Return the number of stored values
1006  return this->node_pt(0)->ntstorage();
1007  } // End of nhistory_values_for_projection
1008 
1009 
1010  /// Number of positional history values (Note: count includes
1011  /// current value!)
1013  {
1014  // Return the number of history values stored by the position timestepper
1015  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
1016  } // End of nhistory_values_for_coordinate_projection
1017 
1018 
1019  /// Return Jacobian of mapping and shape functions of field fld
1020  /// at local coordinate s
1021  double jacobian_and_shape_of_field(const unsigned& fld,
1022  const Vector<double>& s,
1023  Shape& psi)
1024  {
1025 #ifdef PARANOID
1026  // If we're not dealing with the first field
1027  if (fld != 0)
1028  {
1029  // Create a stringstream object to create an error message
1030  std::stringstream error_stream;
1031 
1032  // Create the error string
1033  error_stream << "SpaceTimeUnsteadyHeat elements only store a single "
1034  << "field so fld must be 0 rather than " << fld
1035  << std::endl;
1036 
1037  // Throw an error
1038  throw OomphLibError(
1039  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1040  }
1041 #endif
1042 
1043  // Get the number of dimensions in the element
1044  unsigned n_dim = this->dim();
1045 
1046  // Get the number of nodes in the element
1047  unsigned n_node = this->nnode();
1048 
1049  // Allocate space for the test functions
1050  Shape test(n_node);
1051 
1052  // Allocate space for the derivatives of the shape functions
1053  DShape dpsidx(n_node, n_dim);
1054 
1055  // Allocate space for the derivatives of the test functions
1056  DShape dtestdx(n_node, n_dim);
1057 
1058  // Calculate the shape functions and their derivatives at the local
1059  // coordinate s (and the same for the test functions). On top of this
1060  // calculate the determinant of the Jacobian
1061  double J =
1062  this->dshape_and_dtest_eulerian_ust_heat(s, psi, dpsidx, test, dtestdx);
1063 
1064  // Return the determinant of the Jacobian
1065  return J;
1066  } // End of jacobian_and_shape_of_field
1067 
1068 
1069  /// Return interpolated field fld at local coordinate s, at time
1070  /// level t (t=0: present; t>0: history values)
1071  double get_field(const unsigned& t,
1072  const unsigned& fld,
1073  const Vector<double>& s)
1074  {
1075 #ifdef PARANOID
1076  // If we're not dealing with the first field
1077  if (fld != 0)
1078  {
1079  // Create a stringstream object to create an error message
1080  std::stringstream error_stream;
1081 
1082  // Create the error string
1083  error_stream << "SpaceTimeUnsteadyHeat elements only store a single "
1084  << "field so fld must be 0 rather than " << fld
1085  << std::endl;
1086 
1087  // Throw an error
1088  throw OomphLibError(
1089  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1090  }
1091 #endif
1092 
1093  // Find the index at which the variable is stored
1094  unsigned u_nodal_index = this->u_index_ust_heat();
1095 
1096  // Get the number of nodes in the element
1097  unsigned n_node = this->nnode();
1098 
1099  // Local shape function
1100  Shape psi(n_node);
1101 
1102  // Find values of shape function
1103  this->shape(s, psi);
1104 
1105  // Initialise value of u
1106  double interpolated_u = 0.0;
1107 
1108  // Loop over the local nodes
1109  for (unsigned l = 0; l < n_node; l++)
1110  {
1111  // Update the interpolated solution value
1112  interpolated_u += this->nodal_value(t, l, u_nodal_index) * psi[l];
1113  }
1114 
1115  // Return the interpolated solution value
1116  return interpolated_u;
1117  } // End of get_field
1118 
1119 
1120  /// Return number of values in field fld: One per node
1121  unsigned nvalue_of_field(const unsigned& fld)
1122  {
1123 #ifdef PARANOID
1124  // If we're not dealing with the first field
1125  if (fld != 0)
1126  {
1127  // Create a stringstream object to create an error message
1128  std::stringstream error_stream;
1129 
1130  // Create the error string
1131  error_stream << "SpaceTimeUnsteadyHeat elements only store a single "
1132  << "field so fld must be 0 rather than " << fld
1133  << std::endl;
1134 
1135  // Throw an error
1136  throw OomphLibError(
1137  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1138  }
1139 #endif
1140 
1141  // Return the number of nodes in the element
1142  return this->nnode();
1143  } // End of nvalue_of_field
1144 
1145 
1146  /// Return local equation number of value j in field fld.
1147  int local_equation(const unsigned& fld, const unsigned& j)
1148  {
1149 #ifdef PARANOID
1150  // If we're not dealing with the first field
1151  if (fld != 0)
1152  {
1153  // Create a stringstream object to create an error message
1154  std::stringstream error_stream;
1155 
1156  // Create the error string
1157  error_stream << "SpaceTimeUnsteadyHeat elements only store a single "
1158  << "field so fld must be 0 rather than " << fld
1159  << std::endl;
1160 
1161  // Throw an error
1162  throw OomphLibError(
1163  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1164  }
1165 #endif
1166 
1167  // Get the nodal index of the unknown
1168  const unsigned u_nodal_index = this->u_index_ust_heat();
1169 
1170  // Output the local equation number
1171  return this->nodal_local_eqn(j, u_nodal_index);
1172  } // End of local_equation
1173 
1174 
1175  /// Output FE representation of soln: x,t,u or x,y,t,u
1176  /// at n_plot^(SPATIAL_DIM+1) plot points
1177  void output(std::ostream& outfile, const unsigned& nplot)
1178  {
1179  // Get the dimension of the element
1180  unsigned el_dim = this->dim();
1181 
1182  // Vector of local coordinates
1183  Vector<double> s(el_dim, 0.0);
1184 
1185  // Tecplot header info
1186  outfile << this->tecplot_zone_string(nplot);
1187 
1188  // Get the number of plot points
1189  unsigned num_plot_points = this->nplot_points(nplot);
1190 
1191  // Loop over plot points
1192  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
1193  {
1194  // Get local coordinates of plot point
1195  this->get_s_plot(iplot, nplot, s);
1196 
1197  // Loop over the coordinate directions
1198  for (unsigned i = 0; i < el_dim; i++)
1199  {
1200  // Output the interpolated coordinates
1201  outfile << this->interpolated_x(s, i) << " ";
1202  }
1203 
1204  // Output the interpolated value of u(s)
1205  outfile << this->interpolated_u_ust_heat(s) << " ";
1206 
1207  // Output the interpolated value of du/dt(s)
1208  outfile << this->interpolated_du_dt_ust_heat(s) << " ";
1209 
1210  // History values of coordinates
1211  unsigned n_prev =
1213 
1214  // Loop over the previous timesteps
1215  for (unsigned t = 1; t < n_prev; t++)
1216  {
1217  // Loop over the coordinate directions
1218  for (unsigned i = 0; i < el_dim; i++)
1219  {
1220  // Output the coordinates
1221  outfile << this->interpolated_x(t, s, i) << " ";
1222  }
1223  } // for (unsigned t=1;t<n_prev;t++)
1224 
1225  // Number of history values of velocities
1226  n_prev = this->node_pt(0)->time_stepper_pt()->ntstorage();
1227 
1228  // Loop over the previous timesteps
1229  for (unsigned t = 1; t < n_prev; t++)
1230  {
1231  // Output the solution
1232  outfile << this->interpolated_u_ust_heat(t, s) << " ";
1233  }
1234 
1235  // Finish the line
1236  outfile << std::endl;
1237  } // for (unsigned iplot=0;iplot<num_plot_points;iplot++)
1238 
1239  // Write tecplot footer (e.g. FE connectivity lists)
1240  this->write_tecplot_zone_footer(outfile, nplot);
1241  } // End of output
1242  };
1243 
1244 
1245  //=======================================================================
1246  /// Face geometry for element is the same as that for the underlying
1247  /// wrapped element
1248  //=======================================================================
1249  template<class ELEMENT>
1251  : public virtual FaceGeometry<ELEMENT>
1252  {
1253  public:
1254  FaceGeometry() : FaceGeometry<ELEMENT>() {}
1255  };
1256 
1257 
1258  //=======================================================================
1259  /// Face geometry of the Face Geometry for element is the same as
1260  /// that for the underlying wrapped element
1261  //=======================================================================
1262  template<class ELEMENT>
1265  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
1266  {
1267  public:
1269  };
1270 } // End of namespace oomph
1271 
1272 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
Definition: nodes.cc:879
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
A general Finite Element class.
Definition: elements.h:1317
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot. Broken virtual; can b...
Definition: elements.h:2866
virtual void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size (broken virtual)
Definition: elements.h:1846
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
Definition: elements.h:3165
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
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:3992
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
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2615
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
Definition: elements.h:3152
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
Definition: elements.h:3190
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3328
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
Definition: elements.h:3178
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition: elements.h:1769
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
Definition: elements.h:227
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3443
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Definition: projection.h:183
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (Note: count includes current value!...
void output(std::ostream &outfile, const unsigned &nplot)
Output FE representation of soln: x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
ProjectableUnsteadyHeatSpaceTimeElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
QUnsteadyHeatSpaceTimeElement(const QUnsteadyHeatSpaceTimeElement< SPATIAL_DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
double dshape_and_dtest_eulerian_ust_heat(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape/test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points.
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. x,t,u_exact or x,y,t,u_exact at n_plot^(SPATIAL_...
QUnsteadyHeatSpaceTimeElement()
Constructor: Call constructors for QElement and SpaceTimeUnsteadyHeat equations.
void output(FILE *file_pt)
C-style output function: x,t,u or x,y,t,u.
void output(std::ostream &outfile)
Output function: x,t,u or x,y,t,u.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,t,u_exact or x,y,t,u_exact at n_plot^(SPATIAL_DIM+1) plot po...
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
unsigned required_nvalue(const unsigned &n) const
Required number of 'values' (pinned or dofs) at node n.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
Base class so that we don't need to know the dimension just to set the source function!
void(* SpaceTimeUnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Function pointer to source function fct(t,x,f(x,t)) – x is a Vector!
virtual SpaceTimeUnsteadyHeatSourceFctPt & source_fct_pt()=0
Access function: Pointer to source function.
A class for all isoparametric elements that solve the SpaceTimeUnsteadyHeat equations.
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^SPATIAL_DIM plot points.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error and norm against exact solution.
double *& beta_pt()
Pointer to Beta parameter (thermal conductivity)
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute element residual Vector (wrapper)
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
SpaceTimeUnsteadyHeatEquations(const SpaceTimeUnsteadyHeatEquations &dummy)=delete
Broken copy constructor.
void(* SpaceTimeUnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Function pointer to source function fct(t,x,f(x,t)) – x is a Vector! DRAIG: Why is this here?...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute element residual Vector and element Jacobian matrix (wrapper)
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
double interpolated_u_ust_heat(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void output(FILE *file_pt)
C_style output with default number of plot points.
const double & alpha() const
Alpha parameter (thermal inertia)
double * Beta_pt
Pointer to Beta parameter (thermal conductivity)
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
SpaceTimeUnsteadyHeatSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
void compute_norm(double &norm)
Compute norm of FE solution.
double * Alpha_pt
Pointer to Alpha parameter (thermal inertia)
unsigned self_test()
Self-test: Return 0 for OK.
double *& alpha_pt()
Pointer to Alpha parameter (thermal inertia)
static double Default_beta_parameter
Static default value for the Beta parameter (thermal conductivity): One for natural scaling.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
virtual double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void output(std::ostream &outfile)
Output with default number of plot points.
double interpolated_du_dt_ust_heat(const Vector< double > &s) const
Return FE representation of function value du/dt(s) at local coordinate s.
const double & beta() const
Beta parameter (thermal conductivity)
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
virtual double dshape_and_dtest_eulerian_ust_heat(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at local coordinate s; return Jacobian of map...
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
void output_element_paraview(std::ofstream &outfile, const unsigned &nplot)
C-style output FE representation of soln: x,y,u or x,y,z,u at nplot^SPATIAL_DIM plot points.
virtual unsigned u_index_ust_heat() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
SpaceTimeUnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
double interpolated_u_ust_heat(const unsigned &t, const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s at previous time t (t=0: presen...
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
double du_dt_ust_heat(const unsigned &n) const
Calculate du/dt at the n-th local node. Uses suitably interpolated value for hanging nodes.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i]=du/dx_i.
SpaceTimeUnsteadyHeatEquations()
Constructor: Initialises the Source_fct_pt to null and sets flag to use ALE formulation of the equati...
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
SpaceTimeUnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
static double Default_alpha_parameter
Static default value for the Alpha parameter (thermal inertia): One for natural scaling.
virtual void get_source_ust_heat(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at continous time t and (Eulerian) position x. Virtual so it can be overloaded in der...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
Definition: timesteppers.h:601
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...