Toggle navigation
Documentation
Big picture
The finite element method
The data structure
Not-so-quick guide
Optimisation
Order of action functions
Example codes and tutorials
List of example codes and tutorials
Meshing
Solvers
MPI parallel processing
Post-processing/visualisation
Other
Change log
Creating documentation
Coding conventions
Index
FAQ
About
People
Contact/Get involved
Publications
Acknowledgements
Copyright
Picture show
Go
demo_drivers
meshing
mesh_from_triangle
mesh_from_triangle_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// Copyright (C) 2006-2025 Matthias Heil and Andrew Hazel
7
//LIC//
8
//LIC// This library is free software; you can redistribute it and/or
9
//LIC// modify it under the terms of the GNU Lesser General Public
10
//LIC// License as published by the Free Software Foundation; either
11
//LIC// version 2.1 of the License, or (at your option) any later version.
12
//LIC//
13
//LIC// This library is distributed in the hope that it will be useful,
14
//LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15
//LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16
//LIC// Lesser General Public License for more details.
17
//LIC//
18
//LIC// You should have received a copy of the GNU Lesser General Public
19
//LIC// License along with this library; if not, write to the Free Software
20
//LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21
//LIC// 02110-1301 USA.
22
//LIC//
23
//LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24
//LIC//
25
//LIC//====================================================================
26
// Driver code for a simple test poisson problem using a mesh
27
// generated from an input file generated by the triangle mesh generator
28
// Triangle.
29
30
//Generic routines
31
#include "generic.h"
32
33
// Poisson equations
34
#include "poisson.h"
35
36
// The mesh
37
#include "
meshes/triangle_mesh.h
"
38
39
using namespace
std;
40
41
using namespace
oomph
;
42
43
//====================================================================
44
/// Namespace for exact solution for Poisson equation with sharp step
45
//====================================================================
46
namespace
TanhSolnForPoisson
47
{
48
49
/// Parameter for steepness of step
50
double
Alpha
;
51
52
/// Parameter for angle of step
53
double
Beta
;
54
55
56
/// Exact solution as a Vector
57
void
get_exact_u
(
const
Vector<double>
& x,
Vector<double>
&
u
)
58
{
59
u
[0]=
tanh
(1.0-
Alpha
*(
Beta
*x[0]-x[1]));
60
}
61
62
63
/// Exact solution as a scalar
64
void
get_exact_u
(
const
Vector<double>
& x,
double
&
u
)
65
{
66
u
=
tanh
(1.0-
Alpha
*(
Beta
*x[0]-x[1]));
67
}
68
69
70
/// Source function to make it an exact solution
71
void
get_source
(
const
Vector<double>
& x,
double
&
source
)
72
{
73
source
= 2.0*
tanh
(-1.0+
Alpha
*(
Beta
*x[0]-x[1]))*
74
(1.0-
pow
(
tanh
(-1.0+
Alpha
*(
Beta
*x[0]-x[1])),2.0))*
75
Alpha
*
Alpha
*
Beta
*
Beta
+2.0*
tanh
(-1.0+
Alpha
*(
Beta
*x[0]-x[1]))*
76
(1.0-
pow
(
tanh
(-1.0+
Alpha
*(
Beta
*x[0]-x[1])),2.0))*
Alpha
*
Alpha
;
77
}
78
79
}
80
81
82
83
84
85
86
87
//====================================================================
88
/// Micky mouse Poisson problem.
89
//====================================================================
90
91
// Poisson problem
92
template
<
class
ELEMENT>
93
class
PoissonProblem
:
public
Problem
94
{
95
96
public
:
97
98
99
/// Constructor: Pass pointer to source function and names of
100
/// two triangle input files
101
PoissonProblem
(
PoissonEquations<2>::PoissonSourceFctPt
source_fct_pt
,
102
const
string
&
node_file_name
,
103
const
string
&
element_file_name
,
104
const
string
&
poly_file_name
);
105
106
/// Destructor (empty)
107
~PoissonProblem
(){}
108
109
/// Update the problem specs before solve: (Re)set boundary conditions
110
void
actions_before_newton_solve
()
111
{
112
//Loop over the boundaries
113
unsigned
num_bound
=
mesh_pt
()->nboundary();
114
for
(
unsigned
ibound
=0;
ibound
<
num_bound
;
ibound
++)
115
{
116
// Loop over the nodes on boundary
117
unsigned
num_nod
=
mesh_pt
()->nboundary_node(
ibound
);
118
for
(
unsigned
inod
=0;
inod
<
num_nod
;
inod
++)
119
{
120
Node
*
nod_pt
=
mesh_pt
()->boundary_node_pt(
ibound
,
inod
);
121
double
u
;
122
Vector<double>
x(2);
123
x[0]=
nod_pt
->x(0);
124
x[1]=
nod_pt
->x(1);
125
TanhSolnForPoisson::get_exact_u
(x,
u
);
126
nod_pt
->set_value(0,
u
);
127
}
128
}
129
}
130
131
/// Update the problem specs before solve (empty)
132
void
actions_after_newton_solve
()
133
{}
134
135
#ifdef ADAPT
136
/// Actions performed after the adaptation steps
137
void
actions_after_adapt
();
138
#endif
139
140
#ifdef ADAPT
141
/// Access function for the specific mesh
142
RefineableTriangleMesh<ELEMENT>
*
mesh_pt
()
143
{
144
return
dynamic_cast<
RefineableTriangleMesh<ELEMENT>
*
>
(Problem::mesh_pt());
145
}
146
#else
147
/// Access function for the specific mesh
148
TriangleMesh<ELEMENT>
*
mesh_pt
()
149
{
150
return
dynamic_cast<
TriangleMesh<ELEMENT>
*
>
(Problem::mesh_pt());
151
}
152
#endif
153
154
/// Doc the solution
155
void
doc_solution
(
DocInfo
&
doc_info
);
156
157
private
:
158
159
/// Pointer to source function
160
PoissonEquations<2>::PoissonSourceFctPt
Source_fct_pt
;
161
162
#ifdef ADAPT
163
/// Pointer to the bulk mesh
164
RefineableTriangleMesh<ELEMENT>
*
Bulk_mesh_pt
;
165
166
/// Error estimator
167
Z2ErrorEstimator
*
error_estimator_pt
;
168
#endif
169
170
};
171
172
173
174
//========================================================================
175
/// Constructor for Poisson problem
176
//========================================================================
177
template
<
class
ELEMENT>
178
PoissonProblem<ELEMENT>::
179
PoissonProblem
(
PoissonEquations<2>::PoissonSourceFctPt
source_fct_pt
,
180
const
string
&
node_file_name
,
181
const
string
&
element_file_name
,
182
const
string
&
poly_file_name
)
183
: Source_fct_pt(
source_fct_pt
)
184
{
185
186
// Setup parameters for exact tanh solution
187
188
// Steepness of step
189
TanhSolnForPoisson::Alpha
=1.0;
190
191
// Orientation of step
192
TanhSolnForPoisson::Beta
=1.4;
193
194
#ifdef ADAPT
195
//Create mesh
196
Bulk_mesh_pt
=
new
RefineableTriangleMesh<ELEMENT>
(
node_file_name
,
197
element_file_name
,
198
poly_file_name
);
199
200
// Set error estimator for bulk mesh
201
Z2ErrorEstimator
*
error_estimator_pt
=
new
Z2ErrorEstimator
;
202
Bulk_mesh_pt
->spatial_error_estimator_pt() =
error_estimator_pt
;
203
204
// Set the problem mesh pointer to the bulk mesh
205
Problem::mesh_pt() =
Bulk_mesh_pt
;
206
207
#else
208
//Create mesh
209
Problem::mesh_pt() =
new
TriangleMesh<ELEMENT>
(
node_file_name
,
210
element_file_name
,
211
poly_file_name
);
212
#endif
213
214
// Set the boundary conditions for this problem: All nodes are
215
// free by default -- just pin the ones that have Dirichlet conditions
216
// here.
217
unsigned
num_bound
=
mesh_pt
()->nboundary();
218
for
(
unsigned
ibound
=0;
ibound
<
num_bound
;
ibound
++)
219
{
220
unsigned
num_nod
=
mesh_pt
()->nboundary_node(
ibound
);
221
for
(
unsigned
inod
=0;
inod
<
num_nod
;
inod
++)
222
{
223
mesh_pt
()->boundary_node_pt(
ibound
,
inod
)->pin(0);
224
}
225
}
226
227
// Complete the build of all elements so they are fully functional
228
229
//Find number of elements in mesh
230
unsigned
n_element
=
mesh_pt
()->nelement();
231
232
// Loop over the elements to set up element-specific
233
// things that cannot be handled by constructor
234
for
(
unsigned
i
=0;
i
<
n_element
;
i
++)
235
{
236
// Upcast from GeneralElement to the present element
237
ELEMENT
*
el_pt
=
dynamic_cast<
ELEMENT
*
>
(
mesh_pt
()->element_pt(
i
));
238
239
//Set the source function pointer
240
el_pt
->source_fct_pt() =
Source_fct_pt
;
241
}
242
243
// Setup equation numbering scheme
244
cout
<<
"Number of equations: "
<<
assign_eqn_numbers
() << std::endl;
245
246
}
247
248
#ifdef ADAPT
249
//========================================================================
250
/// Actions performed after the adaptation steps
251
//========================================================================
252
template
<
class
ELEMENT>
253
void
PoissonProblem<ELEMENT>::actions_after_adapt
()
254
{
255
// Since the mesh adaptation strategy is based on mesh re-generation
256
// we need to pin the nodes in the new regenerated mesh, and pass the
257
// source function to the elements
258
259
// Set the boundary conditions for this problem: All nodes are free
260
// by default -- just pin the ones that have Dirichlet conditions
261
// here.
262
unsigned
num_bound
= mesh_pt()->nboundary();
263
for
(
unsigned
ibound
=0;
ibound
<
num_bound
;
ibound
++)
264
{
265
unsigned
num_nod
= mesh_pt()->nboundary_node(
ibound
);
266
for
(
unsigned
inod
=0;
inod
<
num_nod
;
inod
++)
267
{
268
mesh_pt()->boundary_node_pt(
ibound
,
inod
)->pin(0);
269
}
270
}
271
272
// Complete the build of all elements so they are fully functional
273
274
//Find number of elements in mesh
275
unsigned
n_element
= mesh_pt()->nelement();
276
277
// Loop over the elements to set up element-specific
278
// things that cannot be handled by constructor
279
for
(
unsigned
i
=0;
i
<
n_element
;
i
++)
280
{
281
// Upcast from GeneralElement to the present element
282
ELEMENT
*
el_pt
=
dynamic_cast<
ELEMENT
*
>
(mesh_pt()->element_pt(
i
));
283
284
//Set the source function pointer
285
el_pt
->source_fct_pt() = Source_fct_pt;
286
}
287
288
// Setup equation numbering scheme
289
cout
<<
"Number of equations: "
<<
assign_eqn_numbers
() << std::endl;
290
}
291
#endif
292
293
294
//========================================================================
295
/// Doc the solution
296
//========================================================================
297
template
<
class
ELEMENT>
298
void
PoissonProblem<ELEMENT>::doc_solution
(
DocInfo
&
doc_info
)
299
{
300
301
ofstream
some_file
;
302
char
filename
[100];
303
304
// Number of plot points
305
unsigned
npts
;
306
npts
=5;
307
308
309
310
// Output boundaries
311
//------------------
312
snprintf
(
filename
,
sizeof
(
filename
),
"%s/boundaries.dat"
,
doc_info
.directory().c_str());
313
some_file
.open(
filename
);
314
mesh_pt()->output_boundaries(
some_file
);
315
some_file
.close();
316
317
318
// Output solution
319
//----------------
320
snprintf
(
filename
,
sizeof
(
filename
),
"%s/soln%i.dat"
,
doc_info
.directory().c_str(),
321
doc_info
.number());
322
some_file
.open(
filename
);
323
mesh_pt()->output(
some_file
,
npts
);
324
some_file
.close();
325
326
327
// Output exact solution
328
//----------------------
329
snprintf
(
filename
,
sizeof
(
filename
),
"%s/exact_soln%i.dat"
,
doc_info
.directory().c_str(),
330
doc_info
.number());
331
some_file
.open(
filename
);
332
mesh_pt()->output_fct(
some_file
,
npts
,
TanhSolnForPoisson::get_exact_u
);
333
some_file
.close();
334
335
336
// Doc error
337
//----------
338
double
error
,
norm
;
339
snprintf
(
filename
,
sizeof
(
filename
),
"%s/error%i.dat"
,
doc_info
.directory().c_str(),
340
doc_info
.number());
341
some_file
.open(
filename
);
342
mesh_pt()->compute_error(
some_file
,
TanhSolnForPoisson::get_exact_u
,
343
error
,
norm
);
344
some_file
.close();
345
cout
<<
"error: "
<<
sqrt
(
error
) << std::endl;
346
cout
<<
"norm : "
<<
sqrt
(
norm
) << std::endl << std::endl;
347
348
349
// Get norm of solution
350
snprintf
(
filename
,
sizeof
(
filename
),
"%s/norm%i.dat"
,
doc_info
.directory().c_str(),
351
doc_info
.number());
352
some_file
.open(
filename
);
353
double
norm_soln
=0.0;
354
mesh_pt()->compute_norm(
norm_soln
);
355
some_file
<<
sqrt
(
norm_soln
) << std::endl;
356
cout
<<
"Norm of computed solution: "
<<
sqrt
(
norm_soln
) <<
endl
;
357
some_file
.close();
358
359
}
360
361
362
363
364
365
//========================================================================
366
/// Demonstrate how to solve Poisson problem
367
//========================================================================
368
int
main
(
int
argc
,
char
*
argv
[])
369
{
370
// Store command line arguments
371
CommandLineArgs::setup(
argc
,
argv
);
372
373
// Check number of command line arguments: Need exactly two.
374
if
(
argc
!=4)
375
{
376
std::string
error_message
=
377
"Wrong number of command line arguments.\n"
;
378
error_message
+=
379
"Must specify the following file names \n"
;
380
error_message
+=
381
"filename.node then filename.ele then filename.poly\n"
;
382
383
throw
OomphLibError
(
error_message
,
384
OOMPH_CURRENT_FUNCTION
,
385
OOMPH_EXCEPTION_LOCATION
);
386
}
387
388
#ifdef ADAPT
389
const
unsigned
max_adapt
= 3;
390
#endif
391
392
// Convert arguments to strings that specify the input file names
393
string
node_file_name
(
argv
[1]);
394
string
element_file_name
(
argv
[2]);
395
string
poly_file_name
(
argv
[3]);
396
397
// Label for output
398
DocInfo
doc_info
;
399
400
// Output directory
401
doc_info
.set_directory(
"RESLT"
);
402
403
#ifndef ADAPT
404
// Do the problem with cubic elements
405
//-----------------------------------
406
{
407
cout
<< std::endl <<
"Cubic elements"
<< std::endl;
408
cout
<<
"=============="
<< std::endl << std::endl;
409
410
//Set up the problem
411
PoissonProblem<TPoissonElement<2,4>
>
412
problem
(&
TanhSolnForPoisson::get_source
,
node_file_name
,
413
element_file_name
,
414
poly_file_name
);
415
416
// Solve the problem
417
problem
.newton_solve();
418
419
//Output solution
420
problem
.doc_solution(
doc_info
);
421
422
//Increment counter for solutions
423
doc_info
.number()++;
424
}
425
#endif
426
427
428
// Do the problem with quadratic elements
429
//---------------------------------------
430
{
431
cout
<< std::endl <<
"Quadratic elements"
<< std::endl;
432
cout
<<
"==================="
<< std::endl << std::endl;
433
434
#ifdef ADAPT
435
//Set up the problem
436
PoissonProblem<ProjectablePoissonElement<TPoissonElement<2,3>
> >
437
problem
(&
TanhSolnForPoisson::get_source
,
node_file_name
,
438
element_file_name
,
439
poly_file_name
);
440
#else
441
//Set up the problem
442
PoissonProblem<TPoissonElement<2,3>
>
443
problem
(&
TanhSolnForPoisson::get_source
,
node_file_name
,
444
element_file_name
,
445
poly_file_name
);
446
#endif
447
448
#ifdef ADAPT
449
// Solve the problem with adaptation
450
problem
.newton_solve(
max_adapt
);
451
#else
452
// Solve the problem
453
problem
.newton_solve();
454
#endif
455
456
//Output solution
457
problem
.doc_solution(
doc_info
);
458
459
//Increment counter for solutions
460
doc_info
.number()++;
461
}
462
463
464
465
// Do the problem with linear elements
466
//------------------------------------
467
{
468
cout
<< std::endl <<
"Linear elements"
<< std::endl;
469
cout
<<
"==============="
<< std::endl << std::endl;
470
471
#ifdef ADAPT
472
//Set up the problem
473
PoissonProblem<ProjectablePoissonElement<TPoissonElement<2,2>
> >
474
problem
(&
TanhSolnForPoisson::get_source
,
node_file_name
,
475
element_file_name
,
476
poly_file_name
);
477
#else
478
//Set up the problem
479
PoissonProblem<TPoissonElement<2,2>
>
480
problem
(&
TanhSolnForPoisson::get_source
,
node_file_name
,
481
element_file_name
,
482
poly_file_name
);
483
#endif
484
485
#ifdef ADAPT
486
// Solve the problem with adaptation
487
problem
.newton_solve(
max_adapt
);
488
#else
489
// Solve the problem
490
problem
.newton_solve();
491
#endif
492
493
//Output solution
494
problem
.doc_solution(
doc_info
);
495
496
//Increment counter for solutions
497
doc_info
.number()++;
498
}
499
500
}
501
502
503
PoissonProblem
Micky mouse Poisson problem.
Definition
mesh_from_triangle_poisson.cc:94
PoissonProblem::mesh_pt
TriangleMesh< ELEMENT > * mesh_pt()
Access function for the specific mesh.
Definition
mesh_from_triangle_poisson.cc:148
PoissonProblem::mesh_pt
RefineableTriangleMesh< ELEMENT > * mesh_pt()
Access function for the specific mesh.
Definition
mesh_from_triangle_poisson.cc:142
PoissonProblem::Source_fct_pt
PoissonEquations< 2 >::PoissonSourceFctPt Source_fct_pt
Pointer to source function.
Definition
mesh_from_triangle_poisson.cc:160
PoissonProblem::actions_before_newton_solve
void actions_before_newton_solve()
Update the problem specs before solve: (Re)set boundary conditions.
Definition
mesh_from_triangle_poisson.cc:110
PoissonProblem::Bulk_mesh_pt
RefineableTriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the bulk mesh.
Definition
mesh_from_triangle_poisson.cc:164
PoissonProblem::error_estimator_pt
Z2ErrorEstimator * error_estimator_pt
Error estimator.
Definition
mesh_from_triangle_poisson.cc:167
PoissonProblem::actions_after_newton_solve
void actions_after_newton_solve()
Update the problem specs before solve (empty)
Definition
mesh_from_triangle_poisson.cc:132
PoissonProblem::PoissonProblem
PoissonProblem(PoissonEquations< 2 >::PoissonSourceFctPt source_fct_pt, const string &node_file_name, const string &element_file_name, const string &poly_file_name)
Constructor: Pass pointer to source function and names of two triangle input files.
Definition
mesh_from_triangle_poisson.cc:179
PoissonProblem::doc_solution
void doc_solution(DocInfo &doc_info)
Doc the solution.
Definition
mesh_from_triangle_poisson.cc:298
PoissonProblem::~PoissonProblem
~PoissonProblem()
Destructor (empty)
Definition
mesh_from_triangle_poisson.cc:107
PoissonProblem::actions_after_adapt
void actions_after_adapt()
Actions performed after the adaptation steps.
Definition
mesh_from_triangle_poisson.cc:253
oomph::RefineableTriangleMesh
Unstructured refineable Triangle Mesh.
Definition
triangle_mesh.h:2225
main
int main(int argc, char *argv[])
Driver for FlowPastBox test problem.
Definition
mesh_from_triangle_navier_stokes.cc:339
TanhSolnForPoisson
Namespace for exact solution for Poisson equation with sharp step.
Definition
mesh_from_triangle_poisson.cc:47
TanhSolnForPoisson::Beta
double Beta
Parameter for angle of step.
Definition
mesh_from_triangle_poisson.cc:53
TanhSolnForPoisson::get_source
void get_source(const Vector< double > &x, double &source)
Source function to make it an exact solution.
Definition
mesh_from_triangle_poisson.cc:71
TanhSolnForPoisson::Alpha
double Alpha
Parameter for steepness of step.
Definition
mesh_from_triangle_poisson.cc:50
TanhSolnForPoisson::get_exact_u
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
Definition
mesh_from_triangle_poisson.cc:57
oomph
Definition
annular_domain.h:35
triangle_mesh.h