RefineableUnsteadyHeatProblem< ELEMENT > Class Template Reference

Unsteady heat problem in deformable ellipse domain. More...

+ Inheritance diagram for RefineableUnsteadyHeatProblem< ELEMENT >:

Public Member Functions

 RefineableUnsteadyHeatProblem (UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt source_fct_pt)
 Constructor: Pass pointer to source function and timestep. More...
 
 ~RefineableUnsteadyHeatProblem ()
 Destructor: Close trace file. More...
 
void build_mesh ()
 Build mesh function. 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 (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after timestep (empty) More...
 
void actions_before_implicit_timestep ()
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void actions_before_distribute ()
 Actions before distribute: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_distribute ()
 Actions after distribute: Rebuild the mesh of prescribed flux elements. More...
 
double global_temporal_error_norm ()
 Global error norm for adaptive time-stepping. More...
 
void set_initial_condition ()
 
void restart ()
 Restart. More...
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
 
void delete_flux_elements (Mesh *const &surface_mesh_pt)
 Delete UnsteadyHeat flux elements and wipe the surface mesh. More...
 
void doc_solution (const std::string &comment="")
 Doc the solution. More...
 
DocInfodoc_info ()
 Return DocInfo object. More...
 
void dump_it (ofstream &dump_file)
 Dump problem data to allow for later restart. More...
 
void restart (ifstream &restart_file)
 Read problem data for restart. More...
 
RefineableQuarterCircleSectorMesh< ELEMENT > * bulk_mesh_pt ()
 Pointer to bulk mesh. More...
 
doubleepsilon_t ()
 Target error for adaptive timestepping. More...
 
void write_trace_file_header ()
 Write header for trace file. More...
 
doublenext_dt ()
 
 RefineableUnsteadyHeatProblem (UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt source_fct_pt, const bool &use_ale)
 Constructor: Pass pointer to source function and flag to use ALE. More...
 
 ~RefineableUnsteadyHeatProblem ()
 Destructor: Close trace file. 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 (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after timestep (empty) More...
 
void actions_before_implicit_timestep ()
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void set_initial_condition ()
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
 
void delete_flux_elements (Mesh *const &surface_mesh_pt)
 Delete UnsteadyHeat flux elements and wipe the surface mesh. More...
 
void doc_solution ()
 Doc the solution. More...
 
void dump_it (ofstream &dump_file)
 Dump problem data to allow for later restart. More...
 
void restart (ifstream &restart_file)
 Read problem data for restart. More...
 
RefineableQuarterCircleSectorMesh< ELEMENT > * bulk_mesh_pt ()
 Pointer to bulk mesh. More...
 
void possibly_disable_ALE ()
 Switch off ALE terms. More...
 
 RefineableUnsteadyHeatProblem (UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt source_fct_pt)
 Constructor: Pass pointer to source function and timestep. More...
 
 ~RefineableUnsteadyHeatProblem ()
 Destructor: Close trace file. 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 (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after timestep (empty) More...
 
void actions_before_implicit_timestep ()
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
double global_temporal_error_norm ()
 Global error norm for adaptive time-stepping. More...
 
void set_initial_condition ()
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
 
void delete_flux_elements (Mesh *const &surface_mesh_pt)
 Delete UnsteadyHeat flux elements and wipe the surface mesh. More...
 
void doc_solution ()
 Doc the solution. More...
 
void dump_it (ofstream &dump_file)
 Dump problem data to allow for later restart. More...
 
void restart (ifstream &restart_file)
 Read problem data for restart. More...
 
RefineableQuarterCircleSectorMesh< ELEMENT > * bulk_mesh_pt ()
 Pointer to bulk mesh. More...
 
doubleepsilon_t ()
 Target error for adaptive timestepping. More...
 
void write_trace_file_header ()
 Write header for trace file. More...
 
doublenext_dt ()
 
 RefineableUnsteadyHeatProblem (UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt source_fct_pt)
 Constructor: Pass pointer to source function. More...
 
 ~RefineableUnsteadyHeatProblem ()
 Destructor: Close trace file. 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 (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after timestep (empty) More...
 
void actions_before_implicit_timestep ()
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void set_initial_condition ()
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
 
void delete_flux_elements (Mesh *const &surface_mesh_pt)
 Delete UnsteadyHeat flux elements and wipe the surface mesh. More...
 
void doc_solution ()
 Doc the solution. More...
 
void dump_it (ofstream &dump_file)
 Dump problem data to allow for later restart. More...
 
void restart (ifstream &restart_file)
 Read problem data for restart. More...
 
RefineableQuarterCircleSectorMesh< ELEMENT > * bulk_mesh_pt ()
 Pointer to bulk mesh. More...
 
 RefineableUnsteadyHeatProblem (UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt source_fct_pt)
 Constructor: Pass pointer to source function. More...
 
 ~RefineableUnsteadyHeatProblem ()
 Destructor: Close trace file. 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 (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after timestep (empty) More...
 
void actions_before_implicit_timestep ()
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void set_initial_condition ()
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
 
void delete_flux_elements (Mesh *const &surface_mesh_pt)
 Delete UnsteadyHeat flux elements and wipe the surface mesh. More...
 
void doc_solution ()
 Doc the solution. More...
 
void dump_it (ofstream &dump_file)
 Dump problem data to allow for later restart. More...
 
void restart (ifstream &restart_file)
 Read problem data for restart. More...
 
RefineableQuarterCircleSectorMesh< ELEMENT > * bulk_mesh_pt ()
 Pointer to bulk mesh. More...
 
- 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
 
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 Member Functions

void generic_actions_before ()
 Actions before adapt/distribute: Wipe the mesh of prescribed flux elements. More...
 
void generic_actions_after ()
 

Private Attributes

DocInfo Doc_info
 Doc info object. More...
 
double Next_dt
 
GeomObjectBoundary_pt
 Pointer to GeomObject that specifies the domain bondary. More...
 
UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt Source_fct_pt
 Pointer to source function. More...
 
RefineableQuarterCircleSectorMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
MeshSurface_mesh_pt
 Pointer to the "surface" mesh. More...
 
NodeDoc_node_pt
 Pointer to central node (exists at all refinement levels) for doc. More...
 
ofstream Trace_file
 Trace file. More...
 
double Epsilon_t
 Target error for adaptive timestepping. More...
 
bool Use_ALE
 Flag to use ALE. 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_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)
 
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 RefineableUnsteadyHeatProblem< ELEMENT >

Unsteady heat problem in deformable ellipse domain.

Unsteady heat problem in quarter circle domain.

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////

Constructor & Destructor Documentation

◆ RefineableUnsteadyHeatProblem() [1/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::RefineableUnsteadyHeatProblem ( UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt  source_fct_pt)

Constructor: Pass pointer to source function and timestep.

Constructor for UnsteadyHeat problem on deformable ellipse domain. Pass pointer to source function.

Constructor for UnsteadyHeat problem in quarter circle domain. Pass pointer to source function.

370 {
371 
372  // Setup labels for output
373  //------------------------
374 
375  // Output directory
376  Doc_info.set_directory("RESLT");
377 
378  // Output number
379  Doc_info.number()=0;
380 
381 
382  // Setup timestepping
383  //-------------------
384 
385  // Allocate the timestepper -- this constructs the time object as well.
386  // The boolean flag to the constructor enables adaptivity.
387  add_time_stepper_pt(new BDF<2>(true));
388 
389  // Set target error for adaptive timestepping
390  Epsilon_t=1.0e-2;
391 
392  // Setup parameters for tanh solution
393  // ----------------------------------
394 
395  // Steepness of step
397 
398  // Orientation of step
400 
401  // Amplitude for movement of step
403 
404  // Parameter for time-dependence of step movement
406 
407  // Setup mesh
408  //-----------
409 
410  // Build geometric object that forms the curvilinear domain boundary:
411  // an oscillating ellipse
412 
413  // Half axes
414  double a=1.0;
415  double b=1.0;
416 
417  // Variations of half axes
418  double a_hat= 0.1;
419  double b_hat=-0.1;
420 
421  // Period of the oscillation
422  double period=1.0;
423 
424  // Create GeomObject
425  Boundary_pt=new MyEllipse(a,b,a_hat,b_hat,period,Problem::time_pt());
426 
427  // Now build the mesh
428  build_mesh();
429 
430  // Do equation numbering
431  oomph_info <<"Number of equations: " << assign_eqn_numbers() << std::endl;
432 
433 } // end of constructor
Scalar * b
Definition: benchVecAdd.cpp:17
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:51
UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:343
double Epsilon_t
Target error for adaptive timestepping.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:358
GeomObject * Boundary_pt
Pointer to GeomObject that specifies the domain bondary.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:340
DocInfo Doc_info
Doc info object.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:333
void build_mesh()
Build mesh function.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:453
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
void add_time_stepper_pt(TimeStepper *const &time_stepper_pt)
Definition: problem.cc:1545
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
const Scalar * a
Definition: level2_cplx_impl.h:32
double TanPhi
Parameter for angle of step.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:147
double Beta
Parameter for amplitude of step translation.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:141
double Gamma
Parameter for timescale of step translation.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:144
double Alpha
Parameter for steepness of step.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:138
GeneralisedAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: gen_axisym_advection_diffusion_elements.h:229
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References a, oomph::Problem::add_time_stepper_pt(), GlobalParameters::Alpha, oomph::Problem::assign_eqn_numbers(), b, GlobalParameters::Beta, RefineableUnsteadyHeatProblem< ELEMENT >::Boundary_pt, RefineableUnsteadyHeatProblem< ELEMENT >::build_mesh(), RefineableUnsteadyHeatProblem< ELEMENT >::Doc_info, RefineableUnsteadyHeatProblem< ELEMENT >::Epsilon_t, GlobalParameters::Gamma, oomph::DocInfo::number(), oomph::oomph_info, oomph::DocInfo::set_directory(), and GlobalParameters::TanPhi.

◆ ~RefineableUnsteadyHeatProblem() [1/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::~RefineableUnsteadyHeatProblem

Destructor: Close trace file.

441 {
442 
443  // Close trace file
444  Trace_file.close();
445 
446 } // end of destructor
ofstream Trace_file
Trace file.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:355

References oomph::Problem_Parameter::Trace_file.

◆ RefineableUnsteadyHeatProblem() [2/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::RefineableUnsteadyHeatProblem ( UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt  source_fct_pt,
const bool use_ale 
)

Constructor: Pass pointer to source function and flag to use ALE.

Constructor for UnsteadyHeat problem in quarter circle domain. Pass pointer to source function and flag to use ALE.

312  :
314 {
315 
316 
317  // Setup labels for output
318  //------------------------
319  // Output directory
320  if (!use_ale)
321  {
322  Doc_info.set_directory("RESLT_adapt");
323  }
324  else
325  {
326  Doc_info.set_directory("RESLT_adapt_ALE");
327  }
328 
329  // Output number
330  Doc_info.number()=0;
331 
332  // Open trace file
333  char filename[100];
334  sprintf(filename,"%s/trace.dat",Doc_info.directory().c_str());
335  Trace_file.open(filename);
336 
337  Trace_file << "VARIABLES=\"time t\",\"u<SUB>FE</SUB>\",\"u<SUB>exact</SUB>\","
338  << "\"A\","
339  << "\"X<SUB>step</SUB>\","
340  << "\"N<SUB>element</SUB>\","
341  << "\"N<SUB>refined</SUB>\","
342  << "\"N<SUB>unrefined</SUB>\","
343  << "\"norm of error\","
344  << "\"norm of solution\""
345  << std::endl;
346 
347 
348  // Setup parameters for tanh solution
349  // ----------------------------------
350 
351  // Steepness of step
353 
354  // Orientation of step
356 
357  // Amplitude for movement of step
359 
360  // Parameter for time-dependence of step movement
362 
363 
364  //Allocate the timestepper -- This constructs the time object as well
366 
367  // Setup mesh
368  //-----------
369 
370  // Build geometric object that forms the curvilinear domain boundary:
371  // a unit circle
372 
373  // Create GeomObject
375 
376  // Start and end coordinates of curvilinear domain boundary on circle
377  double xi_lo=0.0;
378  double xi_hi=MathematicalConstants::Pi/2.0;
379 
380  // Now create the bulk mesh. Separating line between the two
381  // elements next to the curvilinear boundary is located half-way
382  // along the boundary.
383  double fract_mid=0.5;
385  Boundary_pt,xi_lo,fract_mid,xi_hi,time_stepper_pt());
386 
387  // Create the surface mesh as an empty mesh
388  Surface_mesh_pt=new Mesh;
389 
390  // Create prescribed-flux elements from all elements that are
391  // adjacent to boundary 0 (the horizontal lower boundary), and add them
392  // to the (so far empty) surface mesh.
394 
395  // Add the two sub meshes to the problem
398 
399  // Combine all submeshes into a single global Mesh
401 
402  // Set error estimator for bulk mesh
404  Bulk_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
405 
406  // Set the boundary conditions for this problem: All nodes are
407  // free by default -- just pin the ones that have Dirichlet conditions
408  // here.
409  unsigned n_bound = Bulk_mesh_pt->nboundary();
410  for(unsigned b=0;b<n_bound;b++)
411  {
412  // Leave nodes on boundary 0 free -- this is where we apply the flux
413  // boundary condition
414  if (b!=0)
415  {
416  unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
417  for (unsigned n=0;n<n_node;n++)
418  {
419  Bulk_mesh_pt->boundary_node_pt(b,n)->pin(0);
420  }
421  }
422  }
423 
424  // Extract pointer to the central node (this exists at all refinement levels)
425  // for doc of solution
426  FiniteElement* el0_pt=Bulk_mesh_pt->finite_element_pt(0);
427  unsigned nnod=el0_pt->nnode();
428  Doc_node_pt=el0_pt->node_pt(nnod-1);
429 
430 
431  // Complete the build of all elements so they are fully functional
432  //----------------------------------------------------------------
433 
434  // Find number of elements in mesh
435  unsigned n_element = Bulk_mesh_pt->nelement();
436 
437  // Loop over the elements to set up element-specific
438  // things that cannot be handled by constructor
439  for(unsigned i=0;i<n_element;i++)
440  {
441  // Upcast from FiniteElement to the present element
442  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(i));
443 
444  //Set the source function pointer
445  el_pt->source_fct_pt() = Source_fct_pt;
446  }
447 
448  // Loop over the flux elements to pass pointer to prescribed flux function
449  n_element=Surface_mesh_pt->nelement();
450  for(unsigned e=0;e<n_element;e++)
451  {
452  // Upcast from GeneralisedElement to UnsteadyHeat flux element
454  dynamic_cast<UnsteadyHeatFluxElement<ELEMENT>*>(
456 
457  // Set the pointer to the prescribed flux function
458  el_pt->flux_fct_pt() =
460  }
461 
462 
463  // Switch off ALE if required
465 
466  // Do equation numbering
467  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
468 
469 } // end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:54
Node * Doc_node_pt
Pointer to central node (exists at all refinement levels) for doc.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:352
Mesh * Surface_mesh_pt
Pointer to the "surface" mesh.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:349
void create_flux_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:1007
RefineableQuarterCircleSectorMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:346
bool Use_ALE
Flag to use ALE.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:301
void possibly_disable_ALE()
Switch off ALE terms.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:251
std::string directory() const
Output directory.
Definition: oomph_utilities.h:524
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Definition: mesh.h:67
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
TimeStepper *& time_stepper_pt()
Definition: problem.h:1524
Definition: quarter_circle_sector_mesh.template.h:188
Definition: unsteady_heat_flux_elements.h:56
UnsteadyHeatPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
Definition: unsteady_heat_flux_elements.h:82
Definition: error_estimator.h:266
double Pi
Definition: two_d_biharmonic.cc:235
string filename
Definition: MergeRestartFiles.py:39
Z2ErrorEstimator * error_estimator_pt
Definition: MortaringCantileverCompareToNonMortaring.cpp:190
double Alpha
Parameter for steepness of step.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:103
double Gamma
Parameter for timescale of step translation.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:109
double Beta
Parameter for amplitude of step translation.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:106
void prescribed_flux_on_fixed_y_boundary(const double &time, const Vector< double > &x, double &flux)
Flux required by the exact solution on a boundary on which y is fixed.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:157
double TanPhi
Parameter for angle of step.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:112

References oomph::Problem::add_sub_mesh(), oomph::Problem::add_time_stepper_pt(), TanhSolnForUnsteadyHeat::Alpha, oomph::Problem::assign_eqn_numbers(), b, TanhSolnForUnsteadyHeat::Beta, RefineableUnsteadyHeatProblem< ELEMENT >::Boundary_pt, oomph::Problem::build_global_mesh(), RefineableUnsteadyHeatProblem< ELEMENT >::Bulk_mesh_pt, RefineableUnsteadyHeatProblem< ELEMENT >::create_flux_elements(), oomph::DocInfo::directory(), RefineableUnsteadyHeatProblem< ELEMENT >::Doc_info, RefineableUnsteadyHeatProblem< ELEMENT >::Doc_node_pt, e(), oomph::Mesh::element_pt(), MeshRefinement::error_estimator_pt, MergeRestartFiles::filename, oomph::UnsteadyHeatFluxElement< ELEMENT >::flux_fct_pt(), TanhSolnForUnsteadyHeat::Gamma, i, n, oomph::Mesh::nelement(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::DocInfo::number(), BiharmonicTestFunctions2::Pi, RefineableUnsteadyHeatProblem< ELEMENT >::possibly_disable_ALE(), TanhSolnForUnsteadyHeat::prescribed_flux_on_fixed_y_boundary(), oomph::DocInfo::set_directory(), RefineableUnsteadyHeatProblem< ELEMENT >::Source_fct_pt, RefineableUnsteadyHeatProblem< ELEMENT >::Surface_mesh_pt, TanhSolnForUnsteadyHeat::TanPhi, oomph::Problem::time_stepper_pt(), and RefineableUnsteadyHeatProblem< ELEMENT >::Trace_file.

◆ ~RefineableUnsteadyHeatProblem() [2/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::~RefineableUnsteadyHeatProblem ( )

Destructor: Close trace file.

◆ RefineableUnsteadyHeatProblem() [3/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::RefineableUnsteadyHeatProblem ( UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt  source_fct_pt)

Constructor: Pass pointer to source function and timestep.

◆ ~RefineableUnsteadyHeatProblem() [3/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::~RefineableUnsteadyHeatProblem ( )

Destructor: Close trace file.

◆ RefineableUnsteadyHeatProblem() [4/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::RefineableUnsteadyHeatProblem ( UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt  source_fct_pt)

Constructor: Pass pointer to source function.

◆ ~RefineableUnsteadyHeatProblem() [4/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::~RefineableUnsteadyHeatProblem ( )

Destructor: Close trace file.

◆ RefineableUnsteadyHeatProblem() [5/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::RefineableUnsteadyHeatProblem ( UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt  source_fct_pt)

Constructor: Pass pointer to source function.

◆ ~RefineableUnsteadyHeatProblem() [5/5]

template<class ELEMENT >
RefineableUnsteadyHeatProblem< ELEMENT >::~RefineableUnsteadyHeatProblem ( )

Destructor: Close trace file.

Member Function Documentation

◆ actions_after_adapt() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_adapt
inlinevirtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

260  {
262  }
void generic_actions_after()
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:620

◆ actions_after_adapt() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_after_adapt() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_after_adapt() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_after_adapt() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_after_distribute()

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_distribute ( )
inline

Actions after distribute: Rebuild the mesh of prescribed flux elements.

272  {
274  }

◆ actions_after_implicit_timestep() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after timestep (empty)

Reimplemented from oomph::Problem.

246 {}

◆ actions_after_implicit_timestep() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after timestep (empty)

Reimplemented from oomph::Problem.

207 {}

◆ actions_after_implicit_timestep() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after timestep (empty)

Reimplemented from oomph::Problem.

234 {}

◆ actions_after_implicit_timestep() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after timestep (empty)

Reimplemented from oomph::Problem.

206 {}

◆ actions_after_implicit_timestep() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after timestep (empty)

Reimplemented from oomph::Problem.

234 {}

◆ actions_after_newton_solve() [1/5]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

240 {}

◆ actions_after_newton_solve() [2/5]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

201 {}

◆ actions_after_newton_solve() [3/5]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

228 {}

◆ actions_after_newton_solve() [4/5]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

200 {}

◆ actions_after_newton_solve() [5/5]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

228 {}

◆ actions_before_adapt() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_adapt
inlinevirtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

254  {
256  }
void generic_actions_before()
Actions before adapt/distribute: Wipe the mesh of prescribed flux elements.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:603

◆ actions_before_adapt() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_before_adapt() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_before_adapt() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_before_adapt() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_before_distribute()

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_distribute ( )
inline

Actions before distribute: Wipe the mesh of prescribed flux elements.

266  {
268  }

◆ actions_before_implicit_timestep() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep
virtual

Update the problem specs before next timestep: Set Dirchlet boundary conditions from exact solution.

Actions before timestep: Update the domain shape, then set the boundary conditions for the current time.

Actions before timestep: Set the boundary conditions for the current time.

Reimplemented from oomph::Problem.

567 {
568  // Update the domain shape
569  Bulk_mesh_pt->node_update();
570 
571  // Get current time
572  double time=time_pt()->time();
573 
574  //Loop over all boundaries
575  unsigned num_bound = Bulk_mesh_pt->nboundary();
576  for(unsigned b=0;b<num_bound;b++)
577  {
578  // Exclude flux boundary
579  if (b!=0)
580  {
581  // Loop over the nodes on boundary
582  unsigned num_nod=Bulk_mesh_pt->nboundary_node(b);
583  for (unsigned j=0;j<num_nod;j++)
584  {
585  Node* nod_pt=Bulk_mesh_pt->boundary_node_pt(b,j);
586  double u;
587  Vector<double> x(2);
588  x[0]=nod_pt->x(0);
589  x[1]=nod_pt->x(1);
591  nod_pt->set_value(0,u);
592  }
593  }
594  }
595 } // end of actions_before_implicit_timestep
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
Time *& time_pt()
Return a pointer to the global time object.
Definition: problem.h:1504
double & time()
Return the current value of continuous time.
Definition: problem.cc:11531
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
void get_exact_u(const double &time, const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
Definition: stefan_boltzmann.cc:169
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References b, GlobalParameters::get_exact_u(), j, oomph::Data::set_value(), plotDoE::x, and oomph::Node::x().

◆ actions_before_implicit_timestep() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep ( )
virtual

Update the problem specs before next timestep: Set Dirchlet boundary conditions from exact solution.

Reimplemented from oomph::Problem.

◆ actions_before_implicit_timestep() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep ( )
virtual

Update the problem specs before next timestep: Set Dirchlet boundary conditions from exact solution.

Reimplemented from oomph::Problem.

◆ actions_before_implicit_timestep() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep ( )
virtual

Update the problem specs before next timestep: Set Dirchlet boundary conditions from exact solution.

Reimplemented from oomph::Problem.

◆ actions_before_implicit_timestep() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep ( )
virtual

Update the problem specs before next timestep: Set Dirchlet boundary conditions from exact solution.

Reimplemented from oomph::Problem.

◆ actions_before_newton_solve() [1/5]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

243 {}

◆ actions_before_newton_solve() [2/5]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

204 {}

◆ actions_before_newton_solve() [3/5]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

231 {}

◆ actions_before_newton_solve() [4/5]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

203 {}

◆ actions_before_newton_solve() [5/5]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

231 {}

◆ build_mesh()

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::build_mesh

Build mesh function.

Build meshes (helper fct accessed from constructor and from the load balancing routines)

454 {
455 
456  // Start and end coordinates of curvilinear domain boundary on ellipse
457  double xi_lo=0.0;
458  double xi_hi=MathematicalConstants::Pi/2.0;
459 
460  // Now create the bulk mesh. Separating line between the two
461  // elements next to the curvilinear boundary is located half-way
462  // along the boundary.
463  double fract_mid=0.5;
465  Boundary_pt,xi_lo,fract_mid,xi_hi,time_stepper_pt());
466 
467  // Note: Mesh is completely rebuilt in here
468  // so target errors etc need to be reassigned
469 
470  // Set targets for spatial adaptivity
471  Bulk_mesh_pt->max_permitted_error()=0.001;
472  Bulk_mesh_pt->min_permitted_error()=0.0001;
473 
474  // Need to do this in build_mesh() because it's supposed to
475  // get the problem into the state it was when it was distributed.
476  Bulk_mesh_pt->refine_uniformly();
477  Bulk_mesh_pt->refine_uniformly();
478 
479  // Create the surface mesh as an empty mesh
480  Surface_mesh_pt=new Mesh;
481 
482  // Create prescribed-flux elements from all elements that are
483  // adjacent to boundary 0 (the horizontal lower boundary), and add them
484  // to the (so far empty) surface mesh.
486 
487  // Add the two sub meshes to the problem
490 
491  // Combine all submeshes into a single global Mesh
493 
494  // Set error estimator for bulk mesh
496  Bulk_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
497 
498  // Set the boundary conditions for this problem: All nodes are
499  // free by default -- just pin the ones that have Dirichlet conditions
500  // here.
501  unsigned n_bound = Bulk_mesh_pt->nboundary();
502  for(unsigned b=0;b<n_bound;b++)
503  {
504  // Leave nodes on boundary 0 free -- this is where we apply the flux
505  // boundary condition
506  if (b!=0)
507  {
508  unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
509  for (unsigned n=0;n<n_node;n++)
510  {
511  Bulk_mesh_pt->boundary_node_pt(b,n)->pin(0);
512  }
513  }
514  }
515 
516  // Extract pointer to the central node (this exists at all refinement levels)
517  // for doc of solution
518  FiniteElement* el0_pt=Bulk_mesh_pt->finite_element_pt(0);
519  unsigned nnod=el0_pt->nnode();
520  Doc_node_pt=el0_pt->node_pt(nnod-1);
521 
522  // In parallel make sure that we retain this element on all processors
523  // so the node remains accessible for post-processing
524  el0_pt->must_be_kept_as_halo();
525 
526  // Complete the build of all elements so they are fully functional
527  //----------------------------------------------------------------
528 
529  // Find number of elements in mesh
530  unsigned n_element = Bulk_mesh_pt->nelement();
531 
532  // Loop over the elements to set up element-specific
533  // things that cannot be handled by constructor
534  for(unsigned i=0;i<n_element;i++)
535  {
536  // Upcast from FiniteElement to the present element
537  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(i));
538 
539  //Set the source function pointer
540  el_pt->source_fct_pt() = Source_fct_pt;
541  }
542 
543 
544  // Loop over the flux elements to pass pointer to prescribed flux function
545  n_element=Surface_mesh_pt->nelement();
546  for(unsigned e=0;e<n_element;e++)
547  {
548  // Upcast from GeneralisedElement to UnsteadyHeat flux element
550  dynamic_cast<UnsteadyHeatFluxElement<ELEMENT>*>(
552 
553  // Set the pointer to the prescribed flux function
554  el_pt->flux_fct_pt() =
556  }
557 
558 } // end of build_mesh
void prescribed_flux_on_fixed_y_boundary(const double &time, const Vector< double > &x, double &flux)
Flux required by the exact solution on a boundary on which y is fixed.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:192

References b, e(), MeshRefinement::error_estimator_pt, oomph::UnsteadyHeatFluxElement< ELEMENT >::flux_fct_pt(), i, n, oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), BiharmonicTestFunctions2::Pi, GlobalParameters::prescribed_flux_on_fixed_y_boundary(), and oomph::Source_fct_pt.

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

◆ bulk_mesh_pt() [1/5]

template<class ELEMENT >
RefineableQuarterCircleSectorMesh<ELEMENT>* RefineableUnsteadyHeatProblem< ELEMENT >::bulk_mesh_pt ( )
inline

Pointer to bulk mesh.

309  {
310  return Bulk_mesh_pt;
311  }

◆ bulk_mesh_pt() [2/5]

template<class ELEMENT >
RefineableQuarterCircleSectorMesh<ELEMENT>* RefineableUnsteadyHeatProblem< ELEMENT >::bulk_mesh_pt ( )
inline

Pointer to bulk mesh.

246  {
247  return Bulk_mesh_pt;
248  }

◆ bulk_mesh_pt() [3/5]

template<class ELEMENT >
RefineableQuarterCircleSectorMesh<ELEMENT>* RefineableUnsteadyHeatProblem< ELEMENT >::bulk_mesh_pt ( )
inline

Pointer to bulk mesh.

273  {
274  return Bulk_mesh_pt;
275  }

◆ bulk_mesh_pt() [4/5]

template<class ELEMENT >
RefineableQuarterCircleSectorMesh<ELEMENT>* RefineableUnsteadyHeatProblem< ELEMENT >::bulk_mesh_pt ( )
inline

Pointer to bulk mesh.

245  {
246  return Bulk_mesh_pt;
247  }

◆ bulk_mesh_pt() [5/5]

template<class ELEMENT >
RefineableQuarterCircleSectorMesh<ELEMENT>* RefineableUnsteadyHeatProblem< ELEMENT >::bulk_mesh_pt ( )
inline

Pointer to bulk mesh.

273  {
274  return Bulk_mesh_pt;
275  }

◆ create_flux_elements() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)

Create UnsteadyHeat flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the Mesh object pointed to by surface_mesh_pt

Create UnsteadyHeat Flux Elements on the b-th boundary of the Mesh object pointed to by bulk_mesh_pt and add the elements to the Mesh object pointed to by surface_mesh_pt.

1009 {
1010  // How many bulk elements are adjacent to boundary b?
1011  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
1012 
1013  // Loop over the bulk elements adjacent to boundary b?
1014  for(unsigned e=0;e<n_element;e++)
1015  {
1016  // Get pointer to the bulk element that is adjacent to boundary b
1017  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
1018  bulk_mesh_pt->boundary_element_pt(b,e));
1019 
1020  //What is the face index of element e on boundary b
1021  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
1022 
1023  // Build the corresponding prescribed-flux element
1024  UnsteadyHeatFluxElement<ELEMENT>* flux_element_pt = new
1025  UnsteadyHeatFluxElement<ELEMENT>(bulk_elem_pt,face_index);
1026 
1027  //Add the prescribed-flux element to the surface mesh
1028  surface_mesh_pt->add_element_pt(flux_element_pt);
1029 
1030  } //end of loop over bulk elements adjacent to boundary b
1031 
1032 } // end of create_flux_elements
RefineableQuarterCircleSectorMesh< ELEMENT > * bulk_mesh_pt()
Pointer to bulk mesh.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:308
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617

References oomph::Mesh::add_element_pt(), b, oomph::Mesh::boundary_element_pt(), e(), oomph::Mesh::face_index_at_boundary(), and oomph::Mesh::nboundary_element().

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

◆ create_flux_elements() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)

Create UnsteadyHeat flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the Mesh object pointed to by surface_mesh_pt

◆ create_flux_elements() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)

Create UnsteadyHeat flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the Mesh object pointed to by surface_mesh_pt

◆ create_flux_elements() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)

Create UnsteadyHeat flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the Mesh object pointed to by surface_mesh_pt

◆ create_flux_elements() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)

Create UnsteadyHeat flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the Mesh object pointed to by surface_mesh_pt

◆ delete_flux_elements() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::delete_flux_elements ( Mesh *const &  surface_mesh_pt)

Delete UnsteadyHeat flux elements and wipe the surface mesh.

Delete UnsteadyHeat Flux Elements and wipe the surface mesh.

1041 {
1042  // How many surface elements are in the surface mesh
1043  unsigned n_element = surface_mesh_pt->nelement();
1044 
1045  // Loop over the surface elements
1046  for(unsigned e=0;e<n_element;e++)
1047  {
1048  // Kill surface element
1049  delete surface_mesh_pt->element_pt(e);
1050  }
1051 
1052  // Wipe the mesh
1053  surface_mesh_pt->flush_element_and_node_storage();
1054 
1055 } // end of delete_flux_elements
void flush_element_and_node_storage()
Definition: mesh.h:407

References e(), oomph::Mesh::element_pt(), oomph::Mesh::flush_element_and_node_storage(), and oomph::Mesh::nelement().

◆ delete_flux_elements() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::delete_flux_elements ( Mesh *const &  surface_mesh_pt)

Delete UnsteadyHeat flux elements and wipe the surface mesh.

◆ delete_flux_elements() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::delete_flux_elements ( Mesh *const &  surface_mesh_pt)

Delete UnsteadyHeat flux elements and wipe the surface mesh.

◆ delete_flux_elements() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::delete_flux_elements ( Mesh *const &  surface_mesh_pt)

Delete UnsteadyHeat flux elements and wipe the surface mesh.

◆ delete_flux_elements() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::delete_flux_elements ( Mesh *const &  surface_mesh_pt)

Delete UnsteadyHeat flux elements and wipe the surface mesh.

◆ doc_info()

template<class ELEMENT >
DocInfo& RefineableUnsteadyHeatProblem< ELEMENT >::doc_info ( )
inline

Return DocInfo object.

299 {return Doc_info;}

References GlobalParameters::Doc_info.

◆ doc_solution() [1/5]

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

Doc the solution.

704 {
705  ofstream some_file;
706  char filename[100];
707 
708  // Number of plot points
709  unsigned npts;
710  npts=5;
711 
712  cout << std::endl;
713  cout << "=================================================" << std::endl;
714  cout << "Docing solution for t=" << time_pt()->time() << std::endl;
715  cout << "=================================================" << std::endl;
716 
717  // Output solution
718  //-----------------
719  sprintf(filename,"%s/soln%i.dat",Doc_info.directory().c_str(),
720  Doc_info.number());
721  some_file.open(filename);
722  Bulk_mesh_pt->output(some_file,npts);
723  some_file << "TEXT X=2.5,Y=93.6,F=HELV,HU=POINT,C=BLUE,H=26,T=\"time = "
724  << time_pt()->time() << "\"";
725  some_file << "GEOMETRY X=2.5,Y=98,T=LINE,C=BLUE,LT=0.4" << std::endl;
726  some_file << "1" << std::endl;
727  some_file << "2" << std::endl;
728  some_file << " 0 0" << std::endl;
729  some_file << time_pt()->time()*20.0 << " 0" << std::endl;
730 
731  // Write dummy zones
732  some_file << "ZONE I=2,J=2" << std::endl;
733  some_file << "0.0 0.0 -1.2" << std::endl;
734  some_file << "1.3 0.0 -1.2" << std::endl;
735  some_file << "0.0 1.3 -1.2" << std::endl;
736  some_file << "1.3 1.3 -1.2" << std::endl;
737  some_file << "ZONE I=2,J=2" << std::endl;
738  some_file << "0.0 0.0 1.2" << std::endl;
739  some_file << "1.3 0.0 1.2" << std::endl;
740  some_file << "0.0 1.3 1.2" << std::endl;
741  some_file << "1.3 1.3 1.2" << std::endl;
742 
743  some_file.close();
744 
745 
746  // Output exact solution
747  //----------------------
748  sprintf(filename,"%s/exact_soln%i.dat",Doc_info.directory().c_str(),
749  Doc_info.number());
750  some_file.open(filename);
751  Bulk_mesh_pt->output_fct(some_file,npts,time_pt()->time(),
753 
754  // Write dummy zones
755  some_file << "ZONE I=2,J=2" << std::endl;
756  some_file << "0.0 0.0 -1.2" << std::endl;
757  some_file << "1.3 0.0 -1.2" << std::endl;
758  some_file << "0.0 1.3 -1.2" << std::endl;
759  some_file << "1.3 1.3 -1.2" << std::endl;
760  some_file << "ZONE I=2,J=2" << std::endl;
761  some_file << "0.0 0.0 1.2" << std::endl;
762  some_file << "1.3 0.0 1.2" << std::endl;
763  some_file << "0.0 1.3 1.2" << std::endl;
764  some_file << "1.3 1.3 1.2" << std::endl;
765 
766  some_file.close();
767 
768 
769  // Doc error
770  //----------
771  double error,norm;
772  sprintf(filename,"%s/error%i.dat",Doc_info.directory().c_str(),
773  Doc_info.number());
774  some_file.open(filename);
775  Bulk_mesh_pt->compute_error(some_file,
777  time_pt()->time(),
778  error,norm);
779  some_file.close();
780 
781 
782 
783  // Doc error and write trace file
784  //--------------------------------
785  cout << "error: " << error << std::endl;
786  cout << "norm : " << norm << std::endl << std::endl;
787 
788  Vector<double> x(2);
789  x[0]=Doc_node_pt->x(0);
790  x[1]=Doc_node_pt->x(1);
791  double u_exact;
793  Vector<double > xi_wall(1);
794  Vector<double > r_wall(2);
795  xi_wall[0]=0.0;
796  Boundary_pt->position(xi_wall,r_wall);
797  Trace_file << time_pt()->time()
798  << " " << Doc_node_pt->value(0)
799  << " " << u_exact
800  << " " << r_wall[0]
802  << " " << Bulk_mesh_pt->nelement()
803  << " " << Bulk_mesh_pt->nrefined()
804  << " " << Bulk_mesh_pt->nunrefined()
805  << " " << error
806  << " " << norm << std::endl;
807 
808  // Plot wall posn
809  //---------------
810  sprintf(filename,"%s/Wall%i.dat",Doc_info.directory().c_str(),
811  Doc_info.number());
812  some_file.open(filename);
813 
814  unsigned nplot=100;
815  for (unsigned iplot=0;iplot<nplot;iplot++)
816  {
817  xi_wall[0]=0.5*Pi*double(iplot)/double(nplot-1);
818  Boundary_pt->position(xi_wall,r_wall);
819  some_file << r_wall[0] << " " << r_wall[1] << std::endl;
820  }
821  some_file.close();
822 
823  // Write restart file
824  sprintf(filename,"%s/restart%i.dat",Doc_info.directory().c_str(),
825  Doc_info.number());
826  some_file.open(filename);
827  dump_it(some_file);
828  some_file.close();
829 
830  // Increment number of doc
831  Doc_info.number()++;
832 
833 
834 } // end of doc_solution
void dump_it(ofstream &dump_file)
Dump problem data to allow for later restart.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:1146
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
double value(const unsigned &i) const
Definition: nodes.cc:2408
void get_exact_u(const double &time, const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:121
double step_position(const double &time)
Position of step (x-axis intercept)
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:115
int error
Definition: calibrate.py:297

References oomph::DocInfo::directory(), GlobalParameters::Doc_info, calibrate::error, MergeRestartFiles::filename, TanhSolnForUnsteadyHeat::get_exact_u(), oomph::DocInfo::number(), BiharmonicTestFunctions2::Pi, TanhSolnForUnsteadyHeat::step_position(), oomph::Problem_Parameter::Trace_file, and plotDoE::x.

◆ doc_solution() [2/5]

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

Doc the solution.

◆ doc_solution() [3/5]

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

Doc the solution.

◆ doc_solution() [4/5]

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

Doc the solution.

◆ doc_solution() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::doc_solution ( const std::string &  comment = "")

Doc the solution.

830 {
831  ofstream some_file;
832  char filename[100];
833 
834  // Number of plot points
835  unsigned npts;
836  npts=5;
837 
838  oomph_info << std::endl;
839  oomph_info << "================================================="
840  << std::endl;
841  oomph_info << "Docing solution " << Doc_info.number() << " for t="
842  << time_pt()->time() << " [ndof=" << ndof() << "]\n";
843  oomph_info << "================================================="
844  << std::endl;
845 
846  // Doc halo schemes
847  //-----------------
848  Bulk_mesh_pt->doc_mesh_distribution(Doc_info);
849 
850  // Output solution
851  //-----------------
852  sprintf(filename,"%s/soln%i_on_proc%i.dat",Doc_info.directory().c_str(),
853  Doc_info.number(),this->communicator_pt()->my_rank());
854  some_file.open(filename);
855  Bulk_mesh_pt->output(some_file,npts);
856  some_file << "TEXT X=2.5,Y=93.6,F=HELV,HU=POINT,C=BLUE,H=26,T=\"time = "
857  << time_pt()->time() << "; "<< comment << "\"";
858  some_file << "GEOMETRY X=2.5,Y=98,T=LINE,C=BLUE,LT=0.4" << std::endl;
859  some_file << "1" << std::endl;
860  some_file << "2" << std::endl;
861  some_file << " 0 0" << std::endl;
862  some_file << time_pt()->time()*20.0 << " 0" << std::endl;
863 
864  // Write dummy zones
865  some_file << "ZONE I=2,J=2" << std::endl;
866  some_file << "0.0 0.0 -1.2" << std::endl;
867  some_file << "1.3 0.0 -1.2" << std::endl;
868  some_file << "0.0 1.3 -1.2" << std::endl;
869  some_file << "1.3 1.3 -1.2" << std::endl;
870  some_file << "ZONE I=2,J=2" << std::endl;
871  some_file << "0.0 0.0 1.2" << std::endl;
872  some_file << "1.3 0.0 1.2" << std::endl;
873  some_file << "0.0 1.3 1.2" << std::endl;
874  some_file << "1.3 1.3 1.2" << std::endl;
875 
876  some_file.close();
877 
878 
879  // Output exact solution
880  //----------------------
881  sprintf(filename,"%s/exact_soln%i_on_proc%i.dat",Doc_info.directory().c_str(),
882  Doc_info.number(),this->communicator_pt()->my_rank());
883  some_file.open(filename);
884  Bulk_mesh_pt->output_fct(some_file,npts,time_pt()->time(),
886 
887  // Write dummy zones
888  some_file << "ZONE I=2,J=2" << std::endl;
889  some_file << "0.0 0.0 -1.2" << std::endl;
890  some_file << "1.3 0.0 -1.2" << std::endl;
891  some_file << "0.0 1.3 -1.2" << std::endl;
892  some_file << "1.3 1.3 -1.2" << std::endl;
893  some_file << "ZONE I=2,J=2" << std::endl;
894  some_file << "0.0 0.0 1.2" << std::endl;
895  some_file << "1.3 0.0 1.2" << std::endl;
896  some_file << "0.0 1.3 1.2" << std::endl;
897  some_file << "1.3 1.3 1.2" << std::endl;
898 
899  some_file.close();
900 
901 
902  // Doc error
903  //----------
904  double error,norm;
905 
906  sprintf(filename,"%s/error%i_on_proc%i.dat",Doc_info.directory().c_str(),
907  Doc_info.number(),this->communicator_pt()->my_rank());
908  some_file.open(filename);
909  Bulk_mesh_pt->compute_error(some_file,
911  time_pt()->time(),
912  error,norm);
913  some_file.close();
914 
915 
916  // Collect error and norm from all processors
917  Vector<double> local_stuff(2);
918  local_stuff[0]=error;
919  local_stuff[1]=norm;
920  Vector<double> global_stuff(2);
921 #ifdef OOMPH_HAS_MPI
922  MPI_Allreduce(&local_stuff[0],&global_stuff[0],2,MPI_DOUBLE,MPI_SUM,
923  communicator_pt()->mpi_comm());
924 #else
925  global_stuff[0]=error;
926  global_stuff[1]=norm;
927 #endif
928 
929  // Doc solution and error
930  //-----------------------
931  oomph_info << "error: " << global_stuff[0] << std::endl;
932  oomph_info << "norm : " << global_stuff[1] << std::endl << std::endl;
933 
934  // Get global error norm (on all processors as it involves comms)
935  double temp_error_norm=global_temporal_error_norm();
936 
937  // Gather refinement stats
938  Vector<int> ref_stats(3);
939  ref_stats[0]=Bulk_mesh_pt->nnon_halo_element();
940  ref_stats[1]=Bulk_mesh_pt->nrefined();
941  ref_stats[2]=Bulk_mesh_pt->nunrefined();
942  Vector<int> global_ref_stats(3);
943  MPI_Allreduce(&ref_stats[0],&global_ref_stats[0],3,MPI_INT,MPI_SUM,
944  communicator_pt()->mpi_comm());
945 
946  // Write trace file on root only
947  if (Communicator_pt->my_rank()==0)
948  {
949  Vector<double> x(2,0.0);
950  double u_exact;
952  Vector<double > xi_wall(1);
953  Vector<double > r_wall(2);
954  xi_wall[0]=0.0;
955  Boundary_pt->position(xi_wall,r_wall);
956  Trace_file << time_pt()->time()
957  << " " << u_exact
958  << " " << r_wall[0]
959  << " " << time_pt()->dt()
960  << " " << temp_error_norm
962  time())
963  << " " << global_ref_stats[0]
964  << " " << global_ref_stats[1]
965  << " " << global_ref_stats[2]
966  << " " << error
967  << " " << norm << std::endl;
968 
969  // Plot wall posn
970  //---------------
971  sprintf(filename,"%s/Wall%i.dat",Doc_info.directory().c_str(),
972  Doc_info.number());
973  some_file.open(filename);
974 
975  unsigned nplot=100;
976  for (unsigned iplot=0;iplot<nplot;iplot++)
977  {
978  xi_wall[0]=0.5*Pi*double(iplot)/double(nplot-1);
979  Boundary_pt->position(xi_wall,r_wall);
980  some_file << r_wall[0] << " " << r_wall[1] << std::endl;
981  }
982  some_file.close();
983  }
984 
985 
986  // Write restart file
987  sprintf(filename,"%s/restart%i_on_proc%i.dat",Doc_info.directory().c_str(),
988  Doc_info.number(),this->communicator_pt()->my_rank());
989  some_file.open(filename);
990  some_file.precision(20);
991  dump_it(some_file);
992  some_file.close();
993 
994  // Increment number of doc
995  Doc_info.number()++;
996 
997 } // end of doc_solution
double global_temporal_error_norm()
Global error norm for adaptive time-stepping.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:1063
int my_rank() const
my rank
Definition: communicator.h:176
unsigned long ndof() const
Return the number of dofs.
Definition: problem.h:1674
OomphCommunicator * communicator_pt()
access function to the oomph-lib communicator
Definition: problem.h:1246
OomphCommunicator * Communicator_pt
The communicator for this problem.
Definition: problem.h:1242
double & dt(const unsigned &t=0)
Definition: timesteppers.h:136
double step_position(const double &time)
Position of step (x-axis intercept)
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:150

References oomph::DocInfo::directory(), GlobalParameters::Doc_info, calibrate::error, MergeRestartFiles::filename, GlobalParameters::get_exact_u(), oomph::DocInfo::number(), oomph::oomph_info, BiharmonicTestFunctions2::Pi, GlobalParameters::step_position(), oomph::Problem_Parameter::Trace_file, and plotDoE::x.

◆ dump_it() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::dump_it ( ofstream &  dump_file)

Dump problem data to allow for later restart.

Dump the solution to disk.

1147 {
1148 
1149  // Write step number
1150  dump_file << Doc_info.number() << " # step number" << std::endl;
1151 
1152  // Write suggested next timestep
1153  dump_file << Next_dt << " # suggested next timestep" << std::endl;
1154 
1155  // Dump the refinement pattern and the generic problem data
1156  Problem::dump(dump_file);
1157 
1158 } // end of dump_it
double Next_dt
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:337

References GlobalParameters::Doc_info, and oomph::DocInfo::number().

◆ dump_it() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::dump_it ( ofstream &  dump_file)

Dump problem data to allow for later restart.

◆ dump_it() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::dump_it ( ofstream &  dump_file)

Dump problem data to allow for later restart.

◆ dump_it() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::dump_it ( ofstream &  dump_file)

Dump problem data to allow for later restart.

◆ dump_it() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::dump_it ( ofstream &  dump_file)

Dump problem data to allow for later restart.

◆ epsilon_t() [1/2]

template<class ELEMENT >
double& RefineableUnsteadyHeatProblem< ELEMENT >::epsilon_t ( )
inline

Target error for adaptive timestepping.

314 {return Epsilon_t;}

◆ epsilon_t() [2/2]

template<class ELEMENT >
double& RefineableUnsteadyHeatProblem< ELEMENT >::epsilon_t ( )
inline

Target error for adaptive timestepping.

278 {return Epsilon_t;}

◆ generic_actions_after()

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::generic_actions_after
private

Actions after adapt/distribute: Rebuild the mesh of prescribed flux elements

621 {
622 
623  // Create prescribed-flux elements from all elements that are
624  // adjacent to boundary 0 and add them to surface mesh
626 
627  // Rebuild the global mesh
629 
630  // Loop over the flux elements to pass pointer to prescribed flux function
631  unsigned n_element=Surface_mesh_pt->nelement();
632  for(unsigned e=0;e<n_element;e++)
633  {
634  // Upcast from GeneralisedElement to UnsteadyHeat flux element
636  dynamic_cast<UnsteadyHeatFluxElement<ELEMENT>*>(
638 
639  // Set the pointer to the prescribed flux function
640  el_pt->flux_fct_pt() =
642  }
643 
644 } // end of generic_actions_after
void rebuild_global_mesh()
Definition: problem.cc:1533

References e(), oomph::UnsteadyHeatFluxElement< ELEMENT >::flux_fct_pt(), and GlobalParameters::prescribed_flux_on_fixed_y_boundary().

◆ generic_actions_before()

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::generic_actions_before
private

Actions before adapt/distribute: Wipe the mesh of prescribed flux elements.

Actions before adapt/distribute: Wipe the mesh of prescribed flux elements

604 {
605 
606  // Kill the flux elements and wipe surface mesh
608 
609  // Rebuild the global mesh.
611 
612 } // end of generic_actions_before
void delete_flux_elements(Mesh *const &surface_mesh_pt)
Delete UnsteadyHeat flux elements and wipe the surface mesh.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:1040

◆ global_temporal_error_norm() [1/2]

template<class ELEMENT >
double RefineableUnsteadyHeatProblem< ELEMENT >::global_temporal_error_norm
virtual

Global error norm for adaptive time-stepping.

Global error norm for adaptive timestepping: RMS error, based on difference between predicted and actual value at all nodes.

Reimplemented from oomph::Problem.

1064 {
1065 
1066 #ifdef OOMPH_HAS_MPI
1067 
1068  double global_error = 0.0;
1069 
1070  //Find out how many nodes there are in the problem
1071  unsigned n_node = Bulk_mesh_pt->nnode();
1072 
1073  // Loop over the nodes and calculate the estimated error in the values
1074  // for non-haloes
1075  int count=0;
1076  for(unsigned i=0;i<n_node;i++)
1077  {
1078  Node* nod_pt=Bulk_mesh_pt->node_pt(i);
1079  if (!(nod_pt->is_halo()))
1080  {
1081  // Get error in solution: Difference between predicted and actual
1082  // value for nodal value 0
1083  double error = nod_pt->time_stepper_pt()->
1084  temporal_error_in_value(nod_pt,0);
1085 
1086  //Add the square of the individual error to the global error
1087  global_error += error*error;
1088  count++;
1089  }
1090  }
1091 
1092  // Accumulate
1093  int n_node_local=count;
1094  int n_node_total=0;
1095  MPI_Allreduce(&n_node_local,&n_node_total,1,MPI_INT,MPI_SUM,
1096  communicator_pt()->mpi_comm());
1097 
1098  double global_error_total=0.0;
1099  MPI_Allreduce(&global_error,&global_error_total,1,MPI_DOUBLE,MPI_SUM,
1100  communicator_pt()->mpi_comm());
1101 
1102  // Divide by the number of nodes
1103  global_error_total /= double(n_node_total);
1104 
1105  // Return square root...
1106  oomph_info << "Total error " << n_node_total << " "
1107  << sqrt(global_error_total) << std::endl;
1108  return sqrt(global_error_total);
1109 
1110 #else
1111 
1112  double global_error = 0.0;
1113 
1114  //Find out how many nodes there are in the problem
1115  unsigned n_node = Bulk_mesh_pt->nnode();
1116 
1117  //Loop over the nodes and calculate the errors in the values
1118  for(unsigned i=0;i<n_node;i++)
1119  {
1120  // Get error in solution: Difference between predicted and actual
1121  // value for nodal value 0
1122  double error =
1123  Bulk_mesh_pt->node_pt(i)->time_stepper_pt()->
1124  temporal_error_in_value(Bulk_mesh_pt->node_pt(i),0);
1125 
1126  //Add the square of the individual error to the global error
1127  global_error += error*error;
1128  }
1129 
1130  //Now the global error must be divided by the number of nodes
1131  global_error /= double(n_node);
1132 
1133  //Return the square root of the error
1134  return sqrt(global_error);
1135 
1136 #endif
1137 
1138 
1139 } // end of global_temporal_error_norm
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238

References calibrate::error, i, oomph::oomph_info, sqrt(), and oomph::Data::time_stepper_pt().

◆ global_temporal_error_norm() [2/2]

template<class ELEMENT >
double RefineableUnsteadyHeatProblem< ELEMENT >::global_temporal_error_norm ( )
virtual

Global error norm for adaptive time-stepping.

Reimplemented from oomph::Problem.

◆ next_dt() [1/2]

template<class ELEMENT >
double& RefineableUnsteadyHeatProblem< ELEMENT >::next_dt ( )
inline

Suggestion for next timestep (stored to allow it to be written to (or read from) restart file)

321 {return Next_dt;}

◆ next_dt() [2/2]

template<class ELEMENT >
double& RefineableUnsteadyHeatProblem< ELEMENT >::next_dt ( )
inline

Suggestion for next timestep (stored to allow it to be written to (or read from) restart file)

285 {return Next_dt;}

◆ possibly_disable_ALE()

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::possibly_disable_ALE ( )
inline

Switch off ALE terms.

252  {
253  if (!Use_ALE)
254  {
255  std::cout << "Disabling ALE " << std::endl;
256  }
257  // Loop over the elements
258  unsigned n_element = Bulk_mesh_pt->nelement();
259  for(unsigned i=0;i<n_element;i++)
260  {
261  // Upcast from FiniteElement to the present element
262  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(i));
263 
264  //Disable/enable ALE (just to make sure both versions of the
265  // code to the same amount of setup work...)
266  if (Use_ALE)
267  {
268  el_pt->enable_ALE();
269  }
270  else
271  {
272  el_pt->disable_ALE();
273  }
274  }
275  }

References i.

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

◆ restart() [1/6]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::restart

Restart.

736 {
737 
738  // Pointer to restart file
739  ifstream* restart_file_pt=0;
740 
741  // Open restart file from stem
742  char filename[100];
743  sprintf(filename,"%s_on_proc%i.dat",GlobalParameters::Restart_file.c_str(),
744  this->communicator_pt()->my_rank());
745  restart_file_pt= new ifstream(filename,ios_base::in);
746  if (restart_file_pt!=0)
747  {
748  oomph_info << "Have opened " << filename <<
749  " for restart. " << std::endl;
750  }
751  else
752  {
753  std::ostringstream error_stream;
754  error_stream
755  << "ERROR while trying to open " << filename
756  << " for restart." << std::endl;
757 
758  throw OomphLibError(
759  error_stream.str(),
762  }
763 
764 
765  // Read restart data:
766  //-------------------
767  if (restart_file_pt!=0)
768  {
769  // Read the data from restart file and find out if the restart file
770  // was from an unsteady run
771  restart(*restart_file_pt);
772  }
773 
774 } // end of restart
void restart()
Restart.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:735
Definition: oomph_definitions.h:222
std::string Restart_file
Name of restart file.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:132
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References MergeRestartFiles::filename, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, and GlobalParameters::Restart_file.

◆ restart() [2/6]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::restart ( ifstream &  restart_file)

Read problem data for restart.

Read solution from disk.

1165 {
1166 
1167  double local_next_dt=0.0;
1168  unsigned local_doc_info_number=0;
1169 
1170  if (restart_file.is_open())
1171  {
1172  oomph_info <<"restart file exists"<<endl;
1173 
1174  // Read line up to termination sign
1175  string input_string;
1176  getline(restart_file,input_string,'#');
1177 
1178  // Ignore rest of line
1179  restart_file.ignore(80,'\n');
1180 
1181  // Read in step number
1182  local_doc_info_number=unsigned(atof(input_string.c_str()));
1183 
1184  // Read line up to termination sign
1185  getline(restart_file,input_string,'#');
1186 
1187  // Ignore rest of line
1188  restart_file.ignore(80,'\n');
1189 
1190  // Read suggested next timestep
1191  local_next_dt=double(atof(input_string.c_str()));
1192 
1193  }
1194  else
1195  {
1196  oomph_info <<"restart file does not exist"<<endl;
1197  }
1198 
1199  // Get next dt
1200  double next_dt=0.0;
1201  MPI_Allreduce(&local_next_dt,&next_dt,1,MPI_DOUBLE,MPI_MAX,
1202  communicator_pt()->mpi_comm());
1203  Next_dt=next_dt;
1204  oomph_info << "Next dt: " << Next_dt << std::endl;
1205 
1206  // Get doc info
1207  unsigned doc_info_number=0;
1208  MPI_Allreduce(&local_doc_info_number,&doc_info_number,1,MPI_UNSIGNED,MPI_MAX,
1209  communicator_pt()->mpi_comm());
1210  Doc_info.number()=doc_info_number;
1211  oomph_info << "Restarted doc_info.number(): " << doc_info_number << std::endl;
1212 
1213 
1214  // Refine the mesh and read in the generic problem data
1215  Problem::read(restart_file);
1216 
1217  // Output restarted solution
1218  std::stringstream comment_stream;
1219  comment_stream << "Restarted; Step "
1220  << doc_info().number();
1221  doc_solution(comment_stream.str());
1222 
1223 
1224 } // end of restart
DocInfo & doc_info()
Return DocInfo object.
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:299
void doc_solution()
Doc the solution.
Definition: optimisation/disable_ALE/unsteady_heat/two_d_unsteady_heat_adapt.cc:703
double & next_dt()
Definition: two_d_unsteady_heat_2adapt_load_balance.cc:321
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE std::enable_if_t< PacketLoad, PacketType > read(const TensorMapper &tensorMapper, const StorageIndex &NCIndex, const StorageIndex &CIndex, const StorageIndex &ld)
read, a template function used for loading the data from global memory. This function is used to guar...
Definition: TensorContractionSycl.h:162

References GlobalParameters::Doc_info, oomph::DocInfo::number(), oomph::oomph_info, and Eigen::TensorSycl::internal::read().

◆ restart() [3/6]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::restart ( ifstream &  restart_file)

Read problem data for restart.

◆ restart() [4/6]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::restart ( ifstream &  restart_file)

Read problem data for restart.

◆ restart() [5/6]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::restart ( ifstream &  restart_file)

Read problem data for restart.

◆ restart() [6/6]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::restart ( ifstream &  restart_file)

Read problem data for restart.

◆ set_initial_condition() [1/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::set_initial_condition
virtual

Set initial condition (incl previous timesteps) according to specified function.

Set initial condition: Assign previous and current values from exact solution.

Reimplemented from oomph::Problem.

653 {
654 
655  // Choose initial timestep
656  double dt0=0.005;
657 
658  // Initialise timestep -- also sets the weights for all timesteppers
659  // in the problem.
660  initialise_dt(dt0);
661 
662  // Use this as the first "suggested" timestep
663  Next_dt=dt0;
664 
665  // Backup time in global timestepper
666  double backed_up_time=time_pt()->time();
667 
668  // Past history for velocities must be established for t=time0-deltat, ...
669  // Then provide current values (at t=time0) which will also form
670  // the initial guess for first solve at t=time0+deltat
671 
672  // Vector of exact solution value
673  Vector<double> soln(1);
674  Vector<double> x(2);
675 
676  //Find number of nodes in mesh
677  unsigned num_nod = Bulk_mesh_pt->nnode();
678 
679  // Get continuous times at previous timesteps
680  int nprev_steps=time_stepper_pt()->nprev_values();
681  Vector<double> prev_time(nprev_steps+1);
682  for (int itime=nprev_steps;itime>=0;itime--)
683  {
684  prev_time[itime]=time_pt()->time(unsigned(itime));
685  }
686 
687  // Loop over current & previous timesteps (in outer loop because
688  // the mesh also moves!)
689  for (int itime=nprev_steps;itime>=0;itime--)
690  {
691  double time=prev_time[itime];
692 
693  // Set global time (because this is how the geometric object refers
694  // to continous time
695  time_pt()->time()=time;
696 
697  oomph_info << "setting IC at time =" << time << std::endl;
698 
699  // Update the mesh for this value of the continuous time
700  // (The wall object reads the continous time from the same
701  // global time object)
702  Bulk_mesh_pt->node_update();
703 
704  // Loop over the nodes to set initial guess everywhere
705  for (unsigned jnod=0;jnod<num_nod;jnod++)
706  {
707  // Get nodal coordinates
708  x[0]=Bulk_mesh_pt->node_pt(jnod)->x(0);
709  x[1]=Bulk_mesh_pt->node_pt(jnod)->x(1);
710 
711  // Get intial solution
713 
714  // Assign solution
715  Bulk_mesh_pt->node_pt(jnod)->set_value(itime,0,soln[0]);
716 
717  // Loop over coordinate directions
718  for (unsigned i=0;i<2;i++)
719  {
720  Bulk_mesh_pt->node_pt(jnod)->x(itime,i)=x[i];
721  }
722  }
723  }
724 
725  // Reset backed up time for global timestepper
726  time_pt()->time()=backed_up_time;
727 
728 } // end of set_initial_condition
void initialise_dt(const double &dt)
Definition: problem.cc:13231
virtual unsigned nprev_values() const =0
Number of previous values available: 0 for static, 1 for BDF<1>,...

References GlobalParameters::get_exact_u(), i, oomph::oomph_info, and plotDoE::x.

◆ set_initial_condition() [2/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::set_initial_condition ( )
virtual

Set initial condition (incl previous timesteps) according to specified function. Note that his overloads the virtual function in the Problem base class and is therefore executed automatically to re-assign the initial conditions during the spatially adaptive solution at the first timestep.

Reimplemented from oomph::Problem.

◆ set_initial_condition() [3/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::set_initial_condition ( )
virtual

Set initial condition (incl previous timesteps) according to specified function.

Reimplemented from oomph::Problem.

◆ set_initial_condition() [4/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::set_initial_condition ( )
virtual

Set initial condition (incl previous timesteps) according to specified function. Note that his overloads the virtual function in the Problem base class and is therefore executed automatically to re-assign the initial conditions during the spatially adaptive solution at the first timestep.

Reimplemented from oomph::Problem.

◆ set_initial_condition() [5/5]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::set_initial_condition ( )
virtual

Set initial condition (incl previous timesteps) according to specified function. Note that his overloads the virtual function in the Problem base class and is therefore executed automatically to re-assign the initial conditions during the spatially adaptive solution at the first timestep.

Reimplemented from oomph::Problem.

◆ write_trace_file_header() [1/2]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::write_trace_file_header

Write header for trace file.

Write the tecplot header for the trace file.

784 {
785  // Open trace file on root only
786  if (communicator_pt()->my_rank()==0)
787  {
788  char filename[100];
789  sprintf(filename,"%s/trace.dat",Doc_info.directory().c_str());
790  Trace_file.open(filename);
791 
792  Trace_file << "VARIABLES=\"time t\","
793  << "\"u<SUB>FE</SUB>\","
794  << "\"u<SUB>exact</SUB>\","
795  << "\"A\","
796  << "\"dt\","
797  << "\"global temporal error norm\","
798  << "\"X<SUB>step</SUB>\","
799  << "\"N<SUB>element</SUB>\","
800  << "\"N<SUB>refined</SUB>\","
801  << "\"N<SUB>unrefined</SUB>\","
802  << "\"norm of error\","
803  << "\"norm of solution\""
804  << std::endl;
805  Trace_file << "ZONE T=\"" << time_stepper_pt()->type()
807  << "; initial dt="
808  << time_pt()->dt() << "; ";
809  if (time_stepper_pt()->adaptive_flag())
810  {
811  Trace_file << "adaptive; target error " << Epsilon_t
812  << " \"" << std::endl;
813  }
814  else
815  {
816  Trace_file << "(fixed timestep)\"" << std::endl;
817  }
818  }
819 
820 } // end of write_trace_file_header
std::string type() const
Definition: timesteppers.h:490

References oomph::DocInfo::directory(), GlobalParameters::Doc_info, MergeRestartFiles::filename, and oomph::Problem_Parameter::Trace_file.

◆ write_trace_file_header() [2/2]

template<class ELEMENT >
void RefineableUnsteadyHeatProblem< ELEMENT >::write_trace_file_header ( )

Write header for trace file.

Member Data Documentation

◆ Boundary_pt

template<class ELEMENT >
GeomObject * RefineableUnsteadyHeatProblem< ELEMENT >::Boundary_pt
private

Pointer to GeomObject that specifies the domain bondary.

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

◆ Bulk_mesh_pt

template<class ELEMENT >
RefineableQuarterCircleSectorMesh< ELEMENT > * RefineableUnsteadyHeatProblem< ELEMENT >::Bulk_mesh_pt
private

◆ Doc_info

template<class ELEMENT >
DocInfo RefineableUnsteadyHeatProblem< ELEMENT >::Doc_info
private

◆ Doc_node_pt

template<class ELEMENT >
Node * RefineableUnsteadyHeatProblem< ELEMENT >::Doc_node_pt
private

Pointer to central node (exists at all refinement levels) for doc.

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

◆ Epsilon_t

template<class ELEMENT >
double RefineableUnsteadyHeatProblem< ELEMENT >::Epsilon_t
private

Target error for adaptive timestepping.

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

◆ Next_dt

template<class ELEMENT >
double RefineableUnsteadyHeatProblem< ELEMENT >::Next_dt
private

Suggestion for next timestep (stored to allow it to be written to (or read from) restart file)

◆ Source_fct_pt

template<class ELEMENT >
UnsteadyHeatEquations< 2 >::UnsteadyHeatSourceFctPt RefineableUnsteadyHeatProblem< ELEMENT >::Source_fct_pt
private

◆ Surface_mesh_pt

template<class ELEMENT >
Mesh * RefineableUnsteadyHeatProblem< ELEMENT >::Surface_mesh_pt
private

◆ Trace_file

template<class ELEMENT >
ofstream RefineableUnsteadyHeatProblem< ELEMENT >::Trace_file
private

◆ Use_ALE

template<class ELEMENT >
bool RefineableUnsteadyHeatProblem< ELEMENT >::Use_ALE
private

Flag to use ALE.


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