35 #include "meshes/simple_rectangular_quadmesh.h"
39 using namespace oomph;
54 void get_exact_u(
const Vector<double>& x, Vector<double>& u)
60 void get_source(
const Vector<double>& x,
double& source)
82 template<
class ELEMENT>
92 PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
93 const unsigned& nel_1d,
94 const bool& mess_up_order);
113 PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
126 template<
class ELEMENT>
128 PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
129 const unsigned& nel_1d,
const bool& mess_up_order)
130 : Source_fct_pt(source_fct_pt)
148 Problem::mesh_pt() =
new SimpleRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
156 unsigned n_element=mesh_pt()->nelement();
159 Vector<ELEMENT*> tmp_element_pt(n_element);
160 for (
unsigned e=0;e<n_element;e++)
162 tmp_element_pt[e]=
dynamic_cast<ELEMENT*
>(mesh_pt()->element_pt(e));
166 unsigned e_half=unsigned(0.5*
double(n_element));
167 unsigned e_lo=e_half-3;
168 unsigned e_up=n_element-e_lo;
170 for (
unsigned e=0;e<e_lo;e++)
172 mesh_pt()->element_pt(count)=tmp_element_pt[e];
174 mesh_pt()->element_pt(count)=tmp_element_pt[n_element-e-1];
177 for (
unsigned e=e_lo;e<e_up;e++)
179 mesh_pt()->element_pt(count)=tmp_element_pt[e];
187 unsigned num_bound = mesh_pt()->nboundary();
188 for(
unsigned ibound=0;ibound<num_bound;ibound++)
190 unsigned num_nod= mesh_pt()->nboundary_node(ibound);
191 for (
unsigned inod=0;inod<num_nod;inod++)
193 mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
202 unsigned n_element = mesh_pt()->nelement();
203 for(
unsigned i=0;i<n_element;i++)
206 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(mesh_pt()->element_pt(i));
214 cout <<
"Number of equations: " << assign_eqn_numbers() << std::endl;
225 template<
class ELEMENT>
229 unsigned num_bound = mesh_pt()->nboundary();
232 for(
unsigned ibound=0;ibound<num_bound;ibound++)
235 unsigned num_nod=mesh_pt()->nboundary_node(ibound);
238 for (
unsigned inod=0;inod<num_nod;inod++)
241 Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
253 nod_pt->set_value(0,u[0]);
263 template<
class ELEMENT>
275 sprintf(filename,
"%s/soln%i.dat",doc_info.directory().c_str(),
277 some_file.open(filename);
278 mesh_pt()->output(some_file,npts);
284 sprintf(filename,
"%s/exact_soln%i.dat",doc_info.directory().c_str(),
286 some_file.open(filename);
293 sprintf(filename,
"%s/error%i.dat",doc_info.directory().c_str(),
295 some_file.open(filename);
301 cout <<
"\nNorm of error : " << sqrt(error) << std::endl;
302 cout <<
"Norm of solution: " << sqrt(norm) << std::endl << std::endl;
317 void run(
const string& dir_name,
318 LinearSolver* linear_solver_pt,
319 const unsigned nel_1d,
335 problem.linear_solver_pt() = linear_solver_pt;
343 doc_info.set_directory(dir_name);
350 cout <<
"\n\n\nProblem self-test:";
351 if (problem.self_test()==0)
353 cout <<
"passed: Problem can be solved." << std::endl;
357 throw OomphLibError(
"Self test failed",
358 OOMPH_CURRENT_FUNCTION,
359 OOMPH_EXCEPTION_LOCATION);
370 problem.newton_solve();
389 LinearSolver* linear_solver_pt;
395 clock_t cpu_start,cpu_end;
398 Vector<double> superlu_cr_cpu(2);
399 Vector<double> superlu_cc_cpu(2);
400 Vector<double> hsl_ma42_cpu(2);
401 Vector<double> hsl_ma42_reordered_cpu(2);
402 Vector<double> dense_lu_cpu(2);
403 Vector<double> fd_lu_cpu(2);
413 for (
unsigned mess=0;mess<2;mess++)
430 cout <<
" Use SuperLU with compressed row storage: " << std::endl;
431 cout <<
"========================================= " << std::endl;
438 linear_solver_pt =
new SuperLUSolver;
441 static_cast<SuperLUSolver*
>(linear_solver_pt)
442 ->use_compressed_row_for_superlu_serial();
445 static_cast<SuperLUSolver*
>(linear_solver_pt)->enable_doc_stats();
455 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
464 superlu_cr_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
474 cout <<
" Use SuperLU with compressed column storage: " << std::endl;
475 cout <<
"============================================ " << std::endl;
483 linear_solver_pt =
new SuperLUSolver;
486 static_cast<SuperLUSolver*
>(linear_solver_pt)
487 ->use_compressed_column_for_superlu_serial();
490 static_cast<SuperLUSolver*
>(linear_solver_pt)->enable_doc_stats();
500 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
509 superlu_cc_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
512 #ifdef HAVE_HSL_SOURCES
518 cout <<
" Use HSL frontal solver MA42 without element re-ordering: " << std::endl;
519 cout <<
"========================================================= " << std::endl;
526 linear_solver_pt =
new HSL_MA42;
529 static_cast<HSL_MA42*
>(linear_solver_pt)->enable_doc_stats();
532 dir_name=
"RESLT_frontal";
540 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
550 hsl_ma42_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
558 cout <<
" Use HSL frontal solver MA42 with element re-ordering: " << std::endl;
559 cout <<
"====================================================== " << std::endl;
566 linear_solver_pt =
new HSL_MA42;
569 static_cast<HSL_MA42*
>(linear_solver_pt)->enable_doc_stats();
573 static_cast<HSL_MA42*
>(linear_solver_pt)->enable_reordering();
576 dir_name=
"RESLT_frontal_reordered";
584 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
594 hsl_ma42_reordered_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
603 cout <<
" Use dense matrix LU decomposition: " << std::endl;
604 cout <<
"=================================== " << std::endl;
611 linear_solver_pt =
new DenseLU;
614 dir_name=
"RESLT_dense_LU";
621 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
631 dense_lu_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
640 cout <<
" Use dense FD-ed matrix LU decomposition: " << std::endl;
641 cout <<
"========================================= " << std::endl;
648 linear_solver_pt =
new FD_LU;
651 dir_name=
"RESLT_FD_LU";
658 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
667 fd_lu_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
673 for (
unsigned mess=0;mess<2;mess++)
675 cout << std::endl << std::endl << std::endl ;
678 cout <<
"TIMINGS WITH MESSED UP ORDERING OF ELEMENTS: " << std::endl;
679 cout <<
"============================================ " << std::endl;
684 cout <<
"TIMINGS WITHOUT MESSED UP ORDERING OF ELEMENTS: " << std::endl;
685 cout <<
"=============================================== " << std::endl;
688 cout <<
"CPU time with SuperLU compressed row : "
689 << superlu_cr_cpu[mess] << std::endl;
690 cout <<
"CPU time with SuperLU compressed col : "
691 << superlu_cc_cpu[mess] << std::endl;
692 #ifdef HAVE_HSL_SOURCES
693 cout <<
"CPU time with MA42 frontal solver : "
694 << hsl_ma42_cpu[mess] << std::endl;
695 cout <<
"CPU time with MA42 frontal solver (incl. reordering) : "
696 << hsl_ma42_reordered_cpu[mess] << std::endl;
698 cout <<
"CPU time with dense LU solver (fewer elements!) : "
699 << dense_lu_cpu[mess] << std::endl;
700 cout <<
"CPU time with dense LU solver & FD (fewer elements!) : "
701 << fd_lu_cpu[mess] << std::endl;
702 cout << std::endl << std::endl << std::endl ;
2D Poisson problem on rectangular domain, discretised with 2D QPoisson elements. The specific type of...
PoissonEquations< 2 >::PoissonSourceFctPt Source_fct_pt
Pointer to source function.
void actions_before_newton_solve()
Update the problem specs before solve: Reset boundary conditions to the values from the exact solutio...
void actions_after_newton_solve()
Update the problem after solve (empty)
PoissonProblem(PoissonEquations< 2 >::PoissonSourceFctPt source_fct_pt)
Constructor: Pass pointer to source function.
void doc_solution(DocInfo &doc_info)
Doc the solution. DocInfo object stores flags/labels for where the output gets written to.
~PoissonProblem()
Destructor (empty)
Namespace for exact solution for Poisson equation with "sharp step".
double TanPhi
Parameter for angle Phi of "step".
void get_source(const Vector< double > &x, double &source)
Source function required to make the solution above an exact solution.
double Alpha
Parameter for steepness of "step".
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
void run(const string &dir_name, LinearSolver *linear_solver_pt, const unsigned nel_1d, bool mess_up_order)
Build and run problem with specified linear solver. Also pass flag to specify if order of elements sh...
int main()
Driver code for 2D Poisson problem – compare different solvers.