UnstructuredPoissonProblem< ELEMENT > Class Template Reference

Class definition. More...

+ Inheritance diagram for UnstructuredPoissonProblem< ELEMENT >:

Public Member Functions

 UnstructuredPoissonProblem ()
 Constructor. More...
 
 ~UnstructuredPoissonProblem ()
 Destructor. More...
 
void actions_before_adapt ()
 Actions before adapt. Empty. More...
 
void actions_after_adapt ()
 
void actions_after_newton_solve ()
 Update after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve: Re-apply boundary conditons. More...
 
void doc_solution (const std::string &comment="")
 Doc the solution. More...
 
 UnstructuredPoissonProblem ()
 Constructor. More...
 
 ~UnstructuredPoissonProblem ()
 Destructor. More...
 
void actions_after_newton_solve ()
 Update after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve: Re-apply boundary conditons. More...
 
void doc_solution (const std::string &comment="")
 Doc the solution. More...
 
 UnstructuredPoissonProblem ()
 Constructor. More...
 
 ~UnstructuredPoissonProblem ()
 Destructor. More...
 
void actions_before_adapt ()
 Actions before adapt. Empty. More...
 
void actions_after_adapt ()
 
void actions_after_newton_solve ()
 Update after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve: Re-apply boundary conditons. More...
 
void doc_solution (const std::string &comment="")
 Doc the solution. 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 apply_boundary_conditions ()
 Helper function to apply boundary conditions. More...
 
void complete_problem_setup ()
 
void apply_boundary_conditions ()
 Helper function to apply boundary conditions. More...
 
void complete_problem_setup ()
 
void apply_boundary_conditions ()
 Helper function to apply boundary conditions. More...
 
void complete_problem_setup ()
 

Private Attributes

DocInfo Doc_info
 Doc info object for labeling output. More...
 
RefineableTriangleMesh< ELEMENT > * My_mesh_pt
 Pointers to specific mesh. More...
 
ofstream Trace_file
 Trace file to document norm of solution. More...
 
TriangleMesh< ELEMENT > * My_mesh_pt
 Pointers to specific mesh. More...
 

Additional Inherited Members

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

Detailed Description

template<class ELEMENT>
class UnstructuredPoissonProblem< ELEMENT >

Class definition.

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

Constructor & Destructor Documentation

◆ UnstructuredPoissonProblem() [1/3]

template<class ELEMENT >
UnstructuredPoissonProblem< ELEMENT >::UnstructuredPoissonProblem

Constructor.

155 {
156  // Intrinsic coordinate along GeomObject
157  Vector<double> zeta(1);
158 
159  // Position vector on GeomObject
160  Vector<double> posn(2);
161 
162  // Ellipse defining the outer boundary
163  double x_center = 0.0;
164  double y_center = 0.0;
165  double A = 1.0;
166  double B = 1.0;
167  Ellipse * outer_boundary_ellipse_pt = new Ellipse(A,B);
168 
169  // Pointer to the closed curve that defines the outer boundary
170  TriangleMeshClosedCurve* closed_curve_pt=0;
171 
172  // Build outer boundary as Polygon?
173  //---------------------------------
174  bool polygon_for_outer_boundary=false;
175 #ifdef OUTER_POLYGON
176  polygon_for_outer_boundary=true;
177 #endif
178  if (polygon_for_outer_boundary)
179  {
180  // Number of segments that make up the boundary
181  unsigned n_seg = 5;
182  double unit_zeta = 0.5*MathematicalConstants::Pi/double(n_seg);
183 
184  // The boundary is bounded by two distinct boundaries, each
185  // represented by its own polyline
186  Vector<TriangleMeshCurveSection*> boundary_polyline_pt(2);
187 
188  // Vertex coordinates on boundary
189  Vector<Vector<double> > bound_coords(n_seg+1);
190 
191  // First part of the boundary
192  //---------------------------
193  for(unsigned ipoint=0; ipoint<n_seg+1;ipoint++)
194  {
195  // Resize the vector
196  bound_coords[ipoint].resize(2);
197 
198  // Get the coordinates
199  zeta[0]=unit_zeta*double(ipoint);
200  outer_boundary_ellipse_pt->position(zeta,posn);
201  bound_coords[ipoint][0]=posn[0]+x_center;
202  bound_coords[ipoint][1]=posn[1]+y_center;
203  }
204 
205  // Build the 1st boundary polyline
206  unsigned boundary_id=0;
207  boundary_polyline_pt[0]=new TriangleMeshPolyLine(bound_coords,boundary_id);
208 
209  // Second part of the boundary
210  //----------------------------
211  unit_zeta*=3.0;
212  for(unsigned ipoint=0; ipoint<n_seg+1;ipoint++)
213  {
214  // Resize the vector
215  bound_coords[ipoint].resize(2);
216 
217  // Get the coordinates
218  zeta[0]=(unit_zeta*double(ipoint))+0.5*MathematicalConstants::Pi;
219  outer_boundary_ellipse_pt->position(zeta,posn);
220  bound_coords[ipoint][0]=posn[0]+x_center;
221  bound_coords[ipoint][1]=posn[1]+y_center;
222  }
223 
224  // Build the 2nd boundary polyline
225  boundary_id=1;
226  boundary_polyline_pt[1]=new TriangleMeshPolyLine(bound_coords,boundary_id);
227 
228 
229  // Create the triangle mesh polygon for outer boundary
230  //----------------------------------------------------
231  TriangleMeshPolygon *outer_polygon =
232  new TriangleMeshPolygon(boundary_polyline_pt);
233 
234  // Enable redistribution of polylines
235  outer_polygon->
236  enable_redistribution_of_segments_between_polylines();
237 
238  // Set the pointer
239  closed_curve_pt = outer_polygon;
240 
241  }
242  // Build outer boundary as curvilinear
243  //------------------------------------
244  else
245  {
246 
247  // Provide storage for pointers to the two parts of the curvilinear boundary
248  Vector<TriangleMeshCurveSection*> outer_curvilinear_boundary_pt(2);
249 
250  // First bit
251  //----------
252  double zeta_start=0.0;
253  double zeta_end=MathematicalConstants::Pi;
254  unsigned nsegment=5;
255  unsigned boundary_id=0;
256  outer_curvilinear_boundary_pt[0]=new TriangleMeshCurviLine(
257  outer_boundary_ellipse_pt,zeta_start,zeta_end,nsegment,boundary_id);
258 
259  // Second bit
260  //-----------
261  zeta_start=MathematicalConstants::Pi;
262  zeta_end=2.0*MathematicalConstants::Pi;
263  nsegment=8;
264  boundary_id=1;
265  outer_curvilinear_boundary_pt[1]=new TriangleMeshCurviLine(
266  outer_boundary_ellipse_pt,zeta_start,zeta_end,nsegment,boundary_id);
267 
268  // Combine to curvilinear boundary and define the
269  //--------------------------------
270  // outer boundary
271  //--------------------------------
272  closed_curve_pt=
273  new TriangleMeshClosedCurve(outer_curvilinear_boundary_pt);
274 
275  }
276 
277 
278  // Now build the holes
279  //====================
281 
282  // Build polygonal hole
283  //=====================
284 
285  // Build first hole: A circle
286  x_center = 0.0;
287  y_center = 0.5;
288  A = 0.1;
289  B = 0.1;
290  Ellipse* polygon_ellipse_pt=new Ellipse(A,B);
291 
292  // Number of segments defining upper and lower half of the hole
293  unsigned n_seg = 6;
294  double unit_zeta = MathematicalConstants::Pi/double(n_seg);
295 
296  // This hole is bounded by two distinct boundaries, each
297  // represented by its own polyline
298  Vector<TriangleMeshCurveSection*> hole_polyline_pt(2);
299 
300 
301  // First boundary of polygonal hole
302  //---------------------------------
303 
304  // Vertex coordinates
305  Vector<Vector<double> > bound_hole(n_seg+1);
306  for(unsigned ipoint=0; ipoint<n_seg+1;ipoint++)
307  {
308  // Resize the vector
309  bound_hole[ipoint].resize(2);
310 
311  // Get the coordinates
312  zeta[0]=unit_zeta*double(ipoint);
313  polygon_ellipse_pt->position(zeta,posn);
314  bound_hole[ipoint][0]=posn[0]+x_center;
315  bound_hole[ipoint][1]=posn[1]+y_center;
316  }
317 
318  // Specify the hole boundary id
319  unsigned boundary_id=2;
320 
321  // Build the 1st hole polyline
322  hole_polyline_pt[0] = new TriangleMeshPolyLine(bound_hole,boundary_id);
323 
324 
325  // Second boundary of polygonal hole
326  //----------------------------------
327  for(unsigned ipoint=0; ipoint<n_seg+1;ipoint++)
328  {
329  // Resize the vector
330  bound_hole[ipoint].resize(2);
331 
332  // Get the coordinates
333  zeta[0]=(unit_zeta*double(ipoint))+MathematicalConstants::Pi;
334  polygon_ellipse_pt->position(zeta,posn);
335  bound_hole[ipoint][0]=posn[0]+x_center;
336  bound_hole[ipoint][1]=posn[1]+y_center;
337  }
338 
339  // Specify the hole boundary id
340  boundary_id=3;
341 
342  // Build the 2nd hole polyline
343  hole_polyline_pt[1] = new TriangleMeshPolyLine(bound_hole,boundary_id);
344 
345 
346  // Build the polygonal hole
347  //-------------------------
348 
349  // Inner hole center coordinates
350  Vector<double> hole_center(2);
351  hole_center[0]=x_center;
352  hole_center[1]=y_center;
353 
354  hole_pt[0] = new TriangleMeshPolygon(hole_polyline_pt, hole_center);
355 
356 
357  // Build curvilinear hole
358  //======================
359 
360  // Build second hole: Another ellipse
361  A = 0.2;
362  B = 0.1;
363  Ellipse* ellipse_pt=new Ellipse(A,B);
364 
365  // Build the two parts of the curvilinear boundary
366  Vector<TriangleMeshCurveSection*> curvilinear_boundary_pt(2);
367 
368 
369  // First part of curvilinear boundary
370  //-----------------------------------
371  double zeta_start=0.0;
372  double zeta_end=MathematicalConstants::Pi;
373  unsigned nsegment=10;
374  boundary_id=4;
375  curvilinear_boundary_pt[0]=new TriangleMeshCurviLine(
376  ellipse_pt,zeta_start,zeta_end,
377  nsegment,boundary_id);
378 
379  // Second part of curvilinear boundary
380  //-------------------------------------
381  zeta_start=MathematicalConstants::Pi;
382  zeta_end=2.0*MathematicalConstants::Pi;
383  nsegment=15;
384  boundary_id=5;
385  curvilinear_boundary_pt[1]=new TriangleMeshCurviLine(
386  ellipse_pt,zeta_start,zeta_end,
387  nsegment,boundary_id);
388 
389 
390  // Combine to hole
391  //----------------
392  Vector<double> hole_coords(2);
393  hole_coords[0]=0.0;
394  hole_coords[1]=0.0;
395  Vector<TriangleMeshClosedCurve*> curvilinear_hole_pt(1);
396  hole_pt[1]=
397  new TriangleMeshClosedCurve(curvilinear_boundary_pt,
398  hole_coords);
399 
400  // Uncomment this as an exercise to observe how a
401  // layer of fine elements get left behind near the boundary
402  // once the tanh step has swept past:
403 
404  // closed_curve_pt->disable_polyline_refinement();
405  // closed_curve_pt->disable_polyline_unrefinement();
406 
407  // Now build the mesh
408  //===================
409 
410  // Use the TriangleMeshParameters object for helping on the manage of the
411  // TriangleMesh parameters
412  TriangleMeshParameters triangle_mesh_parameters(closed_curve_pt);
413 
414  // Specify the closed curve using the TriangleMeshParameters object
415  triangle_mesh_parameters.internal_closed_curve_pt() = hole_pt;
416 
417  // Specify the maximum area element
418  double uniform_element_area=0.2;
419  triangle_mesh_parameters.element_area() = uniform_element_area;
420 
421  // Create the mesh
422  My_mesh_pt=new
423  RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
424 
425  // Store as the problem's one and only mesh
426  Problem::mesh_pt()=My_mesh_pt;
427 
428  // Set error estimator for bulk mesh
430  My_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
431 
432  // Set element size limits
433  My_mesh_pt->max_element_size()=0.2;
434  My_mesh_pt->min_element_size()=0.002;
435 
436  // Set boundary condition and complete the build of all elements
438 
439  // Open trace file
440  char filename[100];
441  sprintf(filename,"RESLT/trace.dat");
442  Trace_file.open(filename);
443 
444  // Setup equation numbering scheme
445  oomph_info <<"Number of equations: "
446  << this->assign_eqn_numbers() << std::endl;
447 
448 } // end_of_constructor
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
RefineableTriangleMesh< ELEMENT > * My_mesh_pt
Pointers to specific mesh.
Definition: mesh_from_inline_triangle.cc:139
void complete_problem_setup()
Definition: mesh_from_inline_triangle.cc:458
ofstream Trace_file
Trace file to document norm of solution.
Definition: mesh_from_inline_triangle.cc:142
Definition: geom_objects.h:644
void position(const Vector< double > &zeta, Vector< double > &r) const
Position Vector at Lagrangian coordinate zeta.
Definition: geom_objects.h:745
Definition: matrices.h:74
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
Unstructured refineable Triangle Mesh.
Definition: triangle_mesh.template.h:2249
Base class defining a closed curve for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1339
Definition: unstructured_two_d_mesh_geometry_base.h:662
Definition: triangle_mesh.template.h:94
Class defining a polyline for use in Triangle Mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:868
Class defining a closed polygon for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1451
Definition: error_estimator.h:266
double Pi
Definition: two_d_biharmonic.cc:235
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152
string filename
Definition: MergeRestartFiles.py:39
Z2ErrorEstimator * error_estimator_pt
Definition: MortaringCantileverCompareToNonMortaring.cpp:190
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::TriangleMeshParameters::element_area(), MeshRefinement::error_estimator_pt, MergeRestartFiles::filename, oomph::TriangleMeshParameters::internal_closed_curve_pt(), oomph::oomph_info, BiharmonicTestFunctions2::Pi, oomph::Ellipse::position(), oomph::Problem_Parameter::Trace_file, and Eigen::zeta().

◆ ~UnstructuredPoissonProblem() [1/3]

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

Destructor.

99 {};

◆ UnstructuredPoissonProblem() [2/3]

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

Constructor.

◆ ~UnstructuredPoissonProblem() [2/3]

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

Destructor.

99 {};

◆ UnstructuredPoissonProblem() [3/3]

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

Constructor.

◆ ~UnstructuredPoissonProblem() [3/3]

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

Destructor.

100 {};

Member Function Documentation

◆ actions_after_adapt() [1/2]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::actions_after_adapt ( )
inlinevirtual

Actions after adapt: Setup the problem again – remember that the mesh has been completely rebuilt and its element's don't have any pointers to source fcts etc. yet

Reimplemented from oomph::Problem.

109  {
111  }

◆ actions_after_adapt() [2/2]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::actions_after_adapt ( )
inlinevirtual

Actions after adapt: Setup the problem again – remember that the mesh has been completely rebuilt and its element's don't have any pointers to source fcts etc. yet

Reimplemented from oomph::Problem.

110  {
112  }

◆ actions_after_newton_solve() [1/3]

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

Update after solve (empty)

Reimplemented from oomph::Problem.

114 {}

◆ actions_after_newton_solve() [2/3]

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

Update after solve (empty)

Reimplemented from oomph::Problem.

102 {}

◆ actions_after_newton_solve() [3/3]

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

Update after solve (empty)

Reimplemented from oomph::Problem.

115 {}

◆ actions_before_adapt() [1/2]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::actions_before_adapt ( )
inlinevirtual

Actions before adapt. Empty.

Reimplemented from oomph::Problem.

102 {}

◆ actions_before_adapt() [2/2]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::actions_before_adapt ( )
inlinevirtual

Actions before adapt. Empty.

Reimplemented from oomph::Problem.

103 {}

◆ actions_before_newton_solve() [1/3]

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

Update the problem specs before solve: Re-apply boundary conditons.

Reimplemented from oomph::Problem.

118  {
120  }
void apply_boundary_conditions()
Helper function to apply boundary conditions.
Definition: mesh_from_inline_triangle.cc:502

◆ actions_before_newton_solve() [2/3]

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

Update the problem specs before solve: Re-apply boundary conditons.

Reimplemented from oomph::Problem.

106  {
108  }

◆ actions_before_newton_solve() [3/3]

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

Update the problem specs before solve: Re-apply boundary conditons.

Reimplemented from oomph::Problem.

119  {
121  }

◆ apply_boundary_conditions() [1/3]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::apply_boundary_conditions
private

Helper function to apply boundary conditions.

503 {
504 
505  // Loop over all boundary nodes
506  unsigned nbound=this->My_mesh_pt->nboundary();
507  for(unsigned ibound=0;ibound<nbound;ibound++)
508  {
509  unsigned num_nod=this->My_mesh_pt->nboundary_node(ibound);
510  for (unsigned inod=0;inod<num_nod;inod++)
511  {
512  // Get node
513  Node* nod_pt=this->My_mesh_pt->boundary_node_pt(ibound,inod);
514 
515  // Extract nodal coordinates from node:
516  Vector<double> x(2);
517  x[0]=nod_pt->x(0);
518  x[1]=nod_pt->x(1);
519 
520  // Compute the value of the exact solution at the nodal point
521  Vector<double> u(1);
523 
524  // Assign the value to the one (and only) nodal value at this node
525  nod_pt->set_value(0,u[0]);
526  }
527  }
528 
529 } // end set bc
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
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
Definition: extrude_with_macro_element_representation.cc:206
list x
Definition: plotDoE.py:28

References TanhSolnForPoisson::get_exact_u(), oomph::Data::set_value(), plotDoE::x, and oomph::Node::x().

◆ apply_boundary_conditions() [2/3]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::apply_boundary_conditions ( )
private

Helper function to apply boundary conditions.

◆ apply_boundary_conditions() [3/3]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::apply_boundary_conditions ( )
private

Helper function to apply boundary conditions.

◆ complete_problem_setup() [1/3]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::complete_problem_setup
private

Helper function to (re-)set boundary condition and complete the build of all elements

Set boundary condition exactly, and complete the build of all elements

459 {
460 
461  // Set the boundary conditions for problem: All nodes are
462  // free by default -- just pin the ones that have Dirichlet conditions
463  // here.
464  unsigned nbound=My_mesh_pt->nboundary();
465  for(unsigned ibound=0;ibound<nbound;ibound++)
466  {
467  unsigned num_nod=My_mesh_pt->nboundary_node(ibound);
468  for (unsigned inod=0;inod<num_nod;inod++)
469  {
470  // Get node
471  Node* nod_pt=My_mesh_pt->boundary_node_pt(ibound,inod);
472 
473  // Pin one-and-only unknown value
474  nod_pt->pin(0);
475  }
476  } // end loop over boundaries
477 
478 
479  // Complete the build of all elements so they are fully functional
480  unsigned n_element = My_mesh_pt->nelement();
481  for(unsigned e=0;e<n_element;e++)
482  {
483  // Upcast from GeneralisedElement to the present element
484  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(My_mesh_pt->element_pt(e));
485 
486  //Set the source function pointer
487  el_pt->source_fct_pt() = &TanhSolnForPoisson::get_source;
488  }
489 
490  // Re-apply Dirichlet boundary conditions (projection ignores
491  // boundary conditions!)
493 }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void get_source(const Vector< double > &x, double &source)
Source function to make it an exact solution.
Definition: extrude_with_macro_element_representation.cc:224

References e(), TanhSolnForPoisson::get_source(), and oomph::Data::pin().

◆ complete_problem_setup() [2/3]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::complete_problem_setup ( )
private

Helper function to (re-)set boundary condition and complete the build of all elements

◆ complete_problem_setup() [3/3]

template<class ELEMENT >
void UnstructuredPoissonProblem< ELEMENT >::complete_problem_setup ( )
private

Helper function to (re-)set boundary condition and complete the build of all elements

◆ doc_solution() [1/3]

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

Doc the solution.

538 {
539  ofstream some_file;
540  char filename[100];
541 
542  // Number of plot points
543  unsigned npts;
544  npts=5;
545 
546  sprintf(filename,"RESLT/soln%i.dat",Doc_info.number());
547  some_file.open(filename);
548  this->My_mesh_pt->output(some_file,npts);
549  some_file << "TEXT X = 22, Y = 92, CS=FRAME T = \""
550  << comment << "\"\n";
551  some_file.close();
552 
553  // Output exact solution
554  //----------------------
555  sprintf(filename,"RESLT/exact_soln%i.dat",Doc_info.number());
556  some_file.open(filename);
557  My_mesh_pt->output_fct(some_file,npts,TanhSolnForPoisson::get_exact_u);
558  some_file.close();
559 
560  // Output boundaries
561  //------------------
562  sprintf(filename,"RESLT/boundaries%i.dat",Doc_info.number());
563  some_file.open(filename);
564  My_mesh_pt->output_boundaries(some_file);
565  some_file.close();
566 
567 
568  // Doc error and return of the square of the L2 error
569  //---------------------------------------------------
570  double error,norm,dummy_error,zero_norm;
571  sprintf(filename,"RESLT/error%i.dat",Doc_info.number());
572  some_file.open(filename);
573  My_mesh_pt->compute_error(some_file,TanhSolnForPoisson::get_exact_u,
574  error,norm);
575 
576  My_mesh_pt->compute_error(some_file,TanhSolnForPoisson::zero,
577  dummy_error,zero_norm);
578  some_file.close();
579 
580  // Doc L2 error and norm of solution
581  oomph_info << "\nNorm of error : " << sqrt(error) << std::endl;
582  oomph_info << "Norm of exact solution: " << sqrt(norm) << std::endl;
583  oomph_info << "Norm of computed solution: " << sqrt(dummy_error) << std::endl;
584  Trace_file << sqrt(norm) << " " << sqrt(dummy_error) << std::endl;
585 
586  // Increment the doc_info number
587  Doc_info.number()++;
588 
589 } // end of doc
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
DocInfo Doc_info
Doc info object for labeling output.
Definition: mesh_from_inline_triangle.cc:129
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
void zero(const Vector< double > &x, Vector< double > &u)
Definition: mesh_from_inline_triangle.cc:72
int error
Definition: calibrate.py:297

References GlobalParameters::Doc_info, calibrate::error, MergeRestartFiles::filename, TanhSolnForPoisson::get_exact_u(), oomph::DocInfo::number(), oomph::oomph_info, sqrt(), oomph::Problem_Parameter::Trace_file, and TanhSolnForPoisson::zero().

◆ doc_solution() [2/3]

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

Doc the solution.

◆ doc_solution() [3/3]

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

Doc the solution.

Member Data Documentation

◆ Doc_info

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

Doc info object for labeling output.

◆ My_mesh_pt [1/2]

template<class ELEMENT >
RefineableTriangleMesh< ELEMENT > * UnstructuredPoissonProblem< ELEMENT >::My_mesh_pt
private

Pointers to specific mesh.

◆ My_mesh_pt [2/2]

template<class ELEMENT >
TriangleMesh<ELEMENT>* UnstructuredPoissonProblem< ELEMENT >::My_mesh_pt
private

Pointers to specific mesh.

◆ Trace_file

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

Trace file to document norm of solution.


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