37 #include "meshes/fish_mesh.h"
44 using namespace oomph;
79 template<
class ELEMENT>
89 const bool& fix_position,
const string& directory_name,
90 const unsigned& i_case);
99 fish_mesh_pt()->node_update();
108 fish_mesh_pt()->node_update();
120 return *Load_pt->value_pt(0);
128 return static_cast<ElasticallySupportedRingElement*
>(fish_mesh_pt()->
129 fish_back_pt())->y_c();
147 unsigned n_element = fish_mesh_pt()->nelement();
148 for(
unsigned i=0;i<n_element;i++)
151 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(fish_mesh_pt()->element_pt(i));
156 el_pt->evaluate_shape_derivs_by_direct_fd();
157 if (!done) std::cout <<
"\n\n [CR residuals] Direct FD" << std::endl;
160 else if ( (Case_id==1) || (Case_id==2) )
163 bool i_know_what_i_am_doing=
true;
164 el_pt->evaluate_shape_derivs_by_chain_rule(i_know_what_i_am_doing);
167 el_pt->enable_always_evaluate_dresidual_dnodal_coordinates_by_fd();
168 if (!done) std::cout <<
"\n\n [CR residuals] Chain rule and FD"
173 el_pt->disable_always_evaluate_dresidual_dnodal_coordinates_by_fd();
174 if (!done) std::cout <<
"\n\n [CR residuals] Chain rule and analytic"
179 else if ( (Case_id==3) || (Case_id==4) )
182 bool i_know_what_i_am_doing=
true;
183 el_pt->evaluate_shape_derivs_by_fastest_method(i_know_what_i_am_doing);
186 el_pt->enable_always_evaluate_dresidual_dnodal_coordinates_by_fd();
187 if (!done) std::cout <<
"\n\n [CR residuals] Fastest and FD"
192 el_pt->disable_always_evaluate_dresidual_dnodal_coordinates_by_fd();
193 if (!done) std::cout <<
"\n\n [CR residuals] Fastest and analytic"
239 template<
class ELEMENT>
241 const bool& fix_position,
const string& directory_name,
242 const unsigned& i_case) : Fix_position(fix_position), Case_id(i_case)
247 Doc_info.set_directory(directory_name);
254 sprintf(filename,
"%s/trace.dat",directory_name.c_str());
258 <<
"VARIABLES=\"load\",\"y<sub>circle</sub>\",\"u<sub>control</sub>\""
267 GeomObject* fish_back_pt=
new ElasticallySupportedRingElement(x_c,
y_c,r_back);
270 Fish_mesh_pt=
new AlgebraicRefineableFishMesh<ELEMENT>(fish_back_pt);
292 fish_mesh_pt()->spatial_error_estimator_pt()=
new Z2ErrorEstimator;
299 unsigned nnod=
fish_mesh_pt()->finite_element_pt(0)->nnode();
305 cout << std::endl <<
"Control node is located at: "
307 << std::endl << std::endl;
319 dynamic_cast<ElasticallySupportedRingElement*
>(
336 dynamic_cast<ElasticallySupportedRingElement*
>(
Fish_mesh_pt->fish_back_pt())->
344 for(
unsigned ibound=0;ibound<num_bound;ibound++)
346 unsigned num_nod=
fish_mesh_pt()->nboundary_node(ibound);
347 for (
unsigned inod=0;inod<num_nod;inod++)
355 for(
unsigned ibound=0;ibound<num_bound;ibound++)
358 unsigned num_nod=
fish_mesh_pt()->nboundary_node(ibound);
359 for (
unsigned inod=0;inod<num_nod;inod++)
361 fish_mesh_pt()->boundary_node_pt(ibound,inod)->set_value(0,0.0);
367 for(
unsigned i=0;i<n_element;i++)
370 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(
fish_mesh_pt()->element_pt(i));
380 cout <<
"Number of equations: " << assign_eqn_numbers() << std::endl;
389 template<
class ELEMENT>
403 template<
class ELEMENT>
418 sprintf(filename,
"%s/soln_%i_%i.dat",Doc_info.directory().c_str(),
419 Case_id,Doc_info.number());
423 sprintf(filename,
"%s/soln%i.dat",Doc_info.directory().c_str(),
426 some_file.open(filename);
427 fish_mesh_pt()->output(some_file,npts);
433 <<
static_cast<ElasticallySupportedRingElement*
>(fish_mesh_pt()->
434 fish_back_pt())->load()
436 <<
static_cast<ElasticallySupportedRingElement*
>(fish_mesh_pt()->
437 fish_back_pt())->y_c()
438 <<
" " << Doc_node_pt->value(0) << std::endl;
459 template<
class ELEMENT>
464 bool fix_position=
true;
472 problem.refine_uniformly();
473 problem.refine_uniformly();
485 double dyc=0.6/double(nstep-1);
488 if (CommandLineArgs::Argc>1) nstep=1;
490 for (
unsigned istep=0;istep<nstep;istep++)
493 unsigned max_solve=2;
494 problem.newton_solve(max_solve);
510 template<
class ELEMENT>
515 for (
unsigned i_case=0;i_case<5;i_case++)
518 std::cout <<
"[CR residuals] " << std::endl;
519 std::cout <<
"[CR residuals]=================================================="
521 std::cout <<
"[CR residuals] " << std::endl;
523 bool fix_position=
false;
533 problem.refine_uniformly();
534 problem.refine_uniformly();
541 unsigned max_solve=2;
542 problem.newton_solve(max_solve);
558 int main(
int argc,
char* argv[])
561 CommandLineArgs::setup(argc,argv);
564 typedef AlgebraicElement<RefineableQPoissonElement<2,3> > ELEMENT;
567 demo_fish_poisson<ELEMENT>(
"RESLT");
570 demo_elastic_fish_poisson<ELEMENT>(
"RESLT_coupled");
int main(int argc, char *argv[])
Driver for "elastic" fish poisson solver with adaptation. If there are any command line arguments,...
void demo_elastic_fish_poisson(const string &directory_name)
Demonstrate how to solve "elastic" 2D Poisson problem in deformable fish-shaped domain with mesh adap...
void demo_fish_poisson(const string &directory_name)
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
DocInfo & doc_info()
Access to DocInfo object.
void set_shape_deriv_method()
Helper fct to set method for evaluation of shape derivs.
double & y_c()
Return value of the vertical displacement of the ring that represents the fish's back.
double & load()
Return value of the "load" on the elastically supported ring.
RefineableFishPoissonProblem(const bool &fix_position, const string &directory_name, const unsigned &i_case)
Constructor: Bool flag specifies if position of fish back is prescribed or computed from the coupled ...
bool Fix_position
Is the position of the fish back prescribed?
void actions_before_newton_solve()
Update the problem specs before solve: Update mesh.
Node * Doc_node_pt
Node at which the solution of the Poisson equation is documented.
void doc_solution()
Doc the solution.
AlgebraicRefineableFishMesh< ELEMENT > * fish_mesh_pt()
virtual ~RefineableFishPoissonProblem()
Destructor.
void actions_after_newton_solve()
Update the problem specs after solve (empty)
void actions_before_newton_convergence_check()
Update after Newton step: Update mesh in response to possible changes in the wall shape.
Data * Load_pt
Pointer to data item that stores the "load" on the fish back.
ofstream Trace_file
Trace file.
Mesh * Fish_back_mesh_pt
Pointer to single-element mesh that stores the GeneralisedElement that represents the fish back.
DocInfo Doc_info
Doc info object.
AlgebraicRefineableFishMesh< ELEMENT > * Fish_mesh_pt
Pointer to fish mesh.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void get_source(const Vector< double > &x, double &source)
Const source function.
double Strength
Strength of source function: default value 1.0.