mesh_from_geompack_poisson.cc
Go to the documentation of this file.
1 //LIC// ====================================================================
2 //LIC// This file forms part of oomph-lib, the object-oriented,
3 //LIC// multi-physics finite-element library, available
4 //LIC// at http://www.oomph-lib.org.
5 //LIC//
6 //LIC// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 // Driver function for a simple test poisson problem using a mesh
31 // generated from an input file generated by the quadrilateral mesh generator
32 // Geompack.
33 // files .mh2 and .cs2 should be specified in this order
34 
35 //Generic routines
36 #include "generic.h"
37 
38 // Poisson equations
39 #include "poisson.h"
40 
41 // The mesh
42 #include "meshes/geompack_mesh.h"
43 
44 using namespace std;
45 
46 
47 using namespace oomph;
48 
49 //====================================================================
50 /// Namespace for exact solution for Poisson equation with sharp step
51 //====================================================================
53 {
54 
55  /// Parameter for steepness of step
56  double Alpha;
57 
58  /// Parameter for angle of step
59  double Beta;
60 
61 
62  /// Exact solution as a Vector
63  void get_exact_u(const Vector<double>& x, Vector<double>& u)
64  {
65  u[0]=tanh(1.0-Alpha*(Beta*x[0]-x[1]));
66  }
67 
68 
69  /// Exact solution as a scalar
70  void get_exact_u(const Vector<double>& x, double& u)
71  {
72  u=tanh(1.0-Alpha*(Beta*x[0]-x[1]));
73  }
74 
75 
76  /// Source function to make it an exact solution
77  void get_source(const Vector<double>& x, double& source)
78  {
79  source = 2.0*tanh(-1.0+Alpha*(Beta*x[0]-x[1]))*
80  (1.0-pow(tanh(-1.0+Alpha*(Beta*x[0]-x[1])),2.0))*
81  Alpha*Alpha*Beta*Beta+2.0*tanh(-1.0+Alpha*(Beta*x[0]-x[1]))*
82  (1.0-pow(tanh(-1.0+Alpha*(Beta*x[0]-x[1])),2.0))*Alpha*Alpha;
83  }
84 
85 }
86 
87 
88 
89 
90 
91 
92 
93 //====================================================================
94 /// Micky mouse Poisson problem.
95 //====================================================================
96 template<class ELEMENT>
97 class PoissonProblem : public Problem
98 {
99 
100 public:
101 
102 
103  /// Constructor
104  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
105  const string& mesh_file_name,
106  const string& curve_file_name);
107 
108  /// Destructor (empty)
110 
111  /// Update the problem specs before solve: (Re)set boundary conditions
113  {
114  //Loop over the boundaries
115  unsigned num_bound = mesh_pt()->nboundary();
116  for(unsigned ibound=0;ibound<num_bound;ibound++)
117  {
118  // Loop over the nodes on boundary
119  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
120  for (unsigned inod=0;inod<num_nod;inod++)
121  {
122  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
123  double u;
124  Vector<double> x(2);
125  x[0]=nod_pt->x(0);
126  x[1]=nod_pt->x(1);
128  nod_pt->set_value(0,u);
129  }
130  }
131  }
132 
133  /// Update the problem specs before solve (empty)
135  {}
136 
137 
138  //Access function for the specific mesh
140  {
141  return dynamic_cast<GeompackQuadMesh<ELEMENT>*>(Problem::mesh_pt());
142  }
143 
144  /// Doc the solution
145  void doc_solution(DocInfo& doc_info);
146 
147 private:
148 
149  /// Pointer to source function
150  PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
151 
152 };
153 
154 
155 
156 //========================================================================
157 /// Constructor for Poisson problem
158 //========================================================================
159 template<class ELEMENT>
161  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
162  const string& mesh_file_name,
163  const string& curve_file_name)
164  : Source_fct_pt(source_fct_pt)
165 {
166 
167  // Setup parameters for exact tanh solution
168 
169  // Steepness of step
171 
172  // Orientation of step
174 
175  //Create mesh
176  Problem::mesh_pt() = new GeompackQuadMesh<ELEMENT>(mesh_file_name,
177  curve_file_name);
178 
179  // Set the boundary conditions for this problem: All nodes are
180  // free by default -- just pin the ones that have Dirichlet conditions
181  // here.
182  unsigned num_bound = mesh_pt()->nboundary();
183  for(unsigned ibound=0;ibound<num_bound;ibound++)
184  {
185  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
186  for (unsigned inod=0;inod<num_nod;inod++)
187  {
188  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
189  }
190  }
191 
192  // Complete the build of all elements so they are fully functional
193 
194  //Find number of elements in mesh
195  unsigned n_element = mesh_pt()->nelement();
196 
197  // Loop over the elements to set up element-specific
198  // things that cannot be handled by constructor
199  for(unsigned i=0;i<n_element;i++)
200  {
201  // Upcast from GeneralElement to the present element
202  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
203 
204  //Set the source function pointer
205  el_pt->source_fct_pt() = Source_fct_pt;
206  }
207 
208  // Setup equation numbering scheme
209  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
210 
211 }
212 
213 
214 
215 //========================================================================
216 /// Doc the solution
217 //========================================================================
218 template<class ELEMENT>
220 {
221 
222  ofstream some_file;
223  char filename[100];
224 
225  // Number of plot points
226  unsigned npts;
227  npts=5;
228 
229 
230 
231  // Output boundaries
232  //------------------
233  sprintf(filename,"%s/boundaries.dat",doc_info.directory().c_str());
234  some_file.open(filename);
235  mesh_pt()->output_boundaries(some_file);
236  some_file.close();
237 
238 
239  // Output solution
240  //----------------
241  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
242  doc_info.number());
243  some_file.open(filename);
244  mesh_pt()->output(some_file,npts);
245  some_file.close();
246 
247 
248  // Output exact solution
249  //----------------------
250  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
251  doc_info.number());
252  some_file.open(filename);
253  mesh_pt()->output_fct(some_file,npts,TanhSolnForPoisson::get_exact_u);
254  some_file.close();
255 
256 
257 // Doc error
258  //----------
259  double error,norm;
260  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
261  doc_info.number());
262  some_file.open(filename);
263  mesh_pt()->compute_error(some_file,TanhSolnForPoisson::get_exact_u,
264  error,norm);
265  some_file.close();
266  cout << "error: " << sqrt(error) << std::endl;
267  cout << "norm : " << sqrt(norm) << std::endl << std::endl;
268 
269 
270  }
271 
272 
273 
274 
275 
276 
277 
278 ////////////////////////////////////////////////////////////////////////
279 ////////////////////////////////////////////////////////////////////////
280 ////////////////////////////////////////////////////////////////////////
281 
282 
283 
284 //========================================================================
285 /// Demonstrate how to solve Poisson problem
286 //========================================================================
287 int main(int argc, char* argv[])
288 {
289 // Store command line arguments
290  CommandLineArgs::setup(argc,argv);
291 
292  // Check number of command line arguments: Need exactly two.
293  if (argc!=3)
294  {
295  std::string error_message =
296  "Wrong number of command line arguments.\n";
297  error_message +=
298  "Must specify the following file names \n";
299  error_message +=
300  "filename.mh2 then filename.cs2\n";
301 
302  throw OomphLibError(error_message,
303  OOMPH_CURRENT_FUNCTION,
304  OOMPH_EXCEPTION_LOCATION);
305  }
306 
307  // Convert arguments to strings that specify the input file names
308  string mesh_file_name(argv[1]);
309  string curve_file_name(argv[2]);
310 
311  //Set up the problem
314  mesh_file_name,curve_file_name);
315 
316  /// Label for output
317  DocInfo doc_info;
318 
319  // Output directory
320  doc_info.set_directory("RESLT");
321 
322  // Solve the problem
323  problem.newton_solve();
324 
325  //Output solution
326  problem.doc_solution(doc_info);
327 
328 }
329 
330 
331 
int main(int argc, char *argv[])
Demonstrate how to solve Poisson problem.
~PoissonProblem()
Destructor (empty)
void actions_after_newton_solve()
Update the problem specs before solve (empty)
void doc_solution(DocInfo &doc_info)
Doc the solution.
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
GeompackQuadMesh< ELEMENT > * mesh_pt()
double Alpha
Parameter for steepness of step.
void actions_before_newton_solve()
Update the problem specs before solve: (Re)set boundary conditions.
PoissonEquations< 2 >::PoissonSourceFctPt Source_fct_pt
Pointer to source function.
PoissonProblem(PoissonEquations< 2 >::PoissonSourceFctPt source_fct_pt, const string &mesh_file_name, const string &curve_file_name)
Constructor.
void get_exact_u(const Vector< double > &x, double &u)
Exact solution as a scalar.
Namespace for exact solution for Poisson equation with sharp step.
double Beta
Parameter for angle of step.
Micky mouse Poisson problem.
void get_source(const Vector< double > &x, double &source)
Source function to make it an exact solution.