poisson_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 Poisson elements
27 #ifndef OOMPH_POISSON_ELEMENTS_HEADER
28 #define OOMPH_POISSON_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 <sstream>
37 
38 // OOMPH-LIB headers
39 #include "../generic/projection.h"
40 #include "../generic/nodes.h"
41 #include "../generic/Qelements.h"
42 #include "../generic/oomph_utilities.h"
43 
44 
45 namespace oomph
46 {
47  //=============================================================
53  //=============================================================
54  template<unsigned DIM>
55  class PoissonEquations : public virtual FiniteElement
56  {
57  public:
60  typedef void (*PoissonSourceFctPt)(const Vector<double>& x, double& f);
61 
62 
66  Vector<double>& gradient);
67 
68 
71 
73  PoissonEquations(const PoissonEquations& dummy) = delete;
74 
76  void operator=(const PoissonEquations&) = delete;
77 
85  virtual inline unsigned u_index_poisson() const
86  {
87  return 0;
88  }
89 
90 
94  {
95  // Dimension
96  unsigned dim = s.size();
97 
98  // Resize data for values
99  data.resize(dim + 1);
100 
101  // Write values in the vector
102  for (unsigned i = 0; i < dim; i++)
103  {
104  data[i] = interpolated_x(s, i);
105  }
106  data[dim] = this->interpolated_u_poisson(s);
107  }
108 
109 
112  unsigned nscalar_paraview() const
113  {
114  return 1;
115  }
116 
119  void scalar_value_paraview(std::ofstream& file_out,
120  const unsigned& i,
121  const unsigned& nplot) const
122  {
123 #ifdef PARANOID
124  if (i != 0)
125  {
126  std::stringstream error_stream;
127  error_stream
128  << "Poisson elements only store a single field so i must be 0 rather"
129  << " than " << i << std::endl;
130  throw OomphLibError(
131  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
132  }
133 #endif
134 
135  unsigned local_loop = this->nplot_points_paraview(nplot);
136  for (unsigned j = 0; j < local_loop; j++)
137  {
138  // Get the local coordinate of the required plot point
140  this->get_s_plot(j, nplot, s);
141 
142  file_out << this->interpolated_u_poisson(s) << std::endl;
143  }
144  }
145 
149  std::string scalar_name_paraview(const unsigned& i) const
150  {
151 #ifdef PARANOID
152  if (i != 0)
153  {
154  std::stringstream error_stream;
155  error_stream
156  << "Poisson elements only store a single field so i must be 0 rather"
157  << " than " << i << std::endl;
158  throw OomphLibError(
159  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
160  }
161 #endif
162 
163  return "Poisson solution";
164  }
165 
166 
170  std::ofstream& file_out,
171  const unsigned& i,
172  const unsigned& nplot,
173  FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
174  {
175 #ifdef PARANOID
176  if (i != 0)
177  {
178  std::stringstream error_stream;
179  error_stream << "Poisson equation has only one field. Can't call "
180  << "this function for value " << i << std::endl;
181  throw OomphLibError(
182  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
183  }
184 #endif
185 
186  // Vector of local coordinates
188 
189  // Vector for coordinates
191 
192  // Exact solution Vector
194 
195  // Loop over plot points
196  unsigned num_plot_points = nplot_points_paraview(nplot);
197  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
198  {
199  // Get local coordinates of plot point
200  get_s_plot(iplot, nplot, s);
201 
202  // Get x position as Vector
203  interpolated_x(s, x);
204 
205  // Get exact solution at this point
206  (*exact_soln_pt)(x, exact_soln);
207 
208  // Write it
209  file_out << exact_soln[0] << std::endl;
210  }
211  }
212 
214  void output(std::ostream& outfile)
215  {
216  const unsigned n_plot = 5;
217  output(outfile, n_plot);
218  }
219 
222  void output(std::ostream& outfile, const unsigned& n_plot);
223 
225  void output(FILE* file_pt)
226  {
227  const unsigned n_plot = 5;
228  output(file_pt, n_plot);
229  }
230 
233  void output(FILE* file_pt, const unsigned& n_plot);
234 
237  void output_fct(std::ostream& outfile,
238  const unsigned& n_plot,
240 
244  virtual void output_fct(
245  std::ostream& outfile,
246  const unsigned& n_plot,
247  const double& time,
249  {
250  throw OomphLibError(
251  "There is no time-dependent output_fct() for Poisson elements ",
254  }
255 
256 
258  void compute_norm(double& norm);
259 
261  void compute_error(std::ostream& outfile,
263  double& error,
264  double& norm);
265 
266 
268  void compute_error(std::ostream& outfile,
270  const double& time,
271  double& error,
272  double& norm)
273  {
274  throw OomphLibError(
275  "There is no time-dependent compute_error() for Poisson elements",
278  }
279 
282  {
283  return Source_fct_pt;
284  }
285 
288  {
289  return Source_fct_pt;
290  }
291 
294  {
295  return Source_fct_gradient_pt;
296  }
297 
300  {
301  return Source_fct_gradient_pt;
302  }
303 
304 
309  inline virtual void get_source_poisson(const unsigned& ipt,
310  const Vector<double>& x,
311  double& source) const
312  {
313  // If no source function has been set, return zero
314  if (Source_fct_pt == 0)
315  {
316  source = 0.0;
317  }
318  else
319  {
320  // Get source strength
321  (*Source_fct_pt)(x, source);
322  }
323  }
324 
325 
331  inline virtual void get_source_gradient_poisson(
332  const unsigned& ipt,
333  const Vector<double>& x,
334  Vector<double>& gradient) const
335  {
336  // If no gradient function has been set, FD it
337  if (Source_fct_gradient_pt == 0)
338  {
339  // Reference value
340  double source = 0.0;
341  get_source_poisson(ipt, x, source);
342 
343  // FD it
345  double source_pls = 0.0;
346  Vector<double> x_pls(x);
347  for (unsigned i = 0; i < DIM; i++)
348  {
349  x_pls[i] += eps_fd;
350  get_source_poisson(ipt, x_pls, source_pls);
351  gradient[i] = (source_pls - source) / eps_fd;
352  x_pls[i] = x[i];
353  }
354  }
355  else
356  {
357  // Get gradient
358  (*Source_fct_gradient_pt)(x, gradient);
359  }
360  }
361 
362 
365  {
366  // Find out how many nodes there are in the element
367  const unsigned n_node = nnode();
368 
369  // Get the index at which the unknown is stored
370  const unsigned u_nodal_index = u_index_poisson();
371 
372  // Set up memory for the shape and test functions
373  Shape psi(n_node);
374  DShape dpsidx(n_node, DIM);
375 
376  // Call the derivatives of the shape and test functions
377  dshape_eulerian(s, psi, dpsidx);
378 
379  // Initialise to zero
380  for (unsigned j = 0; j < DIM; j++)
381  {
382  flux[j] = 0.0;
383  }
384 
385  // Loop over nodes
386  for (unsigned l = 0; l < n_node; l++)
387  {
388  // Loop over derivative directions
389  for (unsigned j = 0; j < DIM; j++)
390  {
391  flux[j] += this->nodal_value(l, u_nodal_index) * dpsidx(l, j);
392  }
393  }
394  }
395 
396 
400  Vector<Vector<double>>& dflux_dnodal_u) const
401  {
402  // Find out how many nodes there are in the element
403  const unsigned n_node = nnode();
404 
405  // Set up memory for the shape and test functions
406  Shape psi(n_node);
407  DShape dpsidx(n_node, DIM);
408 
409  // Call the derivatives of the shape and test functions
410  dshape_eulerian(s, psi, dpsidx);
411 
412  // And here it is...
413  for (unsigned i = 0; i < DIM; i++)
414  {
415  for (unsigned j = 0; j < n_node; j++)
416  {
417  dflux_dnodal_u[i][j] = dpsidx(j, i);
418  }
419  }
420  }
421 
422 
425  {
426  // Call the generic residuals function with flag set to 0
427  // using a dummy matrix argument
429  residuals, GeneralisedElement::Dummy_matrix, 0);
430  }
431 
432 
436  DenseMatrix<double>& jacobian)
437  {
438  // Call the generic routine with the flag set to 1
439  fill_in_generic_residual_contribution_poisson(residuals, jacobian, 1);
440  }
441 
442 
445  virtual inline double interpolated_u_poisson(const Vector<double>& s) const
446  {
447  // Find number of nodes
448  const unsigned n_node = nnode();
449 
450  // Get the index at which the poisson unknown is stored
451  const unsigned u_nodal_index = u_index_poisson();
452 
453  // Local shape function
454  Shape psi(n_node);
455 
456  // Find values of shape function
457  shape(s, psi);
458 
459  // Initialise value of u
460  double interpolated_u = 0.0;
461 
462  // Loop over the local nodes and sum
463  for (unsigned l = 0; l < n_node; l++)
464  {
465  interpolated_u += this->nodal_value(l, u_nodal_index) * psi[l];
466  }
467 
468  return (interpolated_u);
469  }
470 
471 
477  RankThreeTensor<double>& dresidual_dnodal_coordinates);
478 
480  unsigned self_test();
481 
482 
483  protected:
487  Shape& psi,
488  DShape& dpsidx,
489  Shape& test,
490  DShape& dtestdx) const = 0;
491 
492 
496  const unsigned& ipt,
497  Shape& psi,
498  DShape& dpsidx,
499  Shape& test,
500  DShape& dtestdx) const = 0;
501 
506  const unsigned& ipt,
507  Shape& psi,
508  DShape& dpsidx,
509  RankFourTensor<double>& d_dpsidx_dX,
510  Shape& test,
511  DShape& dtestdx,
512  RankFourTensor<double>& d_dtestdx_dX,
513  DenseMatrix<double>& djacobian_dX) const = 0;
514 
518  Vector<double>& residuals,
519  DenseMatrix<double>& jacobian,
520  const unsigned& flag);
521 
524 
527  };
528 
529 
533 
534 
535  //======================================================================
538  //======================================================================
539  template<unsigned DIM, unsigned NNODE_1D>
540  class QPoissonElement : public virtual QElement<DIM, NNODE_1D>,
541  public virtual PoissonEquations<DIM>
542  {
543  private:
546  static const unsigned Initial_Nvalue;
547 
548  public:
552 
555 
558 
561  inline unsigned required_nvalue(const unsigned& n) const
562  {
563  return Initial_Nvalue;
564  }
565 
568  void output(std::ostream& outfile)
569  {
571  }
572 
573 
576  void output(std::ostream& outfile, const unsigned& n_plot)
577  {
578  PoissonEquations<DIM>::output(outfile, n_plot);
579  }
580 
581 
584  void output(FILE* file_pt)
585  {
587  }
588 
589 
592  void output(FILE* file_pt, const unsigned& n_plot)
593  {
594  PoissonEquations<DIM>::output(file_pt, n_plot);
595  }
596 
597 
600  void output_fct(std::ostream& outfile,
601  const unsigned& n_plot,
603  {
604  PoissonEquations<DIM>::output_fct(outfile, n_plot, exact_soln_pt);
605  }
606 
607 
611  void output_fct(std::ostream& outfile,
612  const unsigned& n_plot,
613  const double& time,
615  {
616  PoissonEquations<DIM>::output_fct(outfile, n_plot, time, exact_soln_pt);
617  }
618 
619 
620  protected:
624  Shape& psi,
625  DShape& dpsidx,
626  Shape& test,
627  DShape& dtestdx) const;
628 
629 
633  const unsigned& ipt,
634  Shape& psi,
635  DShape& dpsidx,
636  Shape& test,
637  DShape& dtestdx) const;
638 
643  const unsigned& ipt,
644  Shape& psi,
645  DShape& dpsidx,
646  RankFourTensor<double>& d_dpsidx_dX,
647  Shape& test,
648  DShape& dtestdx,
649  RankFourTensor<double>& d_dtestdx_dX,
650  DenseMatrix<double>& djacobian_dX) const;
651  };
652 
653 
654  // Inline functions:
655 
656 
657  //======================================================================
662  //======================================================================
663  template<unsigned DIM, unsigned NNODE_1D>
665  const Vector<double>& s,
666  Shape& psi,
667  DShape& dpsidx,
668  Shape& test,
669  DShape& dtestdx) const
670  {
671  // Call the geometrical shape functions and derivatives
672  const double J = this->dshape_eulerian(s, psi, dpsidx);
673 
674  // Set the test functions equal to the shape functions
675  test = psi;
676  dtestdx = dpsidx;
677 
678  // Return the jacobian
679  return J;
680  }
681 
682 
683  //======================================================================
688  //======================================================================
689  template<unsigned DIM, unsigned NNODE_1D>
692  Shape& psi,
693  DShape& dpsidx,
694  Shape& test,
695  DShape& dtestdx) const
696  {
697  // Call the geometrical shape functions and derivatives
698  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
699 
700  // Set the pointers of the test functions
701  test = psi;
702  dtestdx = dpsidx;
703 
704  // Return the jacobian
705  return J;
706  }
707 
708 
709  //======================================================================
716  //======================================================================
717  template<unsigned DIM, unsigned NNODE_1D>
720  const unsigned& ipt,
721  Shape& psi,
722  DShape& dpsidx,
723  RankFourTensor<double>& d_dpsidx_dX,
724  Shape& test,
725  DShape& dtestdx,
726  RankFourTensor<double>& d_dtestdx_dX,
727  DenseMatrix<double>& djacobian_dX) const
728  {
729  // Call the geometrical shape functions and derivatives
730  const double J = this->dshape_eulerian_at_knot(
731  ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
732 
733  // Set the pointers of the test functions
734  test = psi;
735  dtestdx = dpsidx;
736  d_dtestdx_dX = d_dpsidx_dX;
737 
738  // Return the jacobian
739  return J;
740  }
741 
742 
746 
747 
748  //=======================================================================
753  //=======================================================================
754  template<unsigned DIM, unsigned NNODE_1D>
755  class FaceGeometry<QPoissonElement<DIM, NNODE_1D>>
756  : public virtual QElement<DIM - 1, NNODE_1D>
757  {
758  public:
761  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
762  };
763 
767 
768 
769  //=======================================================================
771  //=======================================================================
772  template<unsigned NNODE_1D>
773  class FaceGeometry<QPoissonElement<1, NNODE_1D>> : public virtual PointElement
774  {
775  public:
779  };
780 
781 
785 
786 
787  //==========================================================
789  //==========================================================
790  template<class POISSON_ELEMENT>
792  : public virtual ProjectableElement<POISSON_ELEMENT>
793  {
794  public:
799  {
800 #ifdef PARANOID
801  if (fld != 0)
802  {
803  std::stringstream error_stream;
804  error_stream << "Poisson elements only store a single field so fld "
805  "must be 0 rather"
806  << " than " << fld << std::endl;
807  throw OomphLibError(
808  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
809  }
810 #endif
811 
812  // Create the vector
813  unsigned nnod = this->nnode();
814  Vector<std::pair<Data*, unsigned>> data_values(nnod);
815 
816  // Loop over all nodes
817  for (unsigned j = 0; j < nnod; j++)
818  {
819  // Add the data value associated field: The node itself
820  data_values[j] = std::make_pair(this->node_pt(j), fld);
821  }
822 
823  // Return the vector
824  return data_values;
825  }
826 
829  {
830  return 1;
831  }
832 
835  unsigned nhistory_values_for_projection(const unsigned& fld)
836  {
837 #ifdef PARANOID
838  if (fld != 0)
839  {
840  std::stringstream error_stream;
841  error_stream << "Poisson elements only store a single field so fld "
842  "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->node_pt(0)->ntstorage();
849  }
850 
854  {
855  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
856  }
857 
860  double jacobian_and_shape_of_field(const unsigned& fld,
861  const Vector<double>& s,
862  Shape& psi)
863  {
864 #ifdef PARANOID
865  if (fld != 0)
866  {
867  std::stringstream error_stream;
868  error_stream << "Poisson elements only store a single field so fld "
869  "must be 0 rather"
870  << " than " << fld << std::endl;
871  throw OomphLibError(
872  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
873  }
874 #endif
875  unsigned n_dim = this->dim();
876  unsigned n_node = this->nnode();
877  Shape test(n_node);
878  DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
879  double J =
880  this->dshape_and_dtest_eulerian_poisson(s, psi, dpsidx, test, dtestdx);
881  return J;
882  }
883 
884 
887  double get_field(const unsigned& t,
888  const unsigned& fld,
889  const Vector<double>& s)
890  {
891 #ifdef PARANOID
892  if (fld != 0)
893  {
894  std::stringstream error_stream;
895  error_stream << "Poisson elements only store a single field so fld "
896  "must be 0 rather"
897  << " than " << fld << std::endl;
898  throw OomphLibError(
899  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
900  }
901 #endif
902  // Find the index at which the variable is stored
903  unsigned u_nodal_index = this->u_index_poisson();
904 
905  // Local shape function
906  unsigned n_node = this->nnode();
907  Shape psi(n_node);
908 
909  // Find values of shape function
910  this->shape(s, psi);
911 
912  // Initialise value of u
913  double interpolated_u = 0.0;
914 
915  // Sum over the local nodes
916  for (unsigned l = 0; l < n_node; l++)
917  {
918  interpolated_u += this->nodal_value(l, u_nodal_index) * psi[l];
919  }
920  return interpolated_u;
921  }
922 
923 
925  unsigned nvalue_of_field(const unsigned& fld)
926  {
927 #ifdef PARANOID
928  if (fld != 0)
929  {
930  std::stringstream error_stream;
931  error_stream << "Poisson elements only store a single field so fld "
932  "must be 0 rather"
933  << " than " << fld << std::endl;
934  throw OomphLibError(
935  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
936  }
937 #endif
938  return this->nnode();
939  }
940 
941 
943  int local_equation(const unsigned& fld, const unsigned& j)
944  {
945 #ifdef PARANOID
946  if (fld != 0)
947  {
948  std::stringstream error_stream;
949  error_stream << "Poisson elements only store a single field so fld "
950  "must be 0 rather"
951  << " than " << fld << std::endl;
952  throw OomphLibError(
953  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
954  }
955 #endif
956  const unsigned u_nodal_index = this->u_index_poisson();
957  return this->nodal_local_eqn(j, u_nodal_index);
958  }
959  };
960 
961 
962  //=======================================================================
965  //=======================================================================
966  template<class ELEMENT>
968  : public virtual FaceGeometry<ELEMENT>
969  {
970  public:
971  FaceGeometry() : FaceGeometry<ELEMENT>() {}
972  };
973 
974 
975  //=======================================================================
978  //=======================================================================
979  template<class ELEMENT>
981  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
982  {
983  public:
985  };
986 
987 
988 } // namespace oomph
989 
990 #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
int data[]
Definition: Map_placement_new.cpp:1
Definition: shape.h:278
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: poisson_elements.h:778
FaceGeometry()
Definition: poisson_elements.h:761
Definition: elements.h:4998
Definition: elements.h:1313
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862
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
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
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 double Default_fd_jacobian_step
Definition: elements.h:1198
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
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: poisson_elements.h:56
virtual void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: poisson_elements.cc:49
void point_output_data(const Vector< double > &s, Vector< double > &data)
Definition: poisson_elements.h:93
unsigned self_test()
Self-test: Return 0 for OK.
Definition: poisson_elements.cc:298
virtual void get_source_gradient_poisson(const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient) const
Definition: poisson_elements.h:331
virtual void get_source_poisson(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: poisson_elements.h:309
virtual double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const =0
PoissonSourceFctGradientPt source_fct_gradient_pt() const
Access function: Pointer to gradient source function. Const version.
Definition: poisson_elements.h:299
void get_dflux_dnodal_u(const Vector< double > &s, Vector< Vector< double >> &dflux_dnodal_u) const
Definition: poisson_elements.h:399
void(* PoissonSourceFctPt)(const Vector< double > &x, double &f)
Definition: poisson_elements.h:60
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: poisson_elements.h:214
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: poisson_elements.h:424
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: poisson_elements.h:119
PoissonEquations()
Constructor (must initialise the Source_fct_pt to null)
Definition: poisson_elements.h:70
PoissonSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: poisson_elements.h:287
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
Definition: poisson_elements.cc:447
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
Definition: poisson_elements.h:523
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: poisson_elements.h:268
PoissonSourceFctGradientPt Source_fct_gradient_pt
Pointer to gradient of source function.
Definition: poisson_elements.h:526
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: poisson_elements.cc:165
void(* PoissonSourceFctGradientPt)(const Vector< double > &x, Vector< double > &gradient)
Definition: poisson_elements.h:65
PoissonSourceFctGradientPt & source_fct_gradient_pt()
Access function: Pointer to gradient of source function.
Definition: poisson_elements.h:293
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.h:244
void operator=(const PoissonEquations &)=delete
Broken assignment operator.
virtual double dshape_and_dtest_eulerian_poisson(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
unsigned nscalar_paraview() const
Definition: poisson_elements.h:112
virtual double interpolated_u_poisson(const Vector< double > &s) const
Definition: poisson_elements.h:445
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: poisson_elements.h:364
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: poisson_elements.cc:500
std::string scalar_name_paraview(const unsigned &i) const
Definition: poisson_elements.h:149
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: poisson_elements.h:435
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.cc:400
virtual unsigned u_index_poisson() const
Definition: poisson_elements.h:85
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: poisson_elements.h:281
virtual double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
Definition: poisson_elements.h:169
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: poisson_elements.h:225
PoissonEquations(const PoissonEquations &dummy)=delete
Broken copy constructor.
Definition: projection.h:183
Poisson upgraded to become projectable.
Definition: poisson_elements.h:793
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: poisson_elements.h:925
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: poisson_elements.h:798
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
Definition: poisson_elements.h:828
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: poisson_elements.h:860
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: poisson_elements.h:887
unsigned nhistory_values_for_coordinate_projection()
Definition: poisson_elements.h:853
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: poisson_elements.h:835
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: poisson_elements.h:943
Definition: Qelements.h:459
Definition: poisson_elements.h:542
double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const
Definition: poisson_elements.h:719
double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: poisson_elements.h:691
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.h:611
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: poisson_elements.h:576
QPoissonElement(const QPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned required_nvalue(const unsigned &n) const
Definition: poisson_elements.h:561
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.h:600
QPoissonElement()
Definition: poisson_elements.h:551
void output(std::ostream &outfile)
Definition: poisson_elements.h:568
double dshape_and_dtest_eulerian_poisson(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: poisson_elements.h:664
void output(FILE *file_pt)
Definition: poisson_elements.h:584
static const unsigned Initial_Nvalue
Definition: poisson_elements.h:546
void operator=(const QPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void output(FILE *file_pt, const unsigned &n_plot)
Definition: poisson_elements.h:592
A Rank 4 Tensor class.
Definition: matrices.h:1701
A Rank 3 Tensor class.
Definition: matrices.h:1370
Definition: shape.h:76
unsigned ntstorage() const
Definition: timesteppers.h:601
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
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 exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
int error
Definition: calibrate.py:297
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
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