oomph::IMRODEElement Class Reference

#include <ode_element_for_imr.h>

+ Inheritance diagram for oomph::IMRODEElement:

Public Member Functions

 IMRODEElement (TimeStepper *time_stepper_pt, SolutionFunctorBase *exact_solution_pt)
 Constructor. More...
 
virtual ~IMRODEElement ()
 Virtual destructor. More...
 
Vector< doubletime_interpolate_u () const
 
double time_interpolate_time () const
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
- Public Member Functions inherited from oomph::ODEElement
 ODEElement ()
 Default constructor: null any pointers. More...
 
 ODEElement (TimeStepper *time_stepper_pt, SolutionFunctorBase *exact_solution_pt)
 
void build (TimeStepper *time_stepper_pt, SolutionFunctorBase *exact_solution_pt)
 Store pointers, create internal data. More...
 
virtual ~ODEElement ()
 
unsigned nvalue () const
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mm)
 
Vector< doubleexact_solution (const double &t) const
 Exact solution. More...
 
Vector< doublederivative_function (const double &t, const Vector< double > &u)
 Exact solution. More...
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual unsigned self_test ()
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 

Private Member Functions

 IMRODEElement (const IMRODEElement &dummy)
 Broken copy constructor. More...
 

Additional Inherited Members

- Public Attributes inherited from oomph::ODEElement
SolutionFunctorBaseExact_solution_pt
 
bool Use_fd_jacobian
 
- Static Public Attributes inherited from oomph::GeneralisedElement
static bool Suppress_warning_about_repeated_internal_data
 
static bool Suppress_warning_about_repeated_external_data = true
 
static double Default_fd_jacobian_step = 1.0e-8
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

Class for ODE elements when using the implicit midpoint rule time integrator IMR. The residual and Jacobian functions are modified to interpolate in time. Note that the normal element can be used for the IMR implementation in the class IMRByBDF.

Constructor & Destructor Documentation

◆ IMRODEElement() [1/2]

oomph::IMRODEElement::IMRODEElement ( TimeStepper time_stepper_pt,
SolutionFunctorBase exact_solution_pt 
)
inline

Constructor.

45  : ODEElement(time_stepper_pt, exact_solution_pt) {}
ODEElement()
Default constructor: null any pointers.
Definition: ode_elements.h:44

◆ ~IMRODEElement()

virtual oomph::IMRODEElement::~IMRODEElement ( )
inlinevirtual

Virtual destructor.

48 {}

◆ IMRODEElement() [2/2]

oomph::IMRODEElement::IMRODEElement ( const IMRODEElement dummy)
inlineprivate

Broken copy constructor.

193  {BrokenCopy::broken_copy("IMRODEElement");}
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Definition: oomph_utilities.cc:212

References oomph::BrokenCopy::broken_copy().

Member Function Documentation

◆ fill_in_contribution_to_jacobian()

void oomph::IMRODEElement::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Modified get jacobian: two differences from the normal verson here which are required when using the IMR implementation: 1) We have to interpolate in time to calculate the analytical contribution. 2) The interpolation in time within the residual calculation modifies the Jacobian (because the jacobian is the derivative w.r.t. the n+1-th values).

Reimplemented from oomph::ODEElement.

140  {
141  // Get residuals
143 
145  {
146  const unsigned n = nvalue();
147  TimeStepper* ts_pt = internal_data_pt(0)->time_stepper_pt();
148 
149  // Get values **interpolated in time**
150  Vector<double> u = time_interpolate_u();
151 
152  // Get **interpolated** continuous time
153  double t = time_interpolate_time();
154 
155  // get df/du jacobian analytically
156  Vector<double> dummy;
157  Exact_solution_pt->jacobian(t, dummy, u, jacobian);
158 
159  // Multiply by the weight of the history value 0 (i.e. at time n+1) in
160  // the calculation of the time-interpolated "current" value (0th
161  // derivative). Need to do this because of the chain rule, just write
162  // out the derivation of dr/du for nvalue=1 to see it.
163  const double w00 = ts_pt->weight(0,0);
164  for(unsigned i=0; i<n; i++)
165  {
166  for(unsigned j=0; j<n; j++)
167  {
168  jacobian(i, j) *= w00;
169  }
170  }
171 
172  // We actually need jacobian of residual = f(imr_t, imr_u) - dudt so
173  // subtract diagonal (dudt)/du term now.
174  const double w10 = ts_pt->weight(1,0);
175  for(unsigned i=0; i<n; i++)
176  {
177  jacobian(i, i) -= w10;
178  }
179  }
180  else
181  {
182  // Use FD for jacobian
184  (residuals, jacobian, true);
185  }
186  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
void fill_in_jacobian_from_internal_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
Definition: elements.cc:1102
Vector< double > time_interpolate_u() const
Definition: ode_element_for_imr.h:53
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Definition: ode_element_for_imr.h:104
double time_interpolate_time() const
Definition: ode_element_for_imr.h:80
bool Use_fd_jacobian
Definition: ode_elements.h:182
SolutionFunctorBase * Exact_solution_pt
Definition: ode_elements.h:180
unsigned nvalue() const
Definition: ode_elements.h:76
virtual void jacobian(const double &t, const Vector< double > &x, const Vector< double > &u, DenseMatrix< double > &jacobian) const
Definition: oomph_utilities.h:1142
virtual bool have_jacobian() const
Is a jacobian function implemented?
Definition: oomph_utilities.h:1153
t
Definition: plotPSD.py:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::ODEElement::Exact_solution_pt, fill_in_contribution_to_residuals(), oomph::GeneralisedElement::fill_in_jacobian_from_internal_by_fd(), oomph::SolutionFunctorBase::have_jacobian(), i, oomph::GeneralisedElement::internal_data_pt(), j, oomph::SolutionFunctorBase::jacobian(), n, oomph::ODEElement::nvalue(), plotPSD::t, time_interpolate_time(), time_interpolate_u(), oomph::Data::time_stepper_pt(), oomph::ODEElement::Use_fd_jacobian, and oomph::TimeStepper::weight().

◆ fill_in_contribution_to_residuals()

void oomph::IMRODEElement::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Modified get residuals: have to interpolate in time for IMR when using the IMR implementation so we can't use the normal ODE element implementation.

Reimplemented from oomph::ODEElement.

105  {
106 
107  // Get pointer to one-and-only internal data object
108  Data* dat_pt = internal_data_pt(0);
109 
110  // Get time stepper
111  TimeStepper* ts_pt = dat_pt->time_stepper_pt();
112 
113  // Get values **interpolated in time**
114  Vector<double> u = time_interpolate_u();
115 
116  // Get **interpolated** continuous time
117  double t = time_interpolate_time();
118 
119  Vector<double> deriv = derivative_function(t, u);
120  for(unsigned j=0, nj=deriv.size(); j<nj; j++)
121  {
122  // Get dudt approximation from time stepper
123  double dudt = ts_pt->time_derivative(1, dat_pt, j);
124 
125  // Residual is difference between the exact derivative and our
126  // time stepper's derivative estimate.
127  residuals[j] = deriv[j] - dudt;
128  }
129  }
Vector< double > derivative_function(const double &t, const Vector< double > &u)
Exact solution.
Definition: ode_elements.h:166

References oomph::ODEElement::derivative_function(), oomph::GeneralisedElement::internal_data_pt(), j, plotPSD::t, oomph::TimeStepper::time_derivative(), time_interpolate_time(), time_interpolate_u(), and oomph::Data::time_stepper_pt().

Referenced by fill_in_contribution_to_jacobian().

◆ time_interpolate_time()

double oomph::IMRODEElement::time_interpolate_time ( ) const
inline

For IMR implemented using IMR we have to interpolate everything in time (to get the midpoint values). This function does that for the continuous time.

81  {
82 
83  // Get pointer to one-and-only internal data object
84  Data* dat_pt = internal_data_pt(0);
85 
86  // Get time stepper
87  TimeStepper* ts_pt = dat_pt->time_stepper_pt();
88 
89  // Number of history values to interpolate over for value/time calculation
90  const unsigned nj = ts_pt->nprev_values_for_value_at_evaluation_time();
91 
92  double t = 0;
93  for(unsigned j=0; j<nj; j++)
94  {
95  t += ts_pt->time_pt()->time(j) * ts_pt->weight(0, j);
96  }
97 
98  return t;
99  }

References oomph::GeneralisedElement::internal_data_pt(), j, oomph::TimeStepper::nprev_values_for_value_at_evaluation_time(), plotPSD::t, oomph::Time::time(), oomph::TimeStepper::time_pt(), oomph::Data::time_stepper_pt(), and oomph::TimeStepper::weight().

Referenced by fill_in_contribution_to_jacobian(), and fill_in_contribution_to_residuals().

◆ time_interpolate_u()

Vector<double> oomph::IMRODEElement::time_interpolate_u ( ) const
inline

For IMR implemented using IMR we have to interpolate everything in time (to get the midpoint values). This function does that for the u values.

54  {
55 
56  // Get pointer to one-and-only internal data object
57  Data* dat_pt = internal_data_pt(0);
58 
59  // Get time stepper
60  TimeStepper* ts_pt = dat_pt->time_stepper_pt();
61 
62  // Number of history values to interpolate over for value/time calculation
63  const unsigned nj = ts_pt->nprev_values_for_value_at_evaluation_time();
64 
65  Vector<double> u(nvalue(), 0.0);
66  for(unsigned j=0; j<nj; j++)
67  {
68  for(unsigned i=0; i<nvalue(); i++)
69  {
70  u[i] += dat_pt->value(j, i) * ts_pt->weight(0, j);
71  }
72  }
73 
74  return u;
75  }

References i, oomph::GeneralisedElement::internal_data_pt(), j, oomph::TimeStepper::nprev_values_for_value_at_evaluation_time(), oomph::ODEElement::nvalue(), oomph::Data::time_stepper_pt(), oomph::Data::value(), and oomph::TimeStepper::weight().

Referenced by fill_in_contribution_to_jacobian(), and fill_in_contribution_to_residuals().


The documentation for this class was generated from the following file: