39 #include "meshes/simple_rectangular_quadmesh.h" 43 using namespace oomph;
58 void get_exact_u(
const Vector<double>& x, Vector<double>& u)
60 u[0]=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
64 void get_source(
const Vector<double>& x,
double& source)
66 source = 2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
67 (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*
68 Alpha*Alpha*TanPhi*TanPhi+2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
69 (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*Alpha*Alpha;
86 template<
class ELEMENT>
96 PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
97 const unsigned& nel_1d,
98 const bool& mess_up_order);
105 void actions_before_newton_solve();
112 void doc_solution(DocInfo& doc_info);
117 PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
130 template<
class ELEMENT>
133 const unsigned& nel_1d,
const bool& mess_up_order)
134 : Source_fct_pt(source_fct_pt)
152 Problem::mesh_pt() =
new SimpleRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
160 unsigned n_element=mesh_pt()->nelement();
163 Vector<ELEMENT*> tmp_element_pt(n_element);
164 for (
unsigned e=0;e<n_element;e++)
166 tmp_element_pt[e]=
dynamic_cast<ELEMENT*
>(mesh_pt()->element_pt(e));
170 unsigned e_half=unsigned(0.5*
double(n_element));
171 unsigned e_lo=e_half-3;
172 unsigned e_up=n_element-e_lo;
174 for (
unsigned e=0;e<e_lo;e++)
176 mesh_pt()->element_pt(count)=tmp_element_pt[e];
178 mesh_pt()->element_pt(count)=tmp_element_pt[n_element-e-1];
181 for (
unsigned e=e_lo;e<e_up;e++)
183 mesh_pt()->element_pt(count)=tmp_element_pt[e];
191 unsigned num_bound = mesh_pt()->nboundary();
192 for(
unsigned ibound=0;ibound<num_bound;ibound++)
194 unsigned num_nod= mesh_pt()->nboundary_node(ibound);
195 for (
unsigned inod=0;inod<num_nod;inod++)
197 mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
206 unsigned n_element = mesh_pt()->nelement();
207 for(
unsigned i=0;i<n_element;i++)
210 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(mesh_pt()->element_pt(i));
218 cout <<
"Number of equations: " << assign_eqn_numbers() << std::endl;
229 template<
class ELEMENT>
233 unsigned num_bound = mesh_pt()->nboundary();
236 for(
unsigned ibound=0;ibound<num_bound;ibound++)
239 unsigned num_nod=mesh_pt()->nboundary_node(ibound);
242 for (
unsigned inod=0;inod<num_nod;inod++)
245 Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
257 nod_pt->set_value(0,u[0]);
267 template<
class ELEMENT>
279 sprintf(filename,
"%s/soln%i.dat",doc_info.directory().c_str(),
281 some_file.open(filename);
282 mesh_pt()->output(some_file,npts);
288 sprintf(filename,
"%s/exact_soln%i.dat",doc_info.directory().c_str(),
290 some_file.open(filename);
297 sprintf(filename,
"%s/error%i.dat",doc_info.directory().c_str(),
299 some_file.open(filename);
305 cout <<
"\nNorm of error : " << sqrt(error) << std::endl;
306 cout <<
"Norm of solution: " << sqrt(norm) << std::endl << std::endl;
321 void run(
const string& dir_name,
322 LinearSolver* linear_solver_pt,
323 const unsigned nel_1d,
339 problem.linear_solver_pt() = linear_solver_pt;
347 doc_info.set_directory(dir_name);
354 cout <<
"\n\n\nProblem self-test:";
355 if (problem.self_test()==0)
357 cout <<
"passed: Problem can be solved." << std::endl;
361 throw OomphLibError(
"Self test failed",
362 OOMPH_CURRENT_FUNCTION,
363 OOMPH_EXCEPTION_LOCATION);
374 problem.newton_solve();
393 LinearSolver* linear_solver_pt;
399 clock_t cpu_start,cpu_end;
402 Vector<double> superlu_cr_cpu(2);
403 Vector<double> superlu_cc_cpu(2);
404 Vector<double> hsl_ma42_cpu(2);
405 Vector<double> hsl_ma42_reordered_cpu(2);
406 Vector<double> dense_lu_cpu(2);
407 Vector<double> fd_lu_cpu(2);
417 for (
unsigned mess=0;mess<2;mess++)
434 cout <<
" Use SuperLU with compressed row storage: " << std::endl;
435 cout <<
"========================================= " << std::endl;
442 linear_solver_pt =
new SuperLUSolver;
445 static_cast<SuperLUSolver*
>(linear_solver_pt)
446 ->use_compressed_row_for_superlu_serial();
449 static_cast<SuperLUSolver*
>(linear_solver_pt)->enable_doc_stats();
459 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
468 superlu_cr_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
478 cout <<
" Use SuperLU with compressed column storage: " << std::endl;
479 cout <<
"============================================ " << std::endl;
487 linear_solver_pt =
new SuperLUSolver;
490 static_cast<SuperLUSolver*
>(linear_solver_pt)
491 ->use_compressed_column_for_superlu_serial();
494 static_cast<SuperLUSolver*
>(linear_solver_pt)->enable_doc_stats();
504 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
513 superlu_cc_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
516 #ifdef HAVE_HSL_SOURCES 522 cout <<
" Use HSL frontal solver MA42 without element re-ordering: " << std::endl;
523 cout <<
"========================================================= " << std::endl;
530 linear_solver_pt =
new HSL_MA42;
533 static_cast<HSL_MA42*
>(linear_solver_pt)->enable_doc_stats();
536 dir_name=
"RESLT_frontal";
544 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
554 hsl_ma42_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
562 cout <<
" Use HSL frontal solver MA42 with element re-ordering: " << std::endl;
563 cout <<
"====================================================== " << std::endl;
570 linear_solver_pt =
new HSL_MA42;
573 static_cast<HSL_MA42*
>(linear_solver_pt)->enable_doc_stats();
577 static_cast<HSL_MA42*
>(linear_solver_pt)->enable_reordering();
580 dir_name=
"RESLT_frontal_reordered";
588 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
598 hsl_ma42_reordered_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
607 cout <<
" Use dense matrix LU decomposition: " << std::endl;
608 cout <<
"=================================== " << std::endl;
615 linear_solver_pt =
new DenseLU;
618 dir_name=
"RESLT_dense_LU";
625 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
635 dense_lu_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
644 cout <<
" Use dense FD-ed matrix LU decomposition: " << std::endl;
645 cout <<
"========================================= " << std::endl;
652 linear_solver_pt =
new FD_LU;
655 dir_name=
"RESLT_FD_LU";
662 run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
671 fd_lu_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
677 for (
unsigned mess=0;mess<2;mess++)
679 cout << std::endl << std::endl << std::endl ;
682 cout <<
"TIMINGS WITH MESSED UP ORDERING OF ELEMENTS: " << std::endl;
683 cout <<
"============================================ " << std::endl;
688 cout <<
"TIMINGS WITHOUT MESSED UP ORDERING OF ELEMENTS: " << std::endl;
689 cout <<
"=============================================== " << std::endl;
692 cout <<
"CPU time with SuperLU compressed row : " 693 << superlu_cr_cpu[mess] << std::endl;
694 cout <<
"CPU time with SuperLU compressed col : " 695 << superlu_cc_cpu[mess] << std::endl;
696 #ifdef HAVE_HSL_SOURCES 697 cout <<
"CPU time with MA42 frontal solver : " 698 << hsl_ma42_cpu[mess] << std::endl;
699 cout <<
"CPU time with MA42 frontal solver (incl. reordering) : " 700 << hsl_ma42_reordered_cpu[mess] << std::endl;
702 cout <<
"CPU time with dense LU solver (fewer elements!) : " 703 << dense_lu_cpu[mess] << std::endl;
704 cout <<
"CPU time with dense LU solver & FD (fewer elements!) : " 705 << fd_lu_cpu[mess] << std::endl;
706 cout << std::endl << std::endl << std::endl ;
int main()
Driver code for 2D Poisson problem – compare different solvers.
PoissonEquations< 2 >::PoissonSourceFctPt Source_fct_pt
Pointer to source function.
~PoissonProblem()
Destructor (empty)
void actions_after_newton_solve()
Update the problem after solve (empty)
void doc_solution(DocInfo &doc_info)
Doc the solution. DocInfo object stores flags/labels for where the output gets written to...
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)
double Alpha
Parameter for steepness of "step".
void actions_before_newton_solve()
Update the problem specs before solve: Reset boundary conditions to the values from the exact solutio...
double TanPhi
Parameter for angle Phi of "step".
Namespace for exact solution for Poisson equation with "sharp step".
PoissonProblem(PoissonEquations< 2 >::PoissonSourceFctPt source_fct_pt)
Constructor: Pass pointer to source function.
void get_source(const Vector< double > &x, double &source)
Source function required to make the solution above an exact solution.