oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH > Class Template Reference

#include <foeppl_von_karman_volume_constraint_element.h>

+ Inheritance diagram for oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >:

Public Member Functions

 FoepplvonKarmanVolumeConstraintElement (MESH< ELEMENT > *bounding_mesh_pt, const Vector< unsigned > &contributing_region, const double &pressure=0.0)
 
 ~FoepplvonKarmanVolumeConstraintElement ()
 
 FoepplvonKarmanVolumeConstraintElement (const FoepplvonKarmanVolumeConstraintElement &)=delete
 Broken copy constructor. More...
 
void operator= (const FoepplvonKarmanVolumeConstraintElement &)=delete
 Broken assignment operator. More...
 
double get_bounded_volume ()
 
void assign_additional_local_eqn_numbers ()
 Assign the equation number for the new equation. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Fill in contribution to elemental residual and Jacobian. More...
 
void set_prescribed_volume (double *volume_pt)
 Set pointer to target volume. More...
 
Datapressure_data_pt () const
 
- 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
 

Protected Attributes

DataVolume_control_pressure_pt
 
unsigned Pressure_data_index
 
int Volume_control_local_eqn
 Local equation number of volume constraint. More...
 
MESH< ELEMENT > * Bounding_mesh_pt
 
Vector< unsignedContributing_region
 
doublePrescribed_volume_pt
 Pointer to target volume. More...
 

Additional Inherited Members

- 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)
 
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_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
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

template<class ELEMENT, template< class > class MESH>
class oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >

A class which allows the user to specify a prescribed volume (as opposed to a prescribed pressure) for in the region bounded by the membrane. Effectively adds an equation to the system for pressure. There would usually only be a single instance of this element in a problem.

Constructor & Destructor Documentation

◆ FoepplvonKarmanVolumeConstraintElement() [1/2]

template<class ELEMENT , template< class > class MESH>
oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::FoepplvonKarmanVolumeConstraintElement ( MESH< ELEMENT > *  bounding_mesh_pt,
const Vector< unsigned > &  contributing_region,
const double pressure = 0.0 
)
inline

Constructor. Takes pointer to mesh of Foeppl von Karman elements and a vector of unsigneds which identifies the regions within it that contribute to the controlled volume defined as int w dA (i.e. the "volume under the membrane"). The optional final argument specifies the initial value for the pressure that is "traded" in exchange for the volume constraint.

66  : Bounding_mesh_pt(bounding_mesh_pt),
67  Contributing_region(contributing_region),
69  {
70  // Create instance of pressure that is traded for volume constraint
71  Volume_control_pressure_pt = new Data(1);
73 
74  // Add the data object as internal data for this element
76  }
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
MESH< ELEMENT > * Bounding_mesh_pt
Definition: foeppl_von_karman_volume_constraint_element.h:262
Data * Volume_control_pressure_pt
Definition: foeppl_von_karman_volume_constraint_element.h:248
unsigned Pressure_data_index
Definition: foeppl_von_karman_volume_constraint_element.h:252
double * Prescribed_volume_pt
Pointer to target volume.
Definition: foeppl_von_karman_volume_constraint_element.h:269
Vector< unsigned > Contributing_region
Definition: foeppl_von_karman_volume_constraint_element.h:266
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:62

References oomph::GeneralisedElement::add_internal_data(), oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Pressure_data_index, oomph::Data::set_value(), and oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Volume_control_pressure_pt.

◆ ~FoepplvonKarmanVolumeConstraintElement()

template<class ELEMENT , template< class > class MESH>
oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::~FoepplvonKarmanVolumeConstraintElement ( )
inline
79 {}

◆ FoepplvonKarmanVolumeConstraintElement() [2/2]

template<class ELEMENT , template< class > class MESH>
oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::FoepplvonKarmanVolumeConstraintElement ( const FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH > &  )
delete

Broken copy constructor.

Member Function Documentation

◆ assign_additional_local_eqn_numbers()

template<class ELEMENT , template< class > class MESH>
void oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::assign_additional_local_eqn_numbers ( )
inlinevirtual

Assign the equation number for the new equation.

Reimplemented from oomph::GeneralisedElement.

114  {
116  }
int Volume_control_local_eqn
Local equation number of volume constraint.
Definition: foeppl_von_karman_volume_constraint_element.h:255
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Definition: elements.h:267

References oomph::GeneralisedElement::internal_local_eqn(), oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Pressure_data_index, and oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Volume_control_local_eqn.

◆ fill_in_contribution_to_jacobian()

template<class ELEMENT , template< class > class MESH>
void oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Fill in contribution to elemental residual and Jacobian.

Reimplemented from oomph::GeneralisedElement.

132  {
133  if (Volume_control_local_eqn >= 0)
134  {
135  // Only add contribution to residual; Jacobian (derivs w.r.t. to
136  // this element's external data is handled by FvK elements)
137  residuals[Volume_control_local_eqn] -= (*Prescribed_volume_pt);
138 
139 
140  /* // We are in charge... */
141  /* else */
142  /* { */
143  /* // NOTE: This is very slow but keep code alive -- can be */
144  /* // recycled if/when we ever make the Jacobians in the */
145  /* // Foeppl von Karman elements analytical. */
146  /* double t_start=TimingHelpers::timer(); */
147 
148  /* // Setup lookup scheme between local and global data */
149  /* std::map<unsigned,unsigned> local_external_eqn; */
150  /* unsigned next=nexternal_data(); */
151  /* for (unsigned j=0;j<next;j++) */
152  /* { */
153  /* Data* data_pt=external_data_pt(j); */
154  /* unsigned nval=data_pt->nvalue(); */
155  /* for (unsigned k=0;k<nval;k++) */
156  /* { */
157  /* int local_eqn=external_local_eqn(j,k); */
158  /* if (local_eqn>=0) */
159  /* { */
160  /* int global_eqn=data_pt->eqn_number(k); */
161  /* local_external_eqn[global_eqn]=local_eqn; */
162  /* } */
163  /* } */
164  /* } */
165 
166 
167  /* double t_end=TimingHelpers::timer(); */
168  /* oomph_info << "Time for local_external_eqn setup: " */
169  /* << t_end-t_start << std::endl; */
170  /* t_start=TimingHelpers::timer(); */
171 
172 
173  /* // Add initial contribution */
174  /* residuals[Volume_control_local_eqn] -= (*Prescribed_volume_pt); */
175 
176  /* // Initialise total bounded volume */
177  /* double bounded_volume = 0.0; */
178 
179  /* // Loop over the bounded regions */
180  /* unsigned n_contributing_regions = Contributing_region.size(); */
181  /* for(unsigned r = 0; r < n_contributing_regions; r++) */
182  /* { */
183  /* // Loop over the elements in the bounded regions */
184  /* unsigned n_inner_el = */
185  /* Bounding_mesh_pt->nregion_element(Contributing_region[r]); */
186  /* for(unsigned e = 0; e < n_inner_el; e++) */
187  /* { */
188  /* // Get element */
189  /* ELEMENT *el_pt = dynamic_cast<ELEMENT*> */
190  /* (Bounding_mesh_pt->region_element_pt(Contributing_region[r],e));
191  */
192 
193  /* // Get element's contribution to bounded volume and */
194  /* // derivative w.r.t. its unknowns */
195  /* double el_bounded_volume=0.0; */
196 
197  /* std::map<unsigned,double> d_bounded_volume_d_unknown; */
198  /* el_pt->fill_in_d_bounded_volume_d_unknown(el_bounded_volume,
199  */
200  /* d_bounded_volume_d_unknown);
201  */
202 
203  /* // Add contribution to bounded volume */
204  /* bounded_volume += el_bounded_volume; */
205 
206 
207  /* // Add contribution to Jacobian */
208  /* for (std::map<unsigned,double>::iterator it= */
209  /* d_bounded_volume_d_unknown.begin();it!= */
210  /* d_bounded_volume_d_unknown.end();it++) */
211  /* { */
212  /* unsigned global_eqn=(*it).first; */
213  /* unsigned local_eqn=local_external_eqn[global_eqn]; */
214  /* jacobian(Volume_control_local_eqn,local_eqn)+=(*it).second;
215  */
216  /* } */
217  /* } */
218  /* } */
219  /* // Add contribution to residuals */
220  /* residuals[Volume_control_local_eqn]+=bounded_volume; */
221 
222 
223  /* t_end=TimingHelpers::timer(); */
224  /* oomph_info << "Time for second part (actual work) of
225  * fill...jacobian: " */
226  /* << t_end-t_start << std::endl; */
227 
228  /* } */
229  }
230  }

References oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Volume_control_local_eqn.

◆ fill_in_contribution_to_residuals()

template<class ELEMENT , template< class > class MESH>
void oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Fill in residual: Difference between actual and prescribed bounded volume.

Reimplemented from oomph::GeneralisedElement.

121  {
122  if (Volume_control_local_eqn >= 0)
123  {
124  residuals[Volume_control_local_eqn] += -(*Prescribed_volume_pt);
125  }
126  }

References oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Volume_control_local_eqn.

◆ get_bounded_volume()

template<class ELEMENT , template< class > class MESH>
double oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::get_bounded_volume ( )
inline

Returns the volume "under the elements" in the constrained regions

91  {
92  double bounded_volume = 0.0;
93 
94  // Loop over the bounded regions
95  unsigned n_contributing_regions = Contributing_region.size();
96  for (unsigned r = 0; r < n_contributing_regions; r++)
97  {
98  // Loop over the elements in the bounded regions
99  unsigned n_inner_el =
100  Bounding_mesh_pt->nregion_element(Contributing_region[r]);
101  for (unsigned e = 0; e < n_inner_el; e++)
102  {
103  // Add the contribution
104  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(
105  Bounding_mesh_pt->region_element_pt(Contributing_region[r], e));
106  bounded_volume += el_pt->get_bounded_volume();
107  }
108  }
109  return bounded_volume;
110  }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
r
Definition: UniformPSDSelfTest.py:20

References oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Bounding_mesh_pt, oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Contributing_region, e(), and UniformPSDSelfTest::r.

◆ operator=()

template<class ELEMENT , template< class > class MESH>
void oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::operator= ( const FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH > &  )
delete

Broken assignment operator.

◆ pressure_data_pt()

template<class ELEMENT , template< class > class MESH>
Data* oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::pressure_data_pt ( ) const
inline

Access to Data object whose single value contains the pressure that has been "traded" for the volume constraint.

241  {
243  }

References oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Volume_control_pressure_pt.

◆ set_prescribed_volume()

template<class ELEMENT , template< class > class MESH>
void oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::set_prescribed_volume ( double volume_pt)
inline

Set pointer to target volume.

234  {
235  Prescribed_volume_pt = volume_pt;
236  }

References oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Prescribed_volume_pt.

Member Data Documentation

◆ Bounding_mesh_pt

template<class ELEMENT , template< class > class MESH>
MESH<ELEMENT>* oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Bounding_mesh_pt
protected

Pointer to mesh of Foeppl von Karman elements that bound the prescribed volume; NULL if the FvK elements that contribute to the bounding volume are in charge of adding their own contribution to the volume constraint equation (and the associated Jacobian)

Referenced by oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::get_bounded_volume().

◆ Contributing_region

template<class ELEMENT , template< class > class MESH>
Vector<unsigned> oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Contributing_region
protected

Region IDs in the bounding mesh that contribute to the prescribed/controlled volume

Referenced by oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::get_bounded_volume().

◆ Prescribed_volume_pt

template<class ELEMENT , template< class > class MESH>
double* oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Prescribed_volume_pt
protected

◆ Pressure_data_index

template<class ELEMENT , template< class > class MESH>
unsigned oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Pressure_data_index
protected

◆ Volume_control_local_eqn

◆ Volume_control_pressure_pt

template<class ELEMENT , template< class > class MESH>
Data* oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::Volume_control_pressure_pt
protected

Data object whose single value contains the pressure that has been "traded" for the volume constraint.

Referenced by oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::FoepplvonKarmanVolumeConstraintElement(), and oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::pressure_data_pt().


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