HomogenisationProblem< ELEMENT > Class Template Reference

Problem class to simulate viscous inclusion propagating along 2D channel. More...

+ Inheritance diagram for HomogenisationProblem< ELEMENT >:

Public Member Functions

 HomogenisationProblem ()
 Constructor. More...
 
 ~HomogenisationProblem ()
 Destructor. More...
 
void complete_problem_setup ()
 Finish the steup of the problem. More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve: More...
 
void calculate_coefficients ()
 
void sub_solve (const unsigned &n_dofs, DoubleVector &dx, DoubleVector &res, const bool &first_solve)
 
void solve ()
 Make our own solve function. More...
 
void doc_solution ()
 Doc the solution. More...
 
 HomogenisationProblem ()
 Constructor. More...
 
 ~HomogenisationProblem ()
 Destructor. More...
 
void complete_problem_setup ()
 Finish the steup of the problem. More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve: More...
 
void calculate_coefficients ()
 
void sub_solve (const unsigned &n_dofs, DoubleVector &dx, DoubleVector &res, const bool &first_solve)
 
void solve ()
 Make our own solve function. More...
 
void doc_solution ()
 
- Public Member Functions inherited from oomph::Problem
virtual void debug_hook_fct (const unsigned &i)
 
void set_analytic_dparameter (double *const &parameter_pt)
 
void unset_analytic_dparameter (double *const &parameter_pt)
 
bool is_dparameter_calculated_analytically (double *const &parameter_pt)
 
void set_analytic_hessian_products ()
 
void unset_analytic_hessian_products ()
 
bool are_hessian_products_calculated_analytically ()
 
void set_pinned_values_to_zero ()
 
bool distributed () const
 
virtual void actions_before_adapt ()
 
virtual void actions_after_adapt ()
 Actions that are to be performed after a mesh adaptation. More...
 
OomphCommunicatorcommunicator_pt ()
 access function to the oomph-lib communicator More...
 
const OomphCommunicatorcommunicator_pt () const
 access function to the oomph-lib communicator, const version More...
 
 Problem ()
 
 Problem (const Problem &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const Problem &)=delete
 Broken assignment operator. More...
 
virtual ~Problem ()
 Virtual destructor to clean up memory. More...
 
Mesh *& mesh_pt ()
 Return a pointer to the global mesh. More...
 
Mesh *const & mesh_pt () const
 Return a pointer to the global mesh (const version) More...
 
Mesh *& mesh_pt (const unsigned &imesh)
 
Mesh *const & mesh_pt (const unsigned &imesh) const
 Return a pointer to the i-th submesh (const version) More...
 
unsigned nsub_mesh () const
 Return number of submeshes. More...
 
unsigned add_sub_mesh (Mesh *const &mesh_pt)
 
void flush_sub_meshes ()
 
void build_global_mesh ()
 
void rebuild_global_mesh ()
 
LinearSolver *& linear_solver_pt ()
 Return a pointer to the linear solver object. More...
 
LinearSolver *const & linear_solver_pt () const
 Return a pointer to the linear solver object (const version) More...
 
LinearSolver *& mass_matrix_solver_for_explicit_timestepper_pt ()
 
LinearSolvermass_matrix_solver_for_explicit_timestepper_pt () const
 
EigenSolver *& eigen_solver_pt ()
 Return a pointer to the eigen solver object. More...
 
EigenSolver *const & eigen_solver_pt () const
 Return a pointer to the eigen solver object (const version) More...
 
Time *& time_pt ()
 Return a pointer to the global time object. More...
 
Timetime_pt () const
 Return a pointer to the global time object (const version). More...
 
doubletime ()
 Return the current value of continuous time. More...
 
double time () const
 Return the current value of continuous time (const version) More...
 
TimeStepper *& time_stepper_pt ()
 
const TimeSteppertime_stepper_pt () const
 
TimeStepper *& time_stepper_pt (const unsigned &i)
 Return a pointer to the i-th timestepper. More...
 
ExplicitTimeStepper *& explicit_time_stepper_pt ()
 Return a pointer to the explicit timestepper. More...
 
unsigned long set_timestepper_for_all_data (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data=false)
 
virtual void shift_time_values ()
 Shift all values along to prepare for next timestep. More...
 
AssemblyHandler *& assembly_handler_pt ()
 Return a pointer to the assembly handler object. More...
 
AssemblyHandler *const & assembly_handler_pt () const
 Return a pointer to the assembly handler object (const version) More...
 
doubleminimum_dt ()
 Access function to min timestep in adaptive timestepping. More...
 
doublemaximum_dt ()
 Access function to max timestep in adaptive timestepping. More...
 
unsignedmax_newton_iterations ()
 Access function to max Newton iterations before giving up. More...
 
void problem_is_nonlinear (const bool &prob_lin)
 Access function to Problem_is_nonlinear. More...
 
doublemax_residuals ()
 
booltime_adaptive_newton_crash_on_solve_fail ()
 Access function for Time_adaptive_newton_crash_on_solve_fail. More...
 
doublenewton_solver_tolerance ()
 
void add_time_stepper_pt (TimeStepper *const &time_stepper_pt)
 
void set_explicit_time_stepper_pt (ExplicitTimeStepper *const &explicit_time_stepper_pt)
 
void initialise_dt (const double &dt)
 
void initialise_dt (const Vector< double > &dt)
 
Data *& global_data_pt (const unsigned &i)
 Return a pointer to the the i-th global data object. More...
 
void add_global_data (Data *const &global_data_pt)
 
void flush_global_data ()
 
LinearAlgebraDistribution *const & dof_distribution_pt () const
 Return the pointer to the dof distribution (read-only) More...
 
unsigned long ndof () const
 Return the number of dofs. More...
 
unsigned ntime_stepper () const
 Return the number of time steppers. More...
 
unsigned nglobal_data () const
 Return the number of global data values. More...
 
unsigned self_test ()
 Self-test: Check meshes and global data. Return 0 for OK. More...
 
void enable_store_local_dof_pt_in_elements ()
 
void disable_store_local_dof_pt_in_elements ()
 
unsigned long assign_eqn_numbers (const bool &assign_local_eqn_numbers=true)
 
void describe_dofs (std::ostream &out= *(oomph_info.stream_pt())) const
 
void enable_discontinuous_formulation ()
 
void disable_discontinuous_formulation ()
 
void get_dofs (DoubleVector &dofs) const
 
void get_dofs (const unsigned &t, DoubleVector &dofs) const
 Return vector of the t'th history value of all dofs. More...
 
void set_dofs (const DoubleVector &dofs)
 Set the values of the dofs. More...
 
void set_dofs (const unsigned &t, DoubleVector &dofs)
 Set the history values of the dofs. More...
 
void set_dofs (const unsigned &t, Vector< double * > &dof_pt)
 
void add_to_dofs (const double &lambda, const DoubleVector &increment_dofs)
 Add lambda x incremenet_dofs[l] to the l-th dof. More...
 
doubleglobal_dof_pt (const unsigned &i)
 
doubledof (const unsigned &i)
 i-th dof in the problem More...
 
double dof (const unsigned &i) const
 i-th dof in the problem (const version) More...
 
double *& dof_pt (const unsigned &i)
 Pointer to i-th dof in the problem. More...
 
doubledof_pt (const unsigned &i) const
 Pointer to i-th dof in the problem (const version) More...
 
virtual void get_inverse_mass_matrix_times_residuals (DoubleVector &Mres)
 
virtual void get_dvaluesdt (DoubleVector &f)
 
virtual void get_residuals (DoubleVector &residuals)
 Get the total residuals Vector for the problem. More...
 
virtual void get_jacobian (DoubleVector &residuals, DenseDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, CRDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, CCDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, SumOfMatrices &jacobian)
 
void get_fd_jacobian (DoubleVector &residuals, DenseMatrix< double > &jacobian)
 Get the full Jacobian by finite differencing. More...
 
void get_derivative_wrt_global_parameter (double *const &parameter_pt, DoubleVector &result)
 
void get_hessian_vector_products (DoubleVectorWithHaloEntries const &Y, Vector< DoubleVectorWithHaloEntries > const &C, Vector< DoubleVectorWithHaloEntries > &product)
 
void solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &steady=true)
 Solve the eigenproblem. More...
 
void solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &steady=true)
 
virtual void get_eigenproblem_matrices (CRDoubleMatrix &mass_matrix, CRDoubleMatrix &main_matrix, const double &shift=0.0)
 
void assign_eigenvector_to_dofs (DoubleVector &eigenvector)
 Assign the eigenvector passed to the function to the dofs. More...
 
void add_eigenvector_to_dofs (const double &epsilon, const DoubleVector &eigenvector)
 
void store_current_dof_values ()
 Store the current values of the degrees of freedom. More...
 
void restore_dof_values ()
 Restore the stored values of the degrees of freedom. More...
 
void enable_jacobian_reuse ()
 
void disable_jacobian_reuse ()
 Disable recycling of Jacobian in Newton iteration. More...
 
bool jacobian_reuse_is_enabled ()
 Is recycling of Jacobian in Newton iteration enabled? More...
 
booluse_predictor_values_as_initial_guess ()
 
void newton_solve ()
 Use Newton method to solve the problem. More...
 
void enable_globally_convergent_newton_method ()
 enable globally convergent Newton method More...
 
void disable_globally_convergent_newton_method ()
 disable globally convergent Newton method More...
 
void newton_solve (unsigned const &max_adapt)
 
void steady_newton_solve (unsigned const &max_adapt=0)
 
void copy (Problem *orig_problem_pt)
 
virtual Problemmake_copy ()
 
virtual void read (std::ifstream &restart_file, bool &unsteady_restart)
 
virtual void read (std::ifstream &restart_file)
 
virtual void dump (std::ofstream &dump_file) const
 
void dump (const std::string &dump_file_name) const
 
void delete_all_external_storage ()
 
virtual void symmetrise_eigenfunction_for_adaptive_pitchfork_tracking ()
 
doublebifurcation_parameter_pt () const
 
void get_bifurcation_eigenfunction (Vector< DoubleVector > &eigenfunction)
 
void activate_fold_tracking (double *const &parameter_pt, const bool &block_solve=true)
 
void activate_bifurcation_tracking (double *const &parameter_pt, const DoubleVector &eigenvector, const bool &block_solve=true)
 
void activate_bifurcation_tracking (double *const &parameter_pt, const DoubleVector &eigenvector, const DoubleVector &normalisation, const bool &block_solve=true)
 
void activate_pitchfork_tracking (double *const &parameter_pt, const DoubleVector &symmetry_vector, const bool &block_solve=true)
 
void activate_hopf_tracking (double *const &parameter_pt, const bool &block_solve=true)
 
void activate_hopf_tracking (double *const &parameter_pt, const double &omega, const DoubleVector &null_real, const DoubleVector &null_imag, const bool &block_solve=true)
 
void deactivate_bifurcation_tracking ()
 
void reset_assembly_handler_to_default ()
 Reset the system to the standard non-augemented state. More...
 
double arc_length_step_solve (double *const &parameter_pt, const double &ds, const unsigned &max_adapt=0)
 
double arc_length_step_solve (Data *const &data_pt, const unsigned &data_index, const double &ds, const unsigned &max_adapt=0)
 
void reset_arc_length_parameters ()
 
intsign_of_jacobian ()
 
void explicit_timestep (const double &dt, const bool &shift_values=true)
 Take an explicit timestep of size dt. More...
 
void unsteady_newton_solve (const double &dt)
 
void unsteady_newton_solve (const double &dt, const bool &shift_values)
 
void unsteady_newton_solve (const double &dt, const unsigned &max_adapt, const bool &first, const bool &shift=true)
 
double doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const bool &first, const bool &shift=true)
 
double doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const unsigned &suppress_resolve_after_spatial_adapt_flag, const bool &first, const bool &shift=true)
 
double adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon)
 
double adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon, const bool &shift_values)
 
void assign_initial_values_impulsive ()
 
void assign_initial_values_impulsive (const double &dt)
 
void calculate_predictions ()
 Calculate predictions. More...
 
void enable_mass_matrix_reuse ()
 
void disable_mass_matrix_reuse ()
 
bool mass_matrix_reuse_is_enabled ()
 Return whether the mass matrix is being reused. More...
 
void refine_uniformly (const Vector< unsigned > &nrefine_for_mesh)
 
void refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh)
 
void refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void refine_uniformly (DocInfo &doc_info)
 
void refine_uniformly_and_prune (DocInfo &doc_info)
 
void refine_uniformly ()
 
void refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform refinement for submesh i_mesh with documentation. More...
 
void refine_uniformly (const unsigned &i_mesh)
 Do uniform refinement for submesh i_mesh without documentation. More...
 
void p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh)
 
void p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh)
 
void p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void p_refine_uniformly (DocInfo &doc_info)
 
void p_refine_uniformly_and_prune (DocInfo &doc_info)
 
void p_refine_uniformly ()
 
void p_refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform p-refinement for submesh i_mesh with documentation. More...
 
void p_refine_uniformly (const unsigned &i_mesh)
 Do uniform p-refinement for submesh i_mesh without documentation. More...
 
void refine_selected_elements (const Vector< unsigned > &elements_to_be_refined)
 
void refine_selected_elements (const Vector< RefineableElement * > &elements_to_be_refined_pt)
 
void refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined)
 
void refine_selected_elements (const unsigned &i_mesh, const Vector< RefineableElement * > &elements_to_be_refined_pt)
 
void refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined)
 
void refine_selected_elements (const Vector< Vector< RefineableElement * >> &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const Vector< unsigned > &elements_to_be_refined)
 
void p_refine_selected_elements (const Vector< PRefineableElement * > &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined)
 
void p_refine_selected_elements (const unsigned &i_mesh, const Vector< PRefineableElement * > &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined)
 
void p_refine_selected_elements (const Vector< Vector< PRefineableElement * >> &elements_to_be_refined_pt)
 
unsigned unrefine_uniformly ()
 
unsigned unrefine_uniformly (const unsigned &i_mesh)
 
void p_unrefine_uniformly (DocInfo &doc_info)
 
void p_unrefine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform p-unrefinement for submesh i_mesh without documentation. More...
 
void adapt (unsigned &n_refined, unsigned &n_unrefined)
 
void adapt ()
 
void p_adapt (unsigned &n_refined, unsigned &n_unrefined)
 
void p_adapt ()
 
void adapt_based_on_error_estimates (unsigned &n_refined, unsigned &n_unrefined, Vector< Vector< double >> &elemental_error)
 
void adapt_based_on_error_estimates (Vector< Vector< double >> &elemental_error)
 
void get_all_error_estimates (Vector< Vector< double >> &elemental_error)
 
void doc_errors (DocInfo &doc_info)
 Get max and min error for all elements in submeshes. More...
 
void doc_errors ()
 Get max and min error for all elements in submeshes. More...
 
void enable_info_in_newton_solve ()
 
void disable_info_in_newton_solve ()
 Disable the output of information when in the newton solver. More...
 
- Public Member Functions inherited from oomph::ExplicitTimeSteppableObject
 ExplicitTimeSteppableObject ()
 Empty constructor. More...
 
 ExplicitTimeSteppableObject (const ExplicitTimeSteppableObject &)=delete
 Broken copy constructor. More...
 
void operator= (const ExplicitTimeSteppableObject &)=delete
 Broken assignment operator. More...
 
virtual ~ExplicitTimeSteppableObject ()
 Empty destructor. More...
 
virtual void actions_before_explicit_stage ()
 
virtual void actions_after_explicit_stage ()
 

Private Attributes

unsigned M
 
unsigned P
 
Vector< Vector< DenseMatrix< double > > > C_eff
 
Vector< GeomObject * > Internal_circle_pt
 
TriangleMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to Bulk_mesh. More...
 
Vector< TriangleMeshPolygon * > Inclusion_polygon_pt
 Vector storing pointer to the inclusion polygons. More...
 
TriangleMeshPolygonOuter_boundary_polyline_pt
 Triangle mesh polygon for outer boundary. More...
 

Additional Inherited Members

- Public Types inherited from oomph::Problem
typedef void(* SpatialErrorEstimatorFctPt) (Mesh *&mesh_pt, Vector< double > &elemental_error)
 Function pointer for spatial error estimator. More...
 
typedef void(* SpatialErrorEstimatorWithDocFctPt) (Mesh *&mesh_pt, Vector< double > &elemental_error, DocInfo &doc_info)
 Function pointer for spatial error estimator with doc. More...
 
- Public Attributes inherited from oomph::Problem
bool Shut_up_in_newton_solve
 
- Static Public Attributes inherited from oomph::Problem
static bool Suppress_warning_about_actions_before_read_unstructured_meshes
 
- Protected Types inherited from oomph::Problem
enum  Assembly_method {
  Perform_assembly_using_vectors_of_pairs , Perform_assembly_using_two_vectors , Perform_assembly_using_maps , Perform_assembly_using_lists ,
  Perform_assembly_using_two_arrays
}
 Enumerated flags to determine which sparse assembly method is used. More...
 
- Protected Member Functions inherited from oomph::Problem
unsigned setup_element_count_per_dof ()
 
virtual void sparse_assemble_row_or_column_compressed (Vector< int * > &column_or_row_index, Vector< int * > &row_or_column_start, Vector< double * > &value, Vector< unsigned > &nnz, Vector< double * > &residual, bool compressed_row_flag)
 
virtual void actions_before_newton_convergence_check ()
 
virtual void actions_before_newton_step ()
 
virtual void actions_after_newton_step ()
 
virtual void actions_before_implicit_timestep ()
 
virtual void actions_after_implicit_timestep ()
 
virtual void actions_after_implicit_timestep_and_error_estimation ()
 
virtual void actions_before_explicit_timestep ()
 Actions that should be performed before each explicit time step. More...
 
virtual void actions_after_explicit_timestep ()
 Actions that should be performed after each explicit time step. More...
 
virtual void actions_before_read_unstructured_meshes ()
 
virtual void actions_after_read_unstructured_meshes ()
 
virtual void actions_after_change_in_global_parameter (double *const &parameter_pt)
 
virtual void actions_after_change_in_bifurcation_parameter ()
 
virtual void actions_after_parameter_increase (double *const &parameter_pt)
 
doubledof_derivative (const unsigned &i)
 
doubledof_current (const unsigned &i)
 
virtual void set_initial_condition ()
 
virtual double global_temporal_error_norm ()
 
unsigned newton_solve_continuation (double *const &parameter_pt)
 
unsigned newton_solve_continuation (double *const &parameter_pt, DoubleVector &z)
 
void calculate_continuation_derivatives (double *const &parameter_pt)
 
void calculate_continuation_derivatives (const DoubleVector &z)
 
void calculate_continuation_derivatives_fd (double *const &parameter_pt)
 
bool does_pointer_correspond_to_problem_data (double *const &parameter_pt)
 
void set_consistent_pinned_values_for_continuation ()
 
- Protected Attributes inherited from oomph::Problem
Vector< Problem * > Copy_of_problem_pt
 
std::map< double *, boolCalculate_dparameter_analytic
 
bool Calculate_hessian_products_analytic
 
LinearAlgebraDistributionDof_distribution_pt
 
Vector< double * > Dof_pt
 Vector of pointers to dofs. More...
 
DoubleVectorWithHaloEntries Element_count_per_dof
 
double Relaxation_factor
 
double Newton_solver_tolerance
 
unsigned Max_newton_iterations
 Maximum number of Newton iterations. More...
 
unsigned Nnewton_iter_taken
 
Vector< doubleMax_res
 Maximum residuals at start and after each newton iteration. More...
 
double Max_residuals
 
bool Time_adaptive_newton_crash_on_solve_fail
 
bool Jacobian_reuse_is_enabled
 Is re-use of Jacobian in Newton iteration enabled? Default: false. More...
 
bool Jacobian_has_been_computed
 
bool Problem_is_nonlinear
 
bool Pause_at_end_of_sparse_assembly
 
bool Doc_time_in_distribute
 
unsigned Sparse_assembly_method
 
unsigned Sparse_assemble_with_arrays_initial_allocation
 
unsigned Sparse_assemble_with_arrays_allocation_increment
 
Vector< Vector< unsigned > > Sparse_assemble_with_arrays_previous_allocation
 
double Numerical_zero_for_sparse_assembly
 
double FD_step_used_in_get_hessian_vector_products
 
bool Mass_matrix_reuse_is_enabled
 
bool Mass_matrix_has_been_computed
 
bool Discontinuous_element_formulation
 
double Minimum_dt
 Minimum desired dt: if dt falls below this value, exit. More...
 
double Maximum_dt
 Maximum desired dt. More...
 
double DTSF_max_increase
 
double DTSF_min_decrease
 
double Minimum_dt_but_still_proceed
 
bool Scale_arc_length
 Boolean to control whether arc-length should be scaled. More...
 
double Desired_proportion_of_arc_length
 Proportion of the arc-length to taken by the parameter. More...
 
double Theta_squared
 
int Sign_of_jacobian
 Storage for the sign of the global Jacobian. More...
 
double Continuation_direction
 
double Parameter_derivative
 Storage for the derivative of the global parameter wrt arc-length. More...
 
double Parameter_current
 Storage for the present value of the global parameter. More...
 
bool Use_continuation_timestepper
 Boolean to control original or new storage of dof stuff. More...
 
unsigned Dof_derivative_offset
 
unsigned Dof_current_offset
 
Vector< doubleDof_derivative
 Storage for the derivative of the problem variables wrt arc-length. More...
 
Vector< doubleDof_current
 Storage for the present values of the variables. More...
 
double Ds_current
 Storage for the current step value. More...
 
unsigned Desired_newton_iterations_ds
 
double Minimum_ds
 Minimum desired value of arc-length. More...
 
bool Bifurcation_detection
 Boolean to control bifurcation detection via determinant of Jacobian. More...
 
bool Bisect_to_find_bifurcation
 Boolean to control wheter bisection is used to located bifurcation. More...
 
bool First_jacobian_sign_change
 Boolean to indicate whether a sign change has occured in the Jacobian. More...
 
bool Arc_length_step_taken
 Boolean to indicate whether an arc-length step has been taken. More...
 
bool Use_finite_differences_for_continuation_derivatives
 
OomphCommunicatorCommunicator_pt
 The communicator for this problem. More...
 
bool Always_take_one_newton_step
 
double Timestep_reduction_factor_after_nonconvergence
 
bool Keep_temporal_error_below_tolerance
 
- Static Protected Attributes inherited from oomph::Problem
static ContinuationStorageScheme Continuation_time_stepper
 Storage for the single static continuation timestorage object. More...
 

Detailed Description

template<class ELEMENT>
class HomogenisationProblem< ELEMENT >

Problem class to simulate viscous inclusion propagating along 2D channel.

Constructor & Destructor Documentation

◆ HomogenisationProblem() [1/2]

template<class ELEMENT >
HomogenisationProblem< ELEMENT >::HomogenisationProblem

Constructor.

324  : M(0), P(0)
325 {
326  // Allocate the timestepper(a Steady default)
327  this->add_time_stepper_pt(new Steady<0>);
328 
329  //Set the number of elements on each side
330  unsigned n_element_on_side = 40;
331 
332  //Each side has uniform length so calculate the spacing between vertices
333  double vertex_spacing = 1.0 / (double)(n_element_on_side);
334 
335  //Set the volume fraction of the (single) fibre
336  double volume_fraction = 0.05;
337 
338  // Build the boundary segments for outer boundary, consisting of
339  //--------------------------------------------------------------
340  // four separate polylines
341  //------------------------
342  Vector<TriangleMeshCurveSection*> boundary_polyline_pt(4);
343 
344  //Each polyline has n_element_on_side + 1 vertices
345  unsigned n_vertex= n_element_on_side + 1;
346  Vector<Vector<double> > vertex_coord(n_vertex);
347  for(unsigned i=0;i<n_vertex;i++)
348  {
349  vertex_coord[i].resize(2);
350  }
351 
352  // First polyline: Left-hand edge, set the nodal positions
353  for(unsigned i=0;i<n_vertex;i++)
354  {
355  vertex_coord[i][0]=0.0;
356  vertex_coord[i][1] = i*vertex_spacing;
357  }
358 
359  // Build the 1st boundary polyline
360  boundary_polyline_pt[0] = new TriangleMeshPolyLine(vertex_coord,3);
361 
362  // Second boundary polyline: Top edge
363  for(unsigned i=0;i<n_vertex;i++)
364  {
365  vertex_coord[i][0]= i*vertex_spacing;
366  vertex_coord[i][1]=1.0;
367  }
368 
369  // Build the 2nd boundary polyline
370  boundary_polyline_pt[1] = new TriangleMeshPolyLine(vertex_coord,2);
371 
372  // Third boundary polyline: Right-hand edge
373  for(unsigned i=0;i<n_vertex;i++)
374  {
375  vertex_coord[i][0]=1.0;
376  vertex_coord[i][1]=1.0 - i*vertex_spacing;
377  }
378 
379  // Build the 3rd boundary polyline
380  boundary_polyline_pt[2] = new TriangleMeshPolyLine(vertex_coord,1);
381 
382  // Fourth boundary polyline: Bottom edge
383  for(unsigned i=0;i<n_vertex;i++)
384  {
385  vertex_coord[i][0] = 1.0 - i*vertex_spacing;
386  vertex_coord[i][1] = 0.0;
387  }
388 
389  // Build the 4th boundary polyline
390  boundary_polyline_pt[3] = new TriangleMeshPolyLine(vertex_coord,0);
391 
392  // Create the triangle mesh polygon for outer boundary
393  Outer_boundary_polyline_pt = new TriangleMeshPolygon(boundary_polyline_pt);
394 
395  //Storage for the outer boundary
396  TriangleMeshClosedCurve* outer_boundary_closed_curve_pt =
398 
399  // Now define initial shape of inclusion with a geometric object
400  //---------------------------------------------------------------
401  // Build an inclusion
402  Internal_circle_pt.resize(1);
403  double rad = std::sqrt(volume_fraction/(4.0*atan(1.0)));
404  double Radius[2] = {rad,rad}; //{0.5*rad,rad};
405 
406  Vector<TriangleMeshClosedCurve*> curvilinear_inclusion_pt(1);
407 
408  for(unsigned h=0;h<1;h++)
409  {
410  Internal_circle_pt[h] =
411  new Circle(0.5,0.5,Radius[h],this->time_stepper_pt());
412 
413  // Build the two parts of the curvilinear boundary
414  //Note that there could well be a memory leak here owing to some stupid
415  //choices
416  Vector<TriangleMeshCurveSection*> curvilinear_boundary_pt(2);
417 
418  // First part of curvilinear boundary
419  //-----------------------------------
420  double zeta_start=0.0;
421  double zeta_end=MathematicalConstants::Pi;
422  unsigned nsegment=20;
423  unsigned boundary_id=4+2*h;
424  curvilinear_boundary_pt[0]=new TriangleMeshCurviLine(
425  Internal_circle_pt[h],zeta_start,zeta_end,
426  nsegment,boundary_id);
427 
428  // Second part of curvilinear boundary
429  //-------------------------------------
430  zeta_start=MathematicalConstants::Pi;
431  zeta_end=2.0*MathematicalConstants::Pi;
432  nsegment=20;
433  boundary_id=5 + 2*h;
434  curvilinear_boundary_pt[1]=new TriangleMeshCurviLine(
435  Internal_circle_pt[h],zeta_start,zeta_end,
436  nsegment,boundary_id);
437 
438  curvilinear_inclusion_pt[h]=
439  new TriangleMeshClosedCurve(curvilinear_boundary_pt);
440 
441  //Delete the curvilines not possible because there are some memory
442  //issues to consider
443  //delete curvilinear_boundary_pt[1];
444  //delete curvilinear_boundary_pt[0];
445  }
446 
447  // Use the TriangleMeshParameters object for gathering all
448  // the necessary arguments for the TriangleMesh object
449  TriangleMeshParameters triangle_mesh_parameters(
450  outer_boundary_closed_curve_pt);
451 
452  // Define the holes on the boundary
453  triangle_mesh_parameters.internal_closed_curve_pt() =
454  curvilinear_inclusion_pt;
455 
456  double uniform_element_area = 0.001;
457 
458  // Define the maximum element area
459  triangle_mesh_parameters.element_area() =
460  uniform_element_area;
461 
462  //Coordinates
463  Vector<double> region1(2,0.5);
464 
465  // Define the region (the internal fibre)
466  triangle_mesh_parameters.add_region_coordinates(1, region1);
467 
468  //Vector<double> region2(2);
469  //region2[0] = 0.5;
470  //region2[1] = 0.5 + 0.5*(Radius[1] + Radius[0]);
471 
472  //Set the second region (the annulus)
473  //triangle_mesh_parameters.add_region_coordinates(2, region2);
474 
475  //Prevent refinement on the boundaries
476  triangle_mesh_parameters.disable_boundary_refinement();
477 
478  // Create the mesh
479  Bulk_mesh_pt =
481  triangle_mesh_parameters, this->time_stepper_pt());
482 
483  //Output the regions
484  unsigned n_region = this->Bulk_mesh_pt->nregion();
485  for(unsigned i=0;i<n_region;i++)
486  {
487  unsigned n_element = Bulk_mesh_pt->nregion_element(i);
488 
489  std::stringstream output_string;
490  output_string << "Region" << i << ".dat";
491  std::ofstream output_file(output_string.str().c_str());
492 
493  for(unsigned e=0;e<n_element;e++)
494  {
495  Bulk_mesh_pt->region_element_pt(i,e)->output(output_file,5);
496  }
497  output_file.close();
498  }
499 
500  // Output boundary and mesh initial mesh for information
501  this->Bulk_mesh_pt->output_boundaries("boundaries.dat");
502  this->Bulk_mesh_pt->output("mesh.dat");
503 
504  // Set boundary condition and complete the build of all elements
505  this->complete_problem_setup();
506 
507  // Add meshes to the problem
508  //----------------------------
509 
510  // Add Bulk_mesh_pt sub meshes
511  this->add_sub_mesh(Bulk_mesh_pt);
512 
513  // Build global mesh
514  this->build_global_mesh();
515 
516  // Setup equation numbering scheme
517  cout <<"Number of equations: " << this->assign_eqn_numbers() << std::endl;
518 
519  //Build storage for the effective modulus
520  C_eff.resize(3);
521  for(unsigned i=0;i<3;i++)
522  {
523  C_eff[i].resize(3);
524  for(unsigned j=0;j<3;j++)
525  {
526  C_eff[i][j].resize(3,3,0.0);
527  }
528  }
529 
530  //Solve the problem
531  solve();
532 
533  //The final effective tensor is actually C_pmij because of the way I have
534  //set up the storage. It shouldn't matter because of the symmetries of the
535  //stress tensor, however.
536 
537  std::cout << "Output effective coefficients\n";
538 
539  std::cout << "C11 " << C_eff[0][0](0,0) << "\n";
540  std::cout << "C12 " << C_eff[0][0](1,1) << "\n";
541  std::cout << "C13 " << C_eff[0][0](2,2) << "\n";
542  std::cout << "C16 " << C_eff[0][0](0,1) << "\n";
543  std::cout << "C22 " << C_eff[1][1](1,1) << "\n";
544  std::cout << "C23 " << C_eff[1][1](2,2) << "\n";
545  std::cout << "C33 " << C_eff[2][2](2,2) << "\n";
546  std::cout << "C36 " << C_eff[2][2](0,1) << "\n";
547  std::cout << "C44 " << C_eff[1][2](1,2) << "\n";
548  std::cout << "C45 " << C_eff[1][2](0,2) << "\n";
549  std::cout << "C55 " << C_eff[0][2](0,2) << "\n";
550  std::cout << "C66 " << C_eff[0][1](0,1) << "\n";
551 
552  //Output to a file
553  std::ofstream output("C_eff.dat");
554 
555  output << "C11 " << C_eff[0][0](0,0) << "\n";
556  output << "C12 " << C_eff[0][0](1,1) << "\n";
557  output << "C13 " << C_eff[0][0](2,2) << "\n";
558  output << "C16 " << C_eff[0][0](0,1) << "\n";
559  output << "C22 " << C_eff[1][1](1,1) << "\n";
560  output << "C23 " << C_eff[1][1](2,2) << "\n";
561  output << "C33 " << C_eff[2][2](2,2) << "\n";
562  output << "C36 " << C_eff[2][2](0,1) << "\n";
563  output << "C44 " << C_eff[1][2](1,2) << "\n";
564  output << "C45 " << C_eff[1][2](0,2) << "\n";
565  output << "C55 " << C_eff[0][2](0,2) << "\n";
566  output << "C66 " << C_eff[0][1](0,1) << std::endl;
567 
568  //Close the output file
569  output.close();
570 
571 } // end_of_constructor
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void solve()
Make our own solve function.
Definition: two_dim.cc:265
void complete_problem_setup()
Finish the steup of the problem.
Definition: two_dim.cc:578
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to Bulk_mesh.
Definition: two_dim.cc:309
Vector< GeomObject * > Internal_circle_pt
Definition: two_dim.cc:147
TriangleMeshPolygon * Outer_boundary_polyline_pt
Triangle mesh polygon for outer boundary.
Definition: two_dim.cc:315
Vector< Vector< DenseMatrix< double > > > C_eff
Definition: two_dim.cc:144
unsigned P
Definition: two_dim.cc:141
unsigned M
Definition: two_dim.cc:141
Definition: geom_objects.h:873
void add_time_stepper_pt(TimeStepper *const &time_stepper_pt)
Definition: problem.cc:1545
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
TimeStepper *& time_stepper_pt()
Definition: problem.h:1524
Base class defining a closed curve for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1339
Definition: unstructured_two_d_mesh_geometry_base.h:662
Definition: triangle_mesh.template.h:94
Class defining a polyline for use in Triangle Mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:868
Class defining a closed polygon for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1451
Definition: triangle_mesh.template.h:424
Definition: oomph-lib/src/generic/Vector.h:58
double Pi
Definition: two_d_biharmonic.cc:235
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 atan(const bfloat16 &a)
Definition: BFloat16.h:636
double Radius
Radius of the pipe.
Definition: pipe.cc:55
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::TriangleMeshParameters::add_region_coordinates(), oomph::Problem::add_sub_mesh(), oomph::Problem::add_time_stepper_pt(), oomph::Problem::assign_eqn_numbers(), Eigen::bfloat16_impl::atan(), oomph::Problem::build_global_mesh(), HomogenisationProblem< ELEMENT >::Bulk_mesh_pt, HomogenisationProblem< ELEMENT >::C_eff, HomogenisationProblem< ELEMENT >::complete_problem_setup(), oomph::TriangleMeshParameters::disable_boundary_refinement(), e(), oomph::TriangleMeshParameters::element_area(), i, HomogenisationProblem< ELEMENT >::Internal_circle_pt, oomph::TriangleMeshParameters::internal_closed_curve_pt(), j, HomogenisationProblem< ELEMENT >::Outer_boundary_polyline_pt, output(), BiharmonicTestFunctions2::Pi, Global_Physical_Variables::Radius, HomogenisationProblem< ELEMENT >::solve(), sqrt(), and oomph::Problem::time_stepper_pt().

◆ ~HomogenisationProblem() [1/2]

template<class ELEMENT >
HomogenisationProblem< ELEMENT >::~HomogenisationProblem ( )
inline

Destructor.

156  {
157  // Kill data associated with outer boundary
158  unsigned n= this->Outer_boundary_polyline_pt->npolyline();
159  for (unsigned j=0;j<n;j++)
160  {
162  }
164 
165  //Kill data associated with inclusions
166  unsigned n_inclusion = Inclusion_polygon_pt.size();
167  for(unsigned iinclusion=0;iinclusion<n_inclusion;iinclusion++)
168  {
169  unsigned n=Inclusion_polygon_pt[iinclusion]->npolyline();
170  for (unsigned j=0;j<n;j++)
171  {
172  delete Inclusion_polygon_pt[iinclusion]->polyline_pt(j);
173  }
174  delete Inclusion_polygon_pt[iinclusion];
175  }
176  // Delete solid mesh
177  delete Bulk_mesh_pt;
178  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Vector< TriangleMeshPolygon * > Inclusion_polygon_pt
Vector storing pointer to the inclusion polygons.
Definition: two_dim.cc:312
unsigned npolyline() const
Number of constituent polylines.
Definition: unstructured_two_d_mesh_geometry_base.h:1482
TriangleMeshPolyLine * polyline_pt(const unsigned &i) const
Pointer to i-th constituent polyline.
Definition: unstructured_two_d_mesh_geometry_base.h:1488

References j, and n.

◆ HomogenisationProblem() [2/2]

template<class ELEMENT >
HomogenisationProblem< ELEMENT >::HomogenisationProblem ( )

Constructor.

◆ ~HomogenisationProblem() [2/2]

template<class ELEMENT >
HomogenisationProblem< ELEMENT >::~HomogenisationProblem ( )
inline

Destructor.

155  {
156  // Kill data associated with outer boundary
157  unsigned n= this->Outer_boundary_polyline_pt->npolyline();
158  for (unsigned j=0;j<n;j++)
159  {
161  }
163 
164  //Kill data associated with inclusions
165  unsigned n_inclusion = Inclusion_polygon_pt.size();
166  for(unsigned iinclusion=0;iinclusion<n_inclusion;iinclusion++)
167  {
168  unsigned n=Inclusion_polygon_pt[iinclusion]->npolyline();
169  for (unsigned j=0;j<n;j++)
170  {
171  delete Inclusion_polygon_pt[iinclusion]->polyline_pt(j);
172  }
173  delete Inclusion_polygon_pt[iinclusion];
174  }
175  // Delete solid mesh
176  delete Bulk_mesh_pt;
177  }

References j, and n.

Member Function Documentation

◆ actions_after_newton_solve() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

184 {}

◆ actions_after_newton_solve() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

183 {}

◆ actions_before_newton_solve() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve:

Reimplemented from oomph::Problem.

187 {}

◆ actions_before_newton_solve() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve:

Reimplemented from oomph::Problem.

186 {}

◆ calculate_coefficients() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::calculate_coefficients ( )
inline

Calculate the values of the effective modulus by integrating over each element

192  {
193  //Loop over all elements
194  unsigned n_element = this->Bulk_mesh_pt->nelement();
195  for(unsigned e=0;e<n_element;e++)
196  {
197  ELEMENT* el_pt =
198  dynamic_cast<ELEMENT*>(this->Bulk_mesh_pt->element_pt(e));
199 
200  //Add the contribution to the effective modulus
201  el_pt->calculate_effective_modulus(C_eff[P][M]);
202  }
203  }
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186

References e(), and Global_Physical_Variables::P.

◆ calculate_coefficients() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::calculate_coefficients ( )
inline

Calculate the values of the effective modulus by integrating over each element

191  {
192  //Loop over all elements
193  unsigned n_element = this->Bulk_mesh_pt->nelement();
194  for(unsigned e=0;e<n_element;e++)
195  {
196  ELEMENT* el_pt =
197  dynamic_cast<ELEMENT*>(this->Bulk_mesh_pt->element_pt(e));
198 
199  //Add the contribution to the effective modulus
200  el_pt->calculate_effective_modulus(C_eff[P][M]);
201  }
202  }

References e(), and Global_Physical_Variables::P.

◆ complete_problem_setup() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::complete_problem_setup

Finish the steup of the problem.

Set boundary conditions and complete the build of all elements.

579  {
580  //Collect vectors of the nodes on the mesh boundaries
581  Vector<Vector<Node*> > boundary_nodes_pt(4);
582 
583  for(unsigned b=0;b<4;b++)
584  {
585  unsigned n_node = this->Bulk_mesh_pt->nboundary_node(b);
586  boundary_nodes_pt[b].resize(n_node);
587  for(unsigned n=0;n<n_node;n++)
588  {
589  boundary_nodes_pt[b][n] = this->Bulk_mesh_pt->boundary_node_pt(b,n);
590  }
591  }
592 
593  //Now let's sort each one
594  std::sort(boundary_nodes_pt[0].begin(), boundary_nodes_pt[0].end(),
596  std::sort(boundary_nodes_pt[1].begin(), boundary_nodes_pt[1].end(),
598  std::sort(boundary_nodes_pt[2].begin(), boundary_nodes_pt[2].end(),
600  std::sort(boundary_nodes_pt[3].begin(), boundary_nodes_pt[3].end(),
602 
603 
604  //Now we can make it periodic
605  double tol = 1.0e-7;
606  unsigned n_node = boundary_nodes_pt[0].size();
607  for(unsigned n=1;n<n_node-1;n++)
608  {
609  //If we have different x coordinates complain
610  if(std::abs(boundary_nodes_pt[0][n]->x(0)
611  - boundary_nodes_pt[2][n]->x(0)) > tol)
612  {
613  std::ostringstream error_stream;
614  error_stream <<
615  "Trying to make periodic nodes across the top/bottom boundary,\n"
616  << "but the nodes have x-coordinates "
617  << boundary_nodes_pt[0][n]->x(0) << " and "
618  << boundary_nodes_pt[2][n]->x(0);
619 
620  throw OomphLibError(error_stream.str(),
623  }
624 
625  //Make it periodic
626  boundary_nodes_pt[0][n]->make_periodic(boundary_nodes_pt[2][n]);
627  }
628 
629  n_node = boundary_nodes_pt[1].size();
630  for(unsigned n=1;n<n_node-1;n++)
631  {
632  //If we have different y coordinates complain
633  if(std::abs(boundary_nodes_pt[1][n]->x(1) - boundary_nodes_pt[3][n]->x(1)) > tol)
634  {
635  std::ostringstream error_stream;
636  error_stream <<
637  "Trying to make periodic nodes across the left/right boundary,\n"
638  << "but the nodes have y-coordinates "
639  << boundary_nodes_pt[1][n]->x(1) << " and "
640  << boundary_nodes_pt[3][n]->x(1);
641 
642  throw OomphLibError(error_stream.str(),
645  }
646 
647  //Make it periodic
648  boundary_nodes_pt[1][n]->make_periodic(boundary_nodes_pt[3][n]);
649  }
650 
651 
652  //Vector of corner nodes
653  Vector<Node*> corner_nodes_pt(4);
654  corner_nodes_pt[0] = boundary_nodes_pt[1][0];
655  corner_nodes_pt[1] = boundary_nodes_pt[1][n_node-1];
656  corner_nodes_pt[2] = boundary_nodes_pt[3][0];
657  corner_nodes_pt[3] = boundary_nodes_pt[3][n_node-1];
658 
659  //PUT IN A TEST TO CHECK THAT THE CORNERS ARE REALLY THE CORNERS!
660 
661  //Make these nodes periodic from the bottom node
662  corner_nodes_pt[0]->make_periodic_nodes(corner_nodes_pt);
663 
664  //Pin the corner to suppress rigid body motions
665  corner_nodes_pt[0]->pin(0);
666  corner_nodes_pt[0]->pin(1);
667  corner_nodes_pt[0]->pin(2);
668 
669  // Complete the build of all elements so they are fully functional
670  // Remember that adaptation for triangle meshes involves a complete
671  // regneration of the mesh (rather than splitting as in tree-based
672  // meshes where such parameters can be passed down from the father
673  // element!)
674  unsigned n_element = this->Bulk_mesh_pt->nregion_element(0);
675  for(unsigned e=0;e<n_element;e++)
676  {
677  // Upcast from GeneralisedElement to the present element
678  ELEMENT* el_pt =
679  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->region_element_pt(0,e));
680 
681  // Set pointer to the elasticity tensor function
682  el_pt->elasticity_tensor_fct_pt() =
684 
685  // Set the values of m and p
686  el_pt->m_pt() = &(this->M);
687  el_pt->p_pt() = &(this->P);
688  }
689 
690  //For the elements within the inclusion (regions 1 and 2
691  //set the viscosity ratio
692  for(unsigned r=1;r<3;r++)
693  {
694  n_element = Bulk_mesh_pt->nregion_element(r);
695  for(unsigned e=0;e<n_element;e++)
696  {
697  // Upcast from GeneralisedElement to the present element
698  ELEMENT* el_pt =
699  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->region_element_pt(r,e));
700 
701  el_pt->elasticity_tensor_fct_pt() =
703 
704  // Set the values of m and p
705  el_pt->m_pt() = &(this->M);
706  el_pt->p_pt() = &(this->P);
707  }
708 
709  }
710 
711  }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
Scalar * b
Definition: benchVecAdd.cpp:17
Definition: two_dim.cc:52
Definition: two_dim.cc:73
Definition: oomph_definitions.h:222
static constexpr lastp1_t end
Definition: IndexedViewHelper.h:79
void bulk_elasticity_tensor_pt(const Vector< double > &x, ElasticityTensor *&E_pt)
Definition: two_dim.cc:116
void fibre_elasticity_tensor_pt(const Vector< double > &x, ElasticityTensor *&E_pt)
Definition: two_dim.cc:124
r
Definition: UniformPSDSelfTest.py:20
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References abs(), b, Problem_Parameter::bulk_elasticity_tensor_pt(), e(), Eigen::placeholders::end, Problem_Parameter::fibre_elasticity_tensor_pt(), n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, Global_Physical_Variables::P, UniformPSDSelfTest::r, and plotDoE::x.

Referenced by HomogenisationProblem< ELEMENT >::HomogenisationProblem().

◆ complete_problem_setup() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::complete_problem_setup ( )

Finish the steup of the problem.

◆ doc_solution() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::doc_solution

Doc the solution.

720 {
721  this->Bulk_mesh_pt->output("soln.dat",5);
722 }

◆ doc_solution() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::doc_solution ( )

◆ solve() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::solve ( )
inline

Make our own solve function.

266  {
267  //Enable the resolve
269 
270  //Find total number of dofs
271  unsigned long n_dofs = ndof();
272 
273  //Set up the Vector to hold the solution
274  //Vector<double> dx(n_dofs,0.0);
275  DoubleVector dx;
276  //Set up a vector to hold the residuals
277  //Vector<double> res(n_dofs);
279 
280  //Only need to loop over the upper "half" of the matrix
281  for(P=0;P<3;P++)
282  {
283  for(M=P;M<3;M++)
284  {
285  if(M==0)
286  {
287  //Do the solve once only
288  sub_solve(n_dofs,dx,res,true);
289  }
290  //Otherwise it's resolves
291  else
292  {
293  sub_solve(n_dofs,dx,res,false);
294  }
295  }
296  }
297 
298  //Disable the resolve
300  }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
void sub_solve(const unsigned &n_dofs, DoubleVector &dx, DoubleVector &res, const bool &first_solve)
Definition: two_dim.cc:207
Definition: double_vector.h:58
virtual void enable_resolve()
Definition: linear_solver.h:135
virtual void disable_resolve()
Definition: linear_solver.h:144
unsigned long ndof() const
Return the number of dofs.
Definition: problem.h:1674
LinearSolver *& linear_solver_pt()
Return a pointer to the linear solver object.
Definition: problem.h:1466

References Global_Physical_Variables::P, and res.

Referenced by HomogenisationProblem< ELEMENT >::HomogenisationProblem().

◆ solve() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::solve ( )
inline

Make our own solve function.

265  {
266  //Enable the resolve
268 
269  //Find total number of dofs
270  unsigned long n_dofs = ndof();
271 
272  //Set up the Vector to hold the solution
273  //Vector<double> dx(n_dofs,0.0);
274  DoubleVector dx;
275  //Set up a vector to hold the residuals
276  //Vector<double> res(n_dofs);
278 
279  //Only need to loop over the upper "half" of the matrix
280  for(P=0;P<3;P++)
281  {
282  for(M=P;M<3;M++)
283  {
284  if(M==0)
285  {
286  //Do the solve once only
287  sub_solve(n_dofs,dx,res,true);
288  }
289  //Otherwise it's resolves
290  else
291  {
292  sub_solve(n_dofs,dx,res,false);
293  }
294  }
295  }
296 
297  //Disable the resolve
299  }

References Global_Physical_Variables::P, and res.

◆ sub_solve() [1/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::sub_solve ( const unsigned n_dofs,
DoubleVector dx,
DoubleVector res,
const bool first_solve 
)
inline

Solve the sub-problems This will only solve the problem if the first_solve flag is true

209  {
210  //Update anything that needs updating
212  //Do any updates that are required
214 
215  //Now do the linear solve
216  if(first_solve)
217  {
218  linear_solver_pt()->solve(this,dx);
219  }
220  else
221  {
222  //Get the new residuals
224  //Now do the linear solve
226  }
227 
228  //Subtract the new values from the true dofs
229  for(unsigned l=0;l<n_dofs;l++)
230  {
231  // This is needed during parallel runs when dofs that are not
232  // held on the current processor are nulled out. Can change
233  // this once/if the Dof_pt vector is distributed too.
234  if (Dof_pt[l]!=0) *Dof_pt[l] -= dx[l];
235  }
236 
237 
238 #ifdef OOMPH_HAS_MPI
239  // Synchronise the solution on different processors
240  this->synchronise_all_dofs();
241 #endif
242 
243  // Do any updates that are required
246 
247  // Maximum residuals
248  double maxres=0.0;
249  //Calculate the new residuals
250  get_residuals(dx);
251  //Get the maximum residuals
252  maxres = dx.max();
253 
254  oomph_info << "Final residuals " << maxres << std::endl;
255 
256  //Now update anything that needs updating
258 
259  //Calculate the current contribution to the coefficients
261  }
void actions_before_newton_solve()
Update the problem specs before solve:
Definition: two_dim.cc:187
void actions_after_newton_solve()
Update the problem specs after solve (empty)
Definition: two_dim.cc:184
void calculate_coefficients()
Definition: two_dim.cc:191
double max() const
returns the maximum coefficient
Definition: double_vector.cc:604
virtual void solve(Problem *const &problem_pt, DoubleVector &result)=0
virtual void resolve(const DoubleVector &rhs, DoubleVector &result)
Definition: linear_solver.h:225
virtual void get_residuals(DoubleVector &residuals)
Get the total residuals Vector for the problem.
Definition: problem.cc:3714
virtual void actions_before_newton_step()
Definition: problem.h:1053
Vector< double * > Dof_pt
Vector of pointers to dofs.
Definition: problem.h:554
virtual void actions_before_newton_convergence_check()
Definition: problem.h:1048
virtual void actions_after_newton_step()
Definition: problem.h:1058
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References GlobalFct::get_residuals(), oomph::DoubleVector::max(), oomph::oomph_info, and res.

◆ sub_solve() [2/2]

template<class ELEMENT >
void HomogenisationProblem< ELEMENT >::sub_solve ( const unsigned n_dofs,
DoubleVector dx,
DoubleVector res,
const bool first_solve 
)
inline

Solve the sub-problems This will only solve the problem if the first_solve flag is true

208  {
209  //Update anything that needs updating
211  //Do any updates that are required
213 
214  //Now do the linear solve
215  if(first_solve)
216  {
217  linear_solver_pt()->solve(this,dx);
218  }
219  else
220  {
221  //Get the new residuals
223  //Now do the linear solve
225  }
226 
227  //Subtract the new values from the true dofs
228  for(unsigned l=0;l<n_dofs;l++)
229  {
230  // This is needed during parallel runs when dofs that are not
231  // held on the current processor are nulled out. Can change
232  // this once/if the Dof_pt vector is distributed too.
233  if (Dof_pt[l]!=0) *Dof_pt[l] -= dx[l];
234  }
235 
236 
237 #ifdef OOMPH_HAS_MPI
238  // Synchronise the solution on different processors
239  this->synchronise_all_dofs();
240 #endif
241 
242  // Do any updates that are required
245 
246  // Maximum residuals
247  double maxres=0.0;
248  //Calculate the new residuals
249  get_residuals(dx);
250  //Get the maximum residuals
251  maxres = dx.max();
252 
253  oomph_info << "Final residuals " << maxres << std::endl;
254 
255  //Now update anything that needs updating
257 
258  //Calculate the current contribution to the coefficients
260  }

References GlobalFct::get_residuals(), oomph::DoubleVector::max(), oomph::oomph_info, and res.

Member Data Documentation

◆ Bulk_mesh_pt

template<class ELEMENT >
TriangleMesh< ELEMENT > * HomogenisationProblem< ELEMENT >::Bulk_mesh_pt
private

Pointer to Bulk_mesh.

Referenced by HomogenisationProblem< ELEMENT >::HomogenisationProblem().

◆ C_eff

template<class ELEMENT >
Vector< Vector< DenseMatrix< double > > > HomogenisationProblem< ELEMENT >::C_eff
private

◆ Inclusion_polygon_pt

template<class ELEMENT >
Vector< TriangleMeshPolygon * > HomogenisationProblem< ELEMENT >::Inclusion_polygon_pt
private

Vector storing pointer to the inclusion polygons.

◆ Internal_circle_pt

template<class ELEMENT >
Vector< GeomObject * > HomogenisationProblem< ELEMENT >::Internal_circle_pt
private

◆ M

template<class ELEMENT >
unsigned HomogenisationProblem< ELEMENT >::M
private

◆ Outer_boundary_polyline_pt

template<class ELEMENT >
TriangleMeshPolygon * HomogenisationProblem< ELEMENT >::Outer_boundary_polyline_pt
private

Triangle mesh polygon for outer boundary.

Referenced by HomogenisationProblem< ELEMENT >::HomogenisationProblem().

◆ P

template<class ELEMENT >
unsigned HomogenisationProblem< ELEMENT >::P
private

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