DarcyProblem< ELEMENT > Class Template Reference

Darcy problem. More...

+ Inheritance diagram for DarcyProblem< ELEMENT >:

Public Member Functions

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

Private Member Functions

void create_pressure_elements ()
 Allocate pressure elements on the top surface. More...
 
void complete_problem_setup ()
 Complete problem setup. More...
 

Private Attributes

TriangleMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
MeshSurface_mesh_pt
 Mesh of boundary condition elements. More...
 
std::ofstream Trace_file
 Trace file for results. 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_after_newton_solve ()
 
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 DarcyProblem< ELEMENT >

Darcy problem.

Constructor & Destructor Documentation

◆ DarcyProblem()

template<class ELEMENT >
DarcyProblem< ELEMENT >::DarcyProblem

Constructor.

280 {
281  // Mesh
282  //-----
283 
284  // Closed curve definiting outer boundary
285  TriangleMeshClosedCurve *outer_boundary_pt=0;
286 
287  // Internal boundary to prevent problem with flux B.C.s
288  Vector<TriangleMeshOpenCurve*> inner_open_boundary_pt;
289 
290  Vector<double> zeta(1);
291  Vector<double> posn(2);
292 
293  Ellipse* outer_boundary_ellipse_pt = new Ellipse(1,1);
294 
295  Vector<TriangleMeshCurveSection*> outer_curvilinear_boundary_pt(2);
296 
297  // First bit
298  double zeta_start = 0.0;
299  double zeta_end = MathematicalConstants::Pi;
300  unsigned nsegment =
301  (unsigned)max(
303  outer_curvilinear_boundary_pt[0] =
304  new TriangleMeshCurviLine(outer_boundary_ellipse_pt, zeta_start,
305  zeta_end, nsegment, 0);
306 
307  // Second bit
308  zeta_start = MathematicalConstants::Pi;
309  zeta_end = 2.0*MathematicalConstants::Pi;
310  nsegment =
311  (unsigned)max(
313  outer_curvilinear_boundary_pt[1] =
314  new TriangleMeshCurviLine(outer_boundary_ellipse_pt, zeta_start,
315  zeta_end, nsegment, 1);
316 
317  outer_boundary_pt =
318  new TriangleMeshClosedCurve(outer_curvilinear_boundary_pt);
319 
320  // Target element area for Triangle
321  double uniform_element_area = TestProblem::Element_area;
322 
323  // Use the TriangleMeshParameters object for gathering all
324  // the necessary arguments for the TriangleMesh object
325  TriangleMeshParameters triangle_mesh_parameters(
326  outer_boundary_pt);
327 
328  // Take the maximum element area
329  triangle_mesh_parameters.element_area() =
330  uniform_element_area;
331 
332 
333 #ifdef ADAPTIVE
334 
335  // Build adaptive "bulk" mesh
336  Bulk_mesh_pt=new RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
337 
338  // Create/set error estimator
339  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
340 
341  // Choose error tolerances
343  {
344  Bulk_mesh_pt->min_permitted_error()=0.0001;
345  Bulk_mesh_pt->max_permitted_error()=0.001;
346  }
347  else
348  {
349  Bulk_mesh_pt->min_permitted_error()=1.0e-8;
350  Bulk_mesh_pt->max_permitted_error()=1.0e-7;
351  }
352 
353  // Actions before projection
354  Bulk_mesh_pt->mesh_update_fct_pt()=&TestProblem::edge_sign_setup<ELEMENT>;
355 
356 #else
357 
358  // Create the buk mesh
359  Bulk_mesh_pt = new TriangleMesh<ELEMENT>(triangle_mesh_parameters);
360 
361 #endif
362 
363  // Create mesh of pressure bc elements
364  Surface_mesh_pt = new Mesh;
366 
367  // Combine into global mesh
371 
372  // Complete problem setup
374 
375  // Assign equation numbers
376  oomph_info << "Number of equations: " << assign_eqn_numbers() << std::endl;
377 
378 
379  // Open trace file
380  char filename[50];
381  sprintf(filename,"%s/trace.dat",TestProblem::Directory.c_str());
382  Trace_file.open(filename);
383 
384 } // end of problem constructor
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Mesh * Surface_mesh_pt
Mesh of boundary condition elements.
Definition: unstructured_two_d_circle.cc:264
void create_pressure_elements()
Allocate pressure elements on the top surface.
Definition: unstructured_two_d_circle.cc:392
void complete_problem_setup()
Complete problem setup.
Definition: unstructured_two_d_circle.cc:426
std::ofstream Trace_file
Trace file for results.
Definition: unstructured_two_d_circle.cc:267
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: unstructured_two_d_circle.cc:259
Definition: geom_objects.h:644
Definition: mesh.h:67
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
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
Definition: triangle_mesh.template.h:424
Definition: error_estimator.h:266
#define max(a, b)
Definition: datatypes.h:23
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
double Element_area
Target element area for Triangle.
Definition: unstructured_two_d_circle.cc:126
bool Use_point_source_solution
Definition: unstructured_two_d_circle.cc:43
std::string Directory
The directory in which the solution is output.
Definition: unstructured_two_d_circle.cc:129
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References TestProblem::Directory, TestProblem::Element_area, oomph::TriangleMeshParameters::element_area(), MergeRestartFiles::filename, max, oomph::oomph_info, oomph::MathematicalConstants::Pi, sqrt(), oomph::Problem_Parameter::Trace_file, TestProblem::Use_point_source_solution, and Eigen::zeta().

Member Function Documentation

◆ actions_before_newton_solve()

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

Actions before Newton solve (empty)

Reimplemented from oomph::Problem.

197 {}

◆ complete_problem_setup()

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

Complete problem setup.

427 {
428  // Setup the signs for the fluxes
429  TestProblem::edge_sign_setup<ELEMENT>(Bulk_mesh_pt);
430 
431  // Loop over all elements
432  unsigned n_element = Bulk_mesh_pt->nelement();
433  for(unsigned e=0;e<n_element;e++)
434  {
435  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
436 
437  // Set source fcts
438  el_pt->source_fct_pt()=TestProblem::source;
439  el_pt->mass_source_fct_pt()=TestProblem::mass_source;
440 
441 #ifdef ADAPTIVE
442 
443  // Pin dofs at vertex nodes (because they're only used for projection)
444  el_pt->pin_superfluous_darcy_dofs();
445 
446 #endif
447 
448  } // end of loop over elements
449 
450 
451  // Apply boundary conditions
452  //--------------------------
453 
454  // Get the nodal index at which values representing edge fluxes
455  // at flux interpolation points are stored
456  ELEMENT *el_pt=dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(0));
457 
458  // Where are the values stored?
459  unsigned n=el_pt->nedge_flux_interpolation_point();
460  Vector<unsigned> q_index(n);
461  for (unsigned j=0;j<n;j++)
462  {
463  q_index[j]=el_pt->q_edge_index(j);
464  }
465 
466  // Coordinate vector
467  Vector<double> x(2);
468 
469  // Impose normal component of flux on lower part of outer boundary
470  unsigned ibound=1;
471  {
472  // Get the number of elements along boundary ibound
473  unsigned n_boundary_element=Bulk_mesh_pt->nboundary_element(ibound);
474 
475  // Loop over the elements along boundary ibound
476  for(unsigned e=0;e<n_boundary_element;e++)
477  {
478  // Upcast the current element to the actual type
479  ELEMENT *el_pt=
480  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->boundary_element_pt(ibound,e));
481 
482  // Loop over the edges
483  for(unsigned edge=0;edge<3;edge++)
484  {
485  // Get pointer to node that stores the edge flux dofs for this edge
486  Node* nod_pt=el_pt->edge_flux_node_pt(edge);
487 
488  // Pin/set values for the flux degrees of freedom
489  if (nod_pt->is_on_boundary(ibound))
490  {
491  for (unsigned j=0;j<n;j++)
492  {
493  nod_pt->pin(q_index[j]);
494  }
495 
496  // Get face index of face associated with edge
497  unsigned f=el_pt->face_index_of_edge(edge);
498 
499  // Build a temporary face element from which we'll extract
500  // the outer unit normal
501  DarcyFaceElement<ELEMENT>* face_el_pt=
502  new DarcyFaceElement<ELEMENT>(el_pt,f);
503 
504  // Loop over the flux interpolation points
505  unsigned n_flux_interpolation_points=
506  el_pt->nedge_flux_interpolation_point();
507  for(unsigned g=0;g<n_flux_interpolation_points;g++)
508  {
509  // Get the global coords of the flux_interpolation point
510  el_pt->edge_flux_interpolation_point_global(edge,g,x);
511 
512  // Get the exact solution
513  Vector<double> exact_soln(4,0.0);
515 
516  // Get unit normal at this flux interpolation point
517  Vector<double> s(1);
518  el_pt->face_local_coordinate_of_flux_interpolation_point(edge,g,s);
519  Vector<double> unit_normal(2);
520  face_el_pt->outer_unit_normal(s,unit_normal);
521 
522 #ifdef PARANOID
523 
524  // Sanity check
525  Vector<double> x_face(2);
526  face_el_pt->interpolated_x(s,x_face);
527  if ((x_face[0]-x[0])*(x_face[0]-x[0])+
528  (x_face[1]-x[1])*(x_face[1]-x[1])>1.0e-3)
529  {
530  std::stringstream error;
531  error << "Discrepancy in coordinate of flux interpolation point\n"
532  << "(computed by bulk and face elements) for edge " << e
533  << " and flux int pt " << g << "\n"
534  << "Face thinks node is at: "
535  << x_face[0] << " " << x_face[1] << "\n"
536  << "Bulk thinks node is at: "
537  << x[0] << " " << x[1] << "\n";
538  throw OomphLibError(
539  error.str(),
542  }
543 #endif
544 
545  // Set the boundary flux
546  nod_pt->set_value(q_index[g],
547  exact_soln[0]*unit_normal[0]+
548  exact_soln[1]*unit_normal[1]);
549 
550  } // End of loop over flux interpolation points
551 
552  // Don't need face element on that edge any more
553  delete face_el_pt;
554 
555  } // End if for edge on required boundary
556  } // End of loop over edges
557 
558  } // End of loop over boundary elements
559  } // End of loop over boundaries
560 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Definition: darcy_face_elements.h:73
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
Definition: nodes.h:906
virtual bool is_on_boundary() const
Definition: nodes.h:1373
Definition: oomph_definitions.h:222
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
void mass_source(const Vector< double > &x, double &m)
Mass source.
Definition: unstructured_two_d_circle.cc:68
void exact_soln(const Vector< double > &x, Vector< double > &soln)
Exact solution: q1,q2,div_q,p.
Definition: unstructured_two_d_circle.cc:98
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
int error
Definition: calibrate.py:297
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References e(), calibrate::error, TestProblem::exact_soln(), f(), oomph::FaceElement::interpolated_x(), oomph::Node::is_on_boundary(), j, TestProblem::mass_source(), n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::FaceElement::outer_unit_normal(), oomph::Data::pin(), s, oomph::Data::set_value(), TestProblem::source(), and plotDoE::x.

◆ create_pressure_elements()

template<class ELEMENT >
void DarcyProblem< ELEMENT >::create_pressure_elements
private

Allocate pressure elements on the top surface.

Creat traction elements.

393 {
394 
395  // Upper part of boundary
396  unsigned bound=0;
397 
398  // How many bulk elements are next to that boundary
399  unsigned n_neigh = Bulk_mesh_pt->nboundary_element(bound);
400 
401  // Now loop over bulk elements and create the face elements
402  for(unsigned n=0;n<n_neigh;n++)
403  {
404  // Create the face element
405  DarcyFaceElement<ELEMENT>* pressure_element_pt
407  (Bulk_mesh_pt->boundary_element_pt(bound,n),
408  Bulk_mesh_pt->face_index_at_boundary(bound,n));
409 
410  // Add to mesh
411  Surface_mesh_pt->add_element_pt(pressure_element_pt);
412 
413  // Set function pointer
414  pressure_element_pt->pressure_fct_pt()=&TestProblem::boundary_pressure;
415 
416  }
417 
418 } // end of assign_traction_elements
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, double &pressure) pressure_fct_pt()
Reference to the pressure function pointer.
Definition: darcy_face_elements.h:148
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617
void boundary_pressure(const double &time, const Vector< double > &x, const Vector< double > &n, double &result)
Pressure around the boundary of the domain.
Definition: unstructured_two_d_circle.cc:82

References TestProblem::boundary_pressure(), n, and oomph::DarcyFaceElement< ELEMENT >::pressure_fct_pt.

◆ doc_shape_functions()

template<class ELEMENT >
void DarcyProblem< ELEMENT >::doc_shape_functions

Doc shape functions.

567 {
568  ofstream some_file;
569  char filename[100];
570 
571  // Coarse solution
572  unsigned npts_coarse=2;
573  sprintf(filename,"%s/coarse_soln.dat",TestProblem::Directory.c_str());
574  some_file.open(filename);
575  Bulk_mesh_pt->output(some_file,npts_coarse);
576  some_file.close();
577 
578  // Backup for all values (double counting plenty but I don't care...)
579  //-------------------------------------------------------------------
580  Vector<double> backed_up_value;
581  unsigned nel=Bulk_mesh_pt->nelement();
582  for (unsigned e=0;e<nel;e++)
583  {
584  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
585 
586  // Backup edge-based q values
587  unsigned n_q_basis_edge = el_pt->nq_basis_edge();
588  for (unsigned j=0;j<n_q_basis_edge;j++)
589  {
590  backed_up_value.push_back(el_pt->q_edge(j));
591  }
592 
593  // Backup internal q values
594  unsigned n_q_basis = el_pt->nq_basis();
595  for (unsigned j=n_q_basis_edge;j<n_q_basis;j++)
596  {
597  backed_up_value.push_back(el_pt->q_internal(j-n_q_basis_edge));
598  }
599 
600  // Backup pressures
601  unsigned n_p_basis = el_pt->np_basis();
602  for (unsigned j=0;j<n_p_basis;j++)
603  {
604  backed_up_value.push_back(el_pt->p_value(j));
605  }
606  }
607 
608 
609 
610  // Loop over all elements to wipe
611  //-------------------------------
612  for (unsigned e=0;e<nel;e++)
613  {
614  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
615 
616  // Kill edge-based q values
617  unsigned n_q_basis_edge = el_pt->nq_basis_edge();
618  for (unsigned j=0;j<n_q_basis_edge;j++)
619  {
620  el_pt->set_q_edge(j,0.0);
621  }
622 
623  // Kill internal q values
624  unsigned n_q_basis = el_pt->nq_basis();
625  for (unsigned j=n_q_basis_edge;j<n_q_basis;j++)
626  {
627  el_pt->set_q_internal(j-n_q_basis_edge,0.0);
628  }
629 
630 
631  // Kill pressures
632  unsigned n_p_basis = el_pt->np_basis();
633  for (unsigned j=0;j<n_p_basis;j++)
634  {
635  el_pt->set_p_value(j,0.0);
636  }
637  }
638 
639 
640  // Choose a specific element for which we'll switch on all the dofs
641  //-----------------------------------------------------------------
642  // one by one
643  //-----------
644  MeshAsGeomObject* mesh_geom_obj_pt=new MeshAsGeomObject(Bulk_mesh_pt);
645  Vector<double> x(2,0.0);
646  Vector<double> s(2,0.0);
647 
648  // Get pointer to GeomObject that contains this point and cast
649  GeomObject* geom_obj_pt=0;
650  mesh_geom_obj_pt->locate_zeta(x,geom_obj_pt,s);
651  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(geom_obj_pt);
652 
653 
654  // Number of plot points
655  unsigned npts=10;
656 
657  // Doc edge-based q values
658  unsigned count_shape_doc=0;
659  unsigned n_q_basis_edge = el_pt->nq_basis_edge();
660  for (unsigned j=0;j<n_q_basis_edge;j++)
661  {
662  // Set this dof to one
663  el_pt->set_q_edge(j,1.0);
664 
665  // Get face index of face associated with this flux basis fct
666  // and erect face element
667  unsigned f=el_pt->face_index_of_q_edge_basis_fct(j);
668  DarcyFaceElement<ELEMENT>* face_el_pt=
669  new DarcyFaceElement<ELEMENT>(el_pt,f);
670 
671  // Get outer unit normal halfway along edge
672  Vector<double> s(1,0.5);
673  Vector<double> unit_normal(2);
674  Vector<double> x(2);
675  face_el_pt->outer_unit_normal(s,unit_normal);
676  face_el_pt->interpolated_x(s,x);
677  delete face_el_pt;
678 
679  // Plot
680  sprintf(filename,"%s/q_shape_fct%i.dat",TestProblem::Directory.c_str(),
681  count_shape_doc++);
682  some_file.open(filename);
683  unsigned nel=Bulk_mesh_pt->nelement();
684  for (unsigned e=0;e<nel;e++)
685  {
686  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e))->
687  output(some_file,npts);
688  }
689  some_file.close();
690 
691 
692  sprintf(filename,"%s/outer_unit_normal%i.dat",
693  TestProblem::Directory.c_str(),j);
694  some_file.open(filename);
695  some_file << x[0] << " "
696  << x[1] << " "
697  << unit_normal[0] << " "
698  << unit_normal[1] << " "
699  << std::endl;
700  some_file.close();
701 
702 
703  sprintf(filename,"%s/flux_interpolation_point%i.dat",
704  TestProblem::Directory.c_str(),j);
705  some_file.open(filename);
706  Vector<double> flux_interpolation_point(2);
707  el_pt->edge_flux_interpolation_point_global(j,flux_interpolation_point);
708  some_file << flux_interpolation_point[0] << " "
709  << flux_interpolation_point[1] << " "
710  << 1.0 << " "
711  << std::endl;
712  some_file.close();
713 
714  sprintf(filename,"%s/q_shape_fct_with_project%i.dat",
715  TestProblem::Directory.c_str(),j);
716  some_file.open(filename);
717  for (unsigned e=0;e<nel;e++)
718  {
719  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e))->
720  output_with_projected_flux(some_file,npts,unit_normal);
721  }
722  some_file.close();
723 
724  // Reset
725  el_pt->set_q_edge(j,0.0);
726  }
727 
728  // Doc internal q values
729  unsigned n_q_basis = el_pt->nq_basis();
730  for (unsigned j=n_q_basis_edge;j<n_q_basis;j++)
731  {
732  el_pt->set_q_internal(j-n_q_basis_edge,1.0);
733  sprintf(filename,"%s/q_shape_fct%i.dat",
734  TestProblem::Directory.c_str(),
735  count_shape_doc++);
736  some_file.open(filename);
737  unsigned nel=Bulk_mesh_pt->nelement();
738  for (unsigned e=0;e<nel;e++)
739  {
740  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e))->
741  output(some_file,npts);
742  }
743  some_file.close();
744  el_pt->set_q_internal(j-n_q_basis_edge,0.0);
745  }
746 
747  // Doc pressures
748  count_shape_doc=0;
749  unsigned n_p_basis = el_pt->np_basis();
750  for (unsigned j=0;j<n_p_basis;j++)
751  {
752  el_pt->set_p_value(j,1.0);
753  sprintf(filename,"%s/p_shape_fct%i.dat",
754  TestProblem::Directory.c_str(),
755  count_shape_doc++);
756  some_file.open(filename);
757  unsigned nel=Bulk_mesh_pt->nelement();
758  for (unsigned e=0;e<nel;e++)
759  {
760  dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e))->
761  output(some_file,npts);
762  }
763  some_file.close();
764  el_pt->set_p_value(j,0.0);
765  }
766 
767 
768 
769  // Restore dofs
770  unsigned count=0;
771  for (unsigned e=0;e<nel;e++)
772  {
773  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
774 
775  // Reset edge-based q values
776  unsigned n_q_basis_edge = el_pt->nq_basis_edge();
777  for (unsigned j=0;j<n_q_basis_edge;j++)
778  {
779  el_pt->set_q_edge(j,backed_up_value[count]);
780  count++;
781  }
782 
783  // Reset internal q values
784  unsigned n_q_basis = el_pt->nq_basis();
785  for (unsigned j=n_q_basis_edge;j<n_q_basis;j++)
786  {
787  el_pt->set_q_internal(j-n_q_basis_edge,backed_up_value[count]);
788  count++;
789  }
790 
791 
792  // Reset pressures
793  unsigned n_p_basis = el_pt->np_basis();
794  for (unsigned j=0;j<n_p_basis;j++)
795  {
796  el_pt->set_p_value(j,backed_up_value[count]);
797  count++;
798  }
799  }
800 }
Definition: geom_objects.h:101
Definition: mesh_as_geometric_object.h:93
void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Definition: mesh_as_geometric_object.h:373
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490

References TestProblem::Directory, e(), f(), MergeRestartFiles::filename, oomph::FaceElement::interpolated_x(), j, oomph::MeshAsGeomObject::locate_zeta(), oomph::FaceElement::outer_unit_normal(), oomph::output(), s, and plotDoE::x.

◆ doc_solution()

template<class ELEMENT >
void DarcyProblem< ELEMENT >::doc_solution ( const unsigned label)

Post-process results, using unsigned to label files.

Write the solution and exact solution to file, and calculate the error.

807 {
808  ofstream some_file;
809  char filename[100];
810 
811 
812  // Output computed solution
813  //-------------------------
814  unsigned npts=5;
815  sprintf(filename,"%s/soln%i.dat",TestProblem::Directory.c_str(),label);
816  some_file.open(filename);
817  Bulk_mesh_pt->output(some_file,npts);
818  some_file.close();
819 
820 
821  // Output coarse solution
822  //-----------------------
823  unsigned npts_coarse=2;
824  sprintf(filename,"%s/coarse_soln%i.dat",TestProblem::Directory.c_str(),label);
825  some_file.open(filename);
826  Bulk_mesh_pt->output(some_file,npts_coarse);
827  some_file.close();
828 
829  // Output exact solution
830  //----------------------
831  sprintf(filename,"%s/exact_soln%i.dat",TestProblem::Directory.c_str(),label);
832  some_file.open(filename);
833  Bulk_mesh_pt->output_fct(some_file,
834  npts,
836  some_file.close();
837 
838  // Output boundary condition elements
839  //-----------------------------------
840  sprintf(filename,"%s/bc_elements%i.dat",TestProblem::Directory.c_str(),label);
841  some_file.open(filename);
842  Surface_mesh_pt->output(some_file,npts);
843  some_file.close();
844 
845  // Doc error
846  //----------
847  Vector<double> norm(2,0.0);
848  Vector<double> error(2,0.0);
849  sprintf(filename,"%s/error%i.dat",TestProblem::Directory.c_str(),label);
850  some_file.open(filename);
851  Bulk_mesh_pt->compute_error(some_file,
853  error,norm);
854  some_file.close();
855 
856  // Doc error norm:
857  cout << std::endl;
858  cout << "Norm of exact q : " << sqrt(norm[0]) << std::endl;
859  cout << "Norm of exact p : " << sqrt(norm[1]) << std::endl << std::endl;
860 
861  cout << "Norm of q error : " << sqrt(error[0]) << std::endl;
862  cout << "Norm of p error : " << sqrt(error[1]) << std::endl;
863  cout << std::endl << std::endl;
864 
865  Trace_file << ndof() << " "
866  << sqrt(error[0]) << " "
867  << sqrt(error[1]) << " "
868  << TestProblem::Element_area << " "
869  << sqrt(norm[0]) << " "
870  << sqrt(norm[1]) << std::endl;
871 }
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
unsigned long ndof() const
Return the number of dofs.
Definition: problem.h:1674
label
Definition: UniformPSDSelfTest.py:24

References TestProblem::Directory, TestProblem::Element_area, calibrate::error, TestProblem::exact_soln(), MergeRestartFiles::filename, UniformPSDSelfTest::label, sqrt(), and oomph::Problem_Parameter::Trace_file.

Member Data Documentation

◆ Bulk_mesh_pt

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

Pointer to the "bulk" mesh.

◆ Surface_mesh_pt

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

Mesh of boundary condition elements.

◆ Trace_file

template<class ELEMENT >
std::ofstream DarcyProblem< ELEMENT >::Trace_file
private

Trace file for results.


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