42 #include "meshes/simple_rectangular_quadmesh.template.h" 46 using namespace oomph;
61 void get_exact_u(
const Vector<double>& x, Vector<double>& u)
63 u[0]=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
69 source = 2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
70 (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*
71 Alpha*Alpha*TanPhi*TanPhi+2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
72 (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*Alpha*Alpha;
89 template<
class ELEMENT>
96 PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt);
103 void actions_before_newton_solve();
110 void doc_solution(DocInfo& doc_info);
115 PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
125 template<
class ELEMENT>
127 PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt)
128 : Source_fct_pt(source_fct_pt)
146 Problem::mesh_pt() =
new SimpleRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
151 unsigned num_bound = mesh_pt()->nboundary();
152 for(
unsigned ibound=0;ibound<num_bound;ibound++)
154 unsigned num_nod= mesh_pt()->nboundary_node(ibound);
155 for (
unsigned inod=0;inod<num_nod;inod++)
157 mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
166 unsigned n_element = mesh_pt()->nelement();
167 for(
unsigned i=0;i<n_element;i++)
170 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(mesh_pt()->element_pt(i));
173 el_pt->source_fct_pt() = Source_fct_pt;
178 cout <<
"Number of equations: " << assign_eqn_numbers() << std::endl;
189 template<
class ELEMENT>
193 unsigned num_bound = mesh_pt()->nboundary();
196 for(
unsigned ibound=0;ibound<num_bound;ibound++)
199 unsigned num_nod=mesh_pt()->nboundary_node(ibound);
202 for (
unsigned inod=0;inod<num_nod;inod++)
205 Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
217 nod_pt->set_value(0,u[0]);
227 template<
class ELEMENT>
239 sprintf(filename,
"%s/soln%i.dat",doc_info.directory().c_str(),
241 some_file.open(filename);
242 mesh_pt()->output(some_file,npts);
248 sprintf(filename,
"%s/exact_soln%i.dat",doc_info.directory().c_str(),
250 some_file.open(filename);
257 sprintf(filename,
"%s/error%i.dat",doc_info.directory().c_str(),
259 some_file.open(filename);
265 cout <<
"\nNorm of error : " << sqrt(error) << std::endl;
266 cout <<
"Norm of solution: " << sqrt(norm) << std::endl << std::endl;
294 doc_info.set_directory(
"RESLT");
301 cout <<
"\n\n\nProblem self-test ";
302 if (problem.self_test()==0)
304 cout <<
"passed: Problem can be solved." << std::endl;
308 throw OomphLibError(
"Self test failed",
309 OOMPH_CURRENT_FUNCTION,
310 OOMPH_EXCEPTION_LOCATION);
324 for (
unsigned istep=0;istep<nstep;istep++)
330 cout <<
"\n\nSolving for TanhSolnForPoisson::Alpha=" 334 problem.newton_solve();
~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.
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.
int main()
Driver code for 2D Poisson problem.