CppNoddy  0.92
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
CppNoddy::TwoD_Mapped_Node_Mesh< _Type > Class Template Reference

A two dimensional (mapped) mesh utility object. More...

#include <TwoD_Mapped_Node_Mesh.h>

Inheritance diagram for CppNoddy::TwoD_Mapped_Node_Mesh< _Type >:
Example::My_Mesh< _Type > Example::My_Mesh< _Type >

Public Member Functions

virtual double FnComp_X (const double &zeta) const
 Mapping function that provides a computational X coordinate from a physical coordinate. More...
 
virtual double FnComp_Xd (const double &zeta) const
 Mapping function that provides the first derivative of the computational m_X coordinate as a function of the physical coordinate. More...
 
virtual double FnComp_Xdd (const double &zeta) const
 Mapping function that provides the second derivative of the computational X coordinate as a function of the physical coordinate. More...
 
virtual double FnComp_Y (const double &eta) const
 Mapping function that provides the computational Y coordinate as a function of the physical coordinate. More...
 
virtual double FnComp_Yd (const double &eta) const
 Mapping function that provides the first derivative of the computational Y coordinate as a function of the physical coordinate. More...
 
virtual double FnComp_Ydd (const double &eta) const
 Mapping function that provides the second derivative of the computational Y coordinate as a function of the physical coordinate. More...
 
 TwoD_Mapped_Node_Mesh ()
 
 TwoD_Mapped_Node_Mesh (const double left, const double right, const double bottom, const double top, const std::size_t nx, const std::size_t ny, const std::size_t nvars)
 ctor of a blank mesh More...
 
 TwoD_Mapped_Node_Mesh (std::string filename, const std::size_t nx, const std::size_t ny, const std::size_t nv)
 
void init_mapping ()
 Sometimes a useful mapping is painful to invert analytically. More...
 
std::pair< double, double > get_comp_step_sizes () const
 
virtual ~TwoD_Mapped_Node_Mesh ()
 dtor More...
 
DenseVector< _Type > operator() (const std::size_t nodex, const std::size_t nodey)
 Access operator for a nodal point that returns a vector. More...
 
_Type & operator() (const std::size_t nodex, const std::size_t nodey, const std::size_t var)
 Access operator for a nodal point/variable in the mesh. More...
 
const _Type & operator() (const std::size_t nodex, const std::size_t nodey, const std::size_t var) const
 Const access operator for a nodal point/variable in the mesh. More...
 
std::pair< double, double > coord (const std::size_t nodex, const std::size_t nodey) const
 Access the physical nodal position - as a pair. More...
 
void set_nodes_vars (const std::size_t nodex, const std::size_t nodey, const DenseVector< _Type > &U)
 Set the variables stored at A SPECIFIED node. More...
 
DenseVector< _Type > get_nodes_vars (const std::size_t nodex, const std::size_t nodey) const
 Get the variables stored at A SPECIFIED node – equivalent to mesh(nodex,nodey). More...
 
std::pair< std::size_t, std::size_t > get_nnodes () const
 Get the number of nodes in the two directions of the 2D mesh. More...
 
std::size_t get_nvars () const
 Get the number of variables that are stored at each node. More...
 
DenseVector< double > & xnodes ()
 Access the vector of x-nodal positions. More...
 
DenseVector< double > & ynodes ()
 Access the vector of y-nodal positions. More...
 
OneD_Node_Mesh< _Type > get_xsection_at_xnode (const std::size_t nodex) const
 Get a cross section of the 2D mesh at a specified (constant) x node. More...
 
OneD_Node_Mesh< _Type > get_xsection_at_ynode (const std::size_t nodey) const
 Get a cross section of the 2D mesh at a specified (constant) y node. More...
 
void read (std::string filename, const bool reset=false)
 A simple method for reading data from a file. More...
 
void dump_gnu (std::string filename) const
 A simple method for dumping data to a file for gnuplot surface plotting. More...
 
void normalise (const std::size_t &var)
 Normalise all data in the mesh based on one variable. More...
 
void scale (const _Type &value)
 Rescale all values stored in the mapped mesh by a scalar. More...
 
double max (unsigned var)
 Find the maximum stored absolute value in the mesh for a given variable – no interpolation is used. More...
 
DenseVector< _Type > get_interpolated_vars (const double &x, const double &y)
 Get a bilinearly interpolated value at a specified point. More...
 
void normalise (const std::size_t &var)
 
void normalise (const std::size_t &var)
 
void read (std::string filename, bool reset)
 
void read (std::string filename, bool reset)
 
void dump_gnu (std::string filename) const
 
void dump_gnu (std::string filename) const
 

Protected Attributes

double m_left
 
double m_right
 
double m_bottom
 
double m_top
 
std::size_t m_nx
 
std::size_t m_ny
 
std::size_t m_nv
 
DenseVector< double > m_compX
 
DenseVector< double > m_compY
 
DenseVector< double > m_X
 
DenseVector< double > m_Y
 
DenseVector< _Type > m_vars
 

Detailed Description

template<typename _Type>
class CppNoddy::TwoD_Mapped_Node_Mesh< _Type >

A two dimensional (mapped) mesh utility object.

Definition at line 21 of file TwoD_Mapped_Node_Mesh.h.

Constructor & Destructor Documentation

◆ TwoD_Mapped_Node_Mesh() [1/3]

template<typename _Type >
CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::TwoD_Mapped_Node_Mesh ( )
inline

Definition at line 90 of file TwoD_Mapped_Node_Mesh.h.

91 {}

◆ TwoD_Mapped_Node_Mesh() [2/3]

template<typename _Type >
CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::TwoD_Mapped_Node_Mesh ( const double  left,
const double  right,
const double  bottom,
const double  top,
const std::size_t  nx,
const std::size_t  ny,
const std::size_t  nvars 
)
inline

ctor of a blank mesh

Definition at line 94 of file TwoD_Mapped_Node_Mesh.h.

94 : m_left(left), m_right(right), m_bottom(bottom), m_top(top), m_nx(nx), m_ny(ny), m_nv(nvars) {
95#ifdef DEBUG
96 std::cout << "[DEBUG] making a blank TwoD_Mapped_Node_Mesh\n";
97 std::cout << "[DEBUG] a " << nx << " by " << ny << " mesh with " << nvars << " variables.\n";
98#endif
99 // initialise the storage, but fill these below in init_mapping()
100 m_X = DenseVector<double> (m_nx,0.0);
101 m_Y = DenseVector<double> (m_ny,0.0);
102 m_compX = DenseVector<double> (m_nx,0.0);
103 m_compY = DenseVector<double> (m_ny,0.0);
104 // we'll store the data as ( x, y, v ) -> x * ny * nv + y * nv + v
105 m_vars = DenseVector<_Type>(m_nx * m_ny * m_nv, 0.0);
106 // user needs to call init_mapping() to set up m_compX, m_compY after this
107#ifdef DEBUG
108 std::cout << "[DEBUG] leaving the TwoD_Mapped_Node_Mesh ctor.\n";
109#endif
110 }

References CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compX, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compY, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nv, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nx, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_ny, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_vars, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_X, and CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_Y.

◆ TwoD_Mapped_Node_Mesh() [3/3]

template<typename _Type >
CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::TwoD_Mapped_Node_Mesh ( std::string  filename,
const std::size_t  nx,
const std::size_t  ny,
const std::size_t  nv 
)
inline

Definition at line 114 of file TwoD_Mapped_Node_Mesh.h.

114 :
115 m_nx(nx), m_ny(ny), m_nv(nv) {
116 // need storage for the coordinates
117 m_X = DenseVector<double>(m_nx, 0.0);
118 m_Y = DenseVector<double>(m_ny, 0.0);
119 m_compX = DenseVector<double> (m_nx,0.0);
120 m_compY = DenseVector<double> (m_ny,0.0);
121 // we'll store the data as ( x, y, v ) -> x * ny * nv + y * nv + v
122 m_vars = DenseVector<_Type>(m_nx * m_ny * m_nv, 0.0);
123 // now read the mesh from the given filename -- this also updates the m_X data
124 read(filename, true);
125 // set up the private member data on the box size.
126 m_left = m_X[0];
127 m_right = m_X[m_nx-1];
128 m_bottom = m_Y[0];
129 m_top = m_Y[m_ny-1];
130 // user needs to call init_mapping() to set up m_compX, m_compY after this
131 }
void read(std::string filename, const bool reset=false)
A simple method for reading data from a file.

References CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_bottom, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compX, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compY, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_left, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nv, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nx, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_ny, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_right, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_top, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_vars, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_X, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_Y, and CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::read().

◆ ~TwoD_Mapped_Node_Mesh()

template<typename _Type >
virtual CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::~TwoD_Mapped_Node_Mesh ( )
inlinevirtual

dtor

Definition at line 143 of file TwoD_Mapped_Node_Mesh.h.

144 {}

Member Function Documentation

◆ coord()

template<typename _Type >
std::pair< double, double > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::coord ( const std::size_t  nodex,
const std::size_t  nodey 
) const
inline

Access the physical nodal position - as a pair.

Parameters
nodexThe x nodal position to return
nodeyThe y nodal position to return
Returns
The spatial position of this node as a pair

Definition at line 371 of file TwoD_Mapped_Node_Mesh.h.

371 {
372#ifdef PARANOID
373 if(nodex > m_nx - 1 || nodey > m_ny - 1) {
374 std::string problem;
375 problem = " The TwoD_Mapped_Node_Mesh.coord method is trying to \n";
376 problem += " access a nodal point that is not in the mesh. \n";
377 throw ExceptionRange(problem, m_nx, nodex, m_ny, nodey);
378 }
379#endif
380 std::pair< double, double > pos;
381 pos.first = m_X[ nodex ];
382 pos.second = m_Y[ nodey ];
383 return pos;
384 }

Referenced by main().

◆ dump_gnu() [1/3]

template<typename _Type >
void CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::dump_gnu ( std::string  filename) const

A simple method for dumping data to a file for gnuplot surface plotting.

Parameters
filenameThe filename to write the data to (will overwrite)

Referenced by main().

◆ dump_gnu() [2/3]

void CppNoddy::TwoD_Mapped_Node_Mesh< double >::dump_gnu ( std::string  filename) const

Definition at line 298 of file TwoD_Mapped_Node_Mesh.cpp.

298 {
299 std::ofstream dump;
300 dump.open(filename.c_str());
301 dump.precision(15);
302 dump.setf(std::ios::showpoint);
303 dump.setf(std::ios::showpos);
304 dump.setf(std::ios::scientific);
305 //
306 for(std::size_t i = 0; i < m_nx; ++i) {
307 for(std::size_t j = 0; j < m_ny; ++j) {
308 dump << m_X[ i ] << " " << m_Y[ j ] << " ";
309 for(std::size_t var = 0; var < m_nv; ++var) {
310 dump << m_vars[(i * m_ny + j) * m_nv + var ] << " ";
311 }
312 dump << "\n";
313 }
314 dump << "\n";
315 }
316 dump.close();
317 }

◆ dump_gnu() [3/3]

void CppNoddy::TwoD_Mapped_Node_Mesh< D_complex >::dump_gnu ( std::string  filename) const

Definition at line 320 of file TwoD_Mapped_Node_Mesh.cpp.

320 {
321 std::ofstream dump;
322 dump.open(filename.c_str());
323 dump.precision(15);
324 dump.setf(std::ios::showpoint);
325 dump.setf(std::ios::showpos);
326 dump.setf(std::ios::scientific);
327 //
328 for(std::size_t i = 0; i < m_nx; ++i) {
329 for(std::size_t j = 0; j < m_ny; ++j) {
330 dump << m_X[ i ] << " " << m_Y[ j ] << " ";
331 for(std::size_t var = 0; var < m_nv; ++var) {
332 dump << real(m_vars[(i * m_ny + j) * m_nv + var ]) << " ";
333 dump << imag(m_vars[(i * m_ny + j) * m_nv + var ]) << " ";
334 }
335 dump << "\n";
336 }
337 dump << "\n";
338 }
339 dump.close();
340 }
DenseVector< double > real(const DenseVector< D_complex > &X)
Return a double DENSE vector containing the real part of a complex DENSE vector.
Definition: Utility.cpp:177
DenseVector< double > imag(const DenseVector< D_complex > &X)
Return a double DENSE vector containing the imaginary part of a complex DENSE vector.
Definition: Utility.cpp:185

◆ FnComp_X()

template<typename _Type >
virtual double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_X ( const double &  zeta) const
inlinevirtual

Mapping function that provides a computational X coordinate from a physical coordinate.

Parameters
zetaThe physical coordinate
Returns
The corresponding computational coordinate

Reimplemented in Example::My_Mesh< _Type >, and Example::My_Mesh< _Type >.

Definition at line 28 of file TwoD_Mapped_Node_Mesh.h.

28 {
29 std::string problem;
30 problem = "The TwoD_Mapped_Node_Mesh::FnComp_X method has not been implemented.\n";
31 problem += "You have to implement this method to define the mesh.\n";
32 throw ExceptionRuntime(problem);
33 }

Referenced by CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_interpolated_vars().

◆ FnComp_Xd()

template<typename _Type >
virtual double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_Xd ( const double &  zeta) const
inlinevirtual

Mapping function that provides the first derivative of the computational m_X coordinate as a function of the physical coordinate.

Parameters
zetaThe physical coordinate
Returns
The corresponding derivative of the computational coordinate

Reimplemented in Example::My_Mesh< _Type >, and Example::My_Mesh< _Type >.

Definition at line 39 of file TwoD_Mapped_Node_Mesh.h.

39 {
40 std::string problem;
41 problem = "The TwoD_Mapped_Node_Mesh::FnComp_Xd method has not been implemented.\n";
42 problem += "You have to implement this method to define the mesh.\n";
43 throw ExceptionRuntime(problem);
44 }

◆ FnComp_Xdd()

template<typename _Type >
virtual double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_Xdd ( const double &  zeta) const
inlinevirtual

Mapping function that provides the second derivative of the computational X coordinate as a function of the physical coordinate.

Parameters
zetaThe physical coordinate
Returns
The corresponding derivative of the computational coordinate

Reimplemented in Example::My_Mesh< _Type >, and Example::My_Mesh< _Type >.

Definition at line 50 of file TwoD_Mapped_Node_Mesh.h.

50 {
51 std::string problem;
52 problem = "The TwoD_Mapped_Node_Mesh::FnComp_Xdd method has not been implemented.\n";
53 problem += "You have to implement this method to define the mesh.\n";
54 throw ExceptionRuntime(problem);
55 }

◆ FnComp_Y()

template<typename _Type >
virtual double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_Y ( const double &  eta) const
inlinevirtual

Mapping function that provides the computational Y coordinate as a function of the physical coordinate.

Parameters
etaThe physical coordinate
Returns
The corresponding derivative of the computational coordinate

Reimplemented in Example::My_Mesh< _Type >, and Example::My_Mesh< _Type >.

Definition at line 61 of file TwoD_Mapped_Node_Mesh.h.

61 {
62 std::string problem;
63 problem = "The TwoD_Mapped_Node_Mesh::FnComp_Y method has not been implemented.\n";
64 problem += "You have to implement this method to define the mesh.\n";
65 throw ExceptionRuntime(problem);
66 }

Referenced by CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_interpolated_vars().

◆ FnComp_Yd()

template<typename _Type >
virtual double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_Yd ( const double &  eta) const
inlinevirtual

Mapping function that provides the first derivative of the computational Y coordinate as a function of the physical coordinate.

Parameters
etaThe physical coordinate
Returns
The corresponding derivative of the computational coordinate

Reimplemented in Example::My_Mesh< _Type >, and Example::My_Mesh< _Type >.

Definition at line 72 of file TwoD_Mapped_Node_Mesh.h.

72 {
73 std::string problem;
74 problem = "The TwoD_Mapped_Node_Mesh::FnComp_Yd method has not been implemented.\n";
75 problem += "You have to implement this method to define the mesh.\n";
76 throw ExceptionRuntime(problem);
77 }

◆ FnComp_Ydd()

template<typename _Type >
virtual double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_Ydd ( const double &  eta) const
inlinevirtual

Mapping function that provides the second derivative of the computational Y coordinate as a function of the physical coordinate.

Parameters
etaThe physical coordinate
Returns
The corresponding derivative of the computational coordinate

Reimplemented in Example::My_Mesh< _Type >, and Example::My_Mesh< _Type >.

Definition at line 83 of file TwoD_Mapped_Node_Mesh.h.

83 {
84 std::string problem;
85 problem = "The TwoD_Mapped_Node_Mesh::FnComp_Ydd method has not been implemented.\n";
86 problem += "You have to implement this method to define the mesh.\n";
87 throw ExceptionRuntime(problem);
88 }

◆ get_comp_step_sizes()

template<typename _Type >
std::pair< double, double > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_comp_step_sizes

Definition at line 155 of file TwoD_Mapped_Node_Mesh.cpp.

155 {
156 std::pair<double,double> steps;
157 steps.first = m_compX[1]-m_compX[0];
158 steps.second = m_compY[1]-m_compY[0];
159 return steps;
160 }

Referenced by main().

◆ get_interpolated_vars()

template<typename _Type >
DenseVector< _Type > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_interpolated_vars ( const double &  x,
const double &  y 
)
inline

Get a bilinearly interpolated value at a specified point.

Parameters
xPhysical (unmapped) x-coordinate in the 2D mesh
yPhysical (unmapped) y-coordinate in the 2D mesh
Returns
A vector of bilinearly interpolated values

Definition at line 251 of file TwoD_Mapped_Node_Mesh.h.

252 {
253 double compX = FnComp_X(x);
254 double compY = FnComp_Y(y);
255 const double tol(1.e-8);
256 // check start and end
257 if ((compX < m_compX[0] - tol) || (compX > m_compX[m_nx-1] + tol)) {
258 std::string problem;
259 problem = " The TwoD_Node_Mesh.get_interpolated_vars method has been called with \n";
260 problem += " an x coordinate that lies outside the mesh. \n";
261 throw ExceptionRuntime(problem);
262 }
263 // check start and end
264 if ((compY < m_compY[0] - tol) || (compY > m_compY[m_ny-1] + tol)) {
265 std::string problem;
266 problem = " The TwoD_Node_Mesh.get_interpolated_vars method has been called with \n";
267 problem += " a y coordinate that lies outside the mesh. \n";
268 throw ExceptionRuntime(problem);
269 }
270 int bottom_j(-1);
271 for (unsigned j = 0; j < m_ny-1; ++j) {
272 if ((compY >= m_compY[j] - tol) && (compY <= m_compY[j+1] + tol)) {
273 bottom_j = j;
274 }
275 }
276 if (bottom_j == -1) {
277 std::string problem;
278 problem = " The TwoD_Node_Mesh.get_interpolated_vars method is broken.\n";
279 throw ExceptionRuntime(problem);
280 }
281 //
282 OneD_Node_Mesh<_Type> bottom_row( m_compX, m_nv );
283 OneD_Node_Mesh<_Type> top_row( m_compX, m_nv );
284 for (std::size_t i = 0; i < m_nx; ++i ) {
285 bottom_row.set_nodes_vars(i, this -> get_nodes_vars(i, bottom_j));
286 top_row.set_nodes_vars(i, this -> get_nodes_vars(i, bottom_j+1));
287 }
288 const double compY1 = m_compY[ bottom_j ];
289 const double compY2 = m_compY[ bottom_j+1 ];
290 DenseVector<_Type> result =
291 top_row.get_interpolated_vars(compX)*(compY-compY1)/(compY2-compY1)
292 + bottom_row.get_interpolated_vars(compX)*(compY2-compY)/(compY2-compY1);
293 return result;
294 }
virtual double FnComp_X(const double &zeta) const
Mapping function that provides a computational X coordinate from a physical coordinate.
DenseVector< _Type > get_nodes_vars(const std::size_t nodex, const std::size_t nodey) const
Get the variables stored at A SPECIFIED node – equivalent to mesh(nodex,nodey).
virtual double FnComp_Y(const double &eta) const
Mapping function that provides the computational Y coordinate as a function of the physical coordinat...

References CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_X(), CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::FnComp_Y(), CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_interpolated_vars(), CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_nodes_vars(), CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compX, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compY, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nv, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nx, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_ny, and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_nodes_vars().

◆ get_nnodes()

template<typename _Type >
std::pair< std::size_t, std::size_t > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_nnodes

Get the number of nodes in the two directions of the 2D mesh.

Returns
A pair consisting of the number of nodes in the 2 directions

Definition at line 163 of file TwoD_Mapped_Node_Mesh.cpp.

163 {
164 std::pair< std::size_t, std::size_t > nodes;
165 nodes.first = m_nx;
166 nodes.second = m_ny;
167 return nodes;
168 }

◆ get_nodes_vars()

template<typename _Type >
DenseVector< _Type > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_nodes_vars ( const std::size_t  nodex,
const std::size_t  nodey 
) const

Get the variables stored at A SPECIFIED node – equivalent to mesh(nodex,nodey).

Parameters
nodexThe x nodal index to be returned
nodeyThe y nodal index to be returned
Returns
The vector of VARIABLES stored at this nodal point

Definition at line 37 of file TwoD_Mapped_Node_Mesh.cpp.

37 {
38#ifdef PARANOID
39 if(nodex > m_nx - 1 || nodey > m_ny - 1) {
40 std::string problem;
41 problem = " The TwoD_Mapped_Node_Mesh.get_nodes_vars method is trying to \n";
42 problem += " access a nodal point that is not in the mesh. \n";
43 throw ExceptionRange(problem, m_nx, nodex, m_ny, nodey);
44 }
45#endif
46 // construct a vector with m_nv elements starting from a pointer
47 DenseVector<_Type> nodes_vars(m_nv, &m_vars[(nodex * m_ny + nodey) * m_nv ]);
48 return nodes_vars;
49 }

Referenced by CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_interpolated_vars().

◆ get_nvars()

template<typename _Type >
std::size_t CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_nvars

Get the number of variables that are stored at each node.

Returns
The number of variables that have data stored at each nodal point

Definition at line 171 of file TwoD_Mapped_Node_Mesh.cpp.

171 {
172 return m_nv;
173 }

◆ get_xsection_at_xnode()

template<typename _Type >
OneD_Node_Mesh< _Type > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_xsection_at_xnode ( const std::size_t  nodex) const

Get a cross section of the 2D mesh at a specified (constant) x node.

Parameters
nodexThe x nodal index at which the cross section is to be taken
Returns
A 1D nodal mesh

Definition at line 52 of file TwoD_Mapped_Node_Mesh.cpp.

52 {
53 OneD_Node_Mesh<_Type> xsection(m_Y, m_nv);
54 for(std::size_t nodey = 0; nodey < m_ny; ++nodey) {
55 xsection.set_nodes_vars(nodey, this -> get_nodes_vars(nodex, nodey));
56 }
57 return xsection;
58 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_nodes_vars().

◆ get_xsection_at_ynode()

template<typename _Type >
OneD_Node_Mesh< _Type > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_xsection_at_ynode ( const std::size_t  nodey) const

Get a cross section of the 2D mesh at a specified (constant) y node.

Parameters
nodeyThe y nodal index at which the cross section is to be taken
Returns
A 1D nodal mesh

Definition at line 61 of file TwoD_Mapped_Node_Mesh.cpp.

61 {
62 OneD_Node_Mesh<_Type> xsection(m_X, m_nv);
63 for(std::size_t nodex = 0; nodex < m_nx; ++nodex) {
64 xsection.set_nodes_vars(nodex, this -> get_nodes_vars(nodex, nodey));
65 }
66 return xsection;
67 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_nodes_vars().

◆ init_mapping()

template<typename _Type >
void CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::init_mapping

Sometimes a useful mapping is painful to invert analytically.

Here we construct the physical node distribution by finding the set of points x_i such that the computation nodes X satisfy X_i = FnComp_X(x_i), for example.

Definition at line 75 of file TwoD_Mapped_Node_Mesh.cpp.

75 {
76 #ifdef DEBUG
77 std::cout << "[DEBUG] Physical domain is [" << m_left << "," << m_right
78 << "] x [" << m_bottom << "," << m_top << "]\n";
79 std::cout << "[DEBUG] Computational domain is ["
80 << FnComp_X(m_left) << "," << FnComp_X(m_right)
81 << "] x [" << FnComp_Y(m_bottom) << "," << FnComp_Y(m_top) << "]\n";
82 #endif
83 {
84 // a uniform mesh in the computational coordinates
85 double comp_left(FnComp_X(m_left));
86 double comp_right(FnComp_X(m_right));
87 const double comp_delta_x = (comp_right - comp_left) / (m_nx - 1);
88 for(std::size_t i = 0; i < m_nx; ++i) {
89 m_compX[i] = comp_left + comp_delta_x * i;
90 }
91 }
92 {
93 // a uniform mesh in the computational coordinates
94 double comp_bottom(FnComp_Y(m_bottom));
95 double comp_top(FnComp_Y(m_top));
96 const double comp_delta_y = (comp_top - comp_bottom) / (m_ny - 1);
97 for(std::size_t j = 0; j < m_ny; ++j) {
98 m_compY[j] = comp_bottom + comp_delta_y * j;
99 }
100 }
101
102 // temporary fill to span the domain -- corrected below
105
106 // we now need the corresponding m_Y coordinates in the physical domain
107 {
108 // start and end points are mapped correctly
109 for(unsigned j = 1; j < m_ny-1; ++j) {
110 // // for each node in m_compY
111 // unsigned kmin(0);
112 // double min(99e9);
113 // for(unsigned k = 0; k < m_ny; ++k) {
114 // // find the y value that is closest to it
115 // if(std::abs(FnComp_Y(m_Y[k]) - m_compY[j]) < min) {
116 // min = std::abs(FnComp_Y(m_Y[k]) - m_compY[j]);
117 // kmin = k;
118 // }
119 // }
120 double y = m_Y[j-1];
121 double delta = 1.e-8;
122 double correction = 1.0;
123 do {
124 double newY = FnComp_Y(y + delta) - m_compY[j];
125 double oldY = FnComp_Y(y) - m_compY[j];
126 double deriv = (newY-oldY)/delta;
127 correction = -oldY/deriv;
128 y += correction;
129 } while(fabs(correction) > 1.e-8);
130 m_Y[ j ] = y;
131 }
132 }
133
134 // we now need the corresponding m_X coordinates in the physical domain
135 {
136 // start and end are already mapped correctly
137 for(unsigned i = 1; i < m_nx-1; ++i) {
138 // use the previous point as an approximate guess at current mapping
139 double x = m_X[i-1];
140 double delta = 1.e-8;
141 double correction = 1.0;
142 do {
143 double newX = FnComp_X(x + delta) - m_compX[i];
144 double oldX = FnComp_X(x) - m_compX[i];
145 double deriv = (newX-oldX)/delta;
146 correction = -oldX/deriv;
147 x += correction;
148 } while(fabs(correction) > 1.e-8);
149 m_X[ i ] = x;
150 }
151 }
152 }
const double delta(0.5)
DenseVector< double > uniform_node_vector(const double &lower, const double &upper, const std::size_t &N)
Return a DENSE vector with the nodal points of a uniform mesh distributed between the upper/lower bou...
Definition: Utility.cpp:113

References CppNoddy::Utility::uniform_node_vector().

Referenced by Example::My_Mesh< _Type >::My_Mesh().

◆ max()

template<typename _Type >
double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::max ( unsigned  var)
inline

Find the maximum stored absolute value in the mesh for a given variable – no interpolation is used.

Parameters
varThe variable index whose maximum is being asked for
Returns
The value of the maximum (abs value)

Definition at line 234 of file TwoD_Mapped_Node_Mesh.h.

234 {
235 double max(0.0);
236 // step through the nodes
237 for(unsigned nodex = 0; nodex < m_nx; ++nodex) {
238 for(unsigned nodey = 0; nodey < m_ny; ++nodey) {
239 if(std::abs(m_vars[(nodex * m_ny + nodey) * m_nv + var ]) > max) {
240 max = std::abs(m_vars[(nodex * m_ny + nodey) * m_nv + var ]);
241 }
242 }
243 }
244 return max;
245 }
double max(unsigned var)
Find the maximum stored absolute value in the mesh for a given variable – no interpolation is used.

References CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nv, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nx, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_ny, CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_vars, and CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::max().

Referenced by CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::max().

◆ normalise() [1/3]

template<typename _Type >
void CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::normalise ( const std::size_t &  var)

Normalise all data in the mesh based on one variable.

Parameters
varThis var will have its peak (absolute) value as +/-unity following the normalisation. All other variables will also be rescaled by the same amount.

◆ normalise() [2/3]

void CppNoddy::TwoD_Mapped_Node_Mesh< double >::normalise ( const std::size_t &  var)

Definition at line 186 of file TwoD_Mapped_Node_Mesh.cpp.

186 {
187 double maxval(max(var));
188 m_vars.scale(1./maxval);
189 }

◆ normalise() [3/3]

void CppNoddy::TwoD_Mapped_Node_Mesh< D_complex >::normalise ( const std::size_t &  var)

Definition at line 192 of file TwoD_Mapped_Node_Mesh.cpp.

192 {
193 //std::cout << "[DEBUG] asked to normalise a complex mesh\n";
194 unsigned max_nx(0);
195 unsigned max_ny(0);
196 double max(0.0);
197 // step through the nodes
198 for(unsigned nodex = 0; nodex < m_X.size(); ++nodex) {
199 for(unsigned nodey = 0; nodey < m_Y.size(); ++nodey) {
200 if(std::abs(m_vars[(nodex * m_ny + nodey) * m_nv + var ]) > max) {
201 max = std::abs(m_vars[(nodex * m_ny + nodey) * m_nv + var ]);
202 max_nx = nodex;
203 max_ny = nodey;
204 }
205 }
206 }
207 D_complex factor(m_vars[(max_nx * m_ny + max_ny) * m_nv + var ]);
208 //std::cout << "[DEBUG] MAX |variable| had complex value of " << factor << "\n";
209 m_vars.scale(1./factor);
210 }
std::size_t size() const
A pass-thru definition to get the size of the vector.
Definition: DenseVector.h:330
std::complex< double > D_complex
A complex double precision number using std::complex.
Definition: Types.h:98

◆ operator()() [1/3]

template<typename _Type >
DenseVector< _Type > CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::operator() ( const std::size_t  nodex,
const std::size_t  nodey 
)
inline

Access operator for a nodal point that returns a vector.

Parameters
nodexThe nodal index value in the first direction
nodeyThe nodal index value in the second direction
Returns
The vector of variables stored at the node

Definition at line 319 of file TwoD_Mapped_Node_Mesh.h.

319 {
320#ifdef PARANOID
321 if(nodex > m_nx - 1 || nodey > m_ny - 1) {
322 std::string problem;
323 problem = " The TwoD_Mapped_Node_Mesh.operator() method is trying to \n";
324 problem += " access a nodal point that is not in the mesh. \n";
325 throw ExceptionRange(problem, m_nx, nodex, m_ny, nodey);
326 }
327#endif
328 return get_nodes_vars(nodex,nodey);
329 }

◆ operator()() [2/3]

template<typename _Type >
_Type & CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::operator() ( const std::size_t  nodex,
const std::size_t  nodey,
const std::size_t  var 
)
inline

Access operator for a nodal point/variable in the mesh.

Parameters
nodexThe nodal index value in the first direction
nodeyThe nodal index value in the second direction
varThe variable index to be accessed

Definition at line 333 of file TwoD_Mapped_Node_Mesh.h.

333 {
334#ifdef PARANOID
335 if(nodex > m_nx - 1 || nodey > m_ny - 1) {
336 std::string problem;
337 problem = " The TwoD_Mapped_Node_Mesh.operator() method is trying to \n";
338 problem += " access a nodal point that is not in the mesh. \n";
339 throw ExceptionRange(problem, m_nx, nodex, m_ny, nodey);
340 }
341 if(var > m_nv - 1) {
342 std::string problem;
343 problem = " The TwoD_Mapped_Node_Mesh.operator() method is trying to \n";
344 problem += " access a variable index that is not in the mesh. \n";
345 throw ExceptionRange(problem, m_nv, var);
346 }
347#endif
348 return m_vars[(nodex * m_ny + nodey) * m_nv + var ];
349 }

◆ operator()() [3/3]

template<typename _Type >
const _Type & CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::operator() ( const std::size_t  nodex,
const std::size_t  nodey,
const std::size_t  var 
) const
inline

Const access operator for a nodal point/variable in the mesh.

Parameters
nodexThe nodal index value in the first direction
nodeyThe nodal index value in the second direction
varThe variable index to be accessed

Definition at line 352 of file TwoD_Mapped_Node_Mesh.h.

352 {
353#ifdef PARANOID
354 if(nodex > m_nx - 1 || nodey > m_ny - 1) {
355 std::string problem;
356 problem = " The TwoD_Mapped_Node_Mesh.operator() method is trying to \n";
357 problem += " access a nodal point that is not in the mesh. \n";
358 throw ExceptionRange(problem, m_nx, nodex, m_ny, nodey);
359 }
360 if(var > m_nv - 1) {
361 std::string problem;
362 problem = " The TwoD_Mapped_Node_Mesh.operator() method is trying to \n";
363 problem += " access a variable index that is not in the mesh. \n";
364 throw ExceptionRange(problem, m_nv, var);
365 }
366#endif
367 return m_vars[(nodex * m_ny + nodey) * m_nv + var ];
368 }

◆ read() [1/3]

void CppNoddy::TwoD_Mapped_Node_Mesh< double >::read ( std::string  filename,
bool  reset 
)

Definition at line 213 of file TwoD_Mapped_Node_Mesh.cpp.

213 {
214 std::ifstream dump;
215 dump.open(filename.c_str());
216 if(dump.good() != true) {
217 std::string problem;
218 problem = " The TwoD_Node_Mesh.read method is trying to read a \n";
219 problem += " file (" + filename + ") that doesn't exist.\n";
220 throw ExceptionRuntime(problem);
221 }
222 dump.precision(15);
223 dump.setf(std::ios::showpoint);
224 dump.setf(std::ios::showpos);
225 dump.setf(std::ios::scientific);
226 for(std::size_t i = 0; i < m_nx; ++i) {
227 for(std::size_t j = 0; j < m_ny; ++j) {
228 double x, y;
229 dump >> x;
230 dump >> y;
231 for(std::size_t var = 0; var < m_nv; ++var) {
232 double value;
233 dump >> value;
234 m_vars[(i * m_ny + j) * m_nv + var ] = value;
235 }
236 if(reset != true) {
237 // if not reseting the mesh we should check the node positions
238 if((std::fabs(x - m_X[ i ]) > 1.e-6) || (std::fabs(y - m_Y[ j ]) > 1.e-6)) {
239 std::cout << " Read x = " << x << " Expected x = " << m_X[ i ] << "; Read y = " << y << " Expected y = " << m_Y[ j ] << " \n";
240 std::cout << " Absolute differences are " << fabs(x - m_X[i]) << " and " << fabs(y - m_Y[j]) << "\n";
241 std::string problem;
242 problem = " The TwoD_Node_Mesh.read method is trying to read a \n";
243 problem += " file whose nodal points are in a different position. \n";
244 throw ExceptionRuntime(problem);
245 }
246 } else {
247 m_X[ i ] = x;
248 m_Y[ j ] = y;
249 }
250 }
251 }
252 }

◆ read() [2/3]

void CppNoddy::TwoD_Mapped_Node_Mesh< D_complex >::read ( std::string  filename,
bool  reset 
)

Definition at line 256 of file TwoD_Mapped_Node_Mesh.cpp.

256 {
257 std::ifstream dump;
258 dump.open(filename.c_str());
259 dump.precision(15);
260 dump.setf(std::ios::showpoint);
261 dump.setf(std::ios::showpos);
262 dump.setf(std::ios::scientific);
263 //
264 // 18/06/2017: switched i and j below for consistency with double
265 //
266 for(std::size_t i = 0; i < m_nx; ++i) {
267 for(std::size_t j = 0; j < m_ny; ++j) {
268 double x, y;
269 dump >> x;
270 dump >> y;
271 for(std::size_t var = 0; var < m_nv; ++var) {
272 double value_r, value_i;
273 dump >> value_r;
274 dump >> value_i;
275 m_vars[(i * m_ny + j) * m_nv + var ] = D_complex(value_r, value_i);
276 }
277 if(reset != true) {
278 // if not reseting the mesh we should check the node positions
279 if((std::fabs(x - m_X[ i ]) > 1.e-6) || (std::fabs(y - m_Y[ j ]) > 1.e-6)) {
280 std::cout << " Read x = " << x << " Expected x = " << m_X[ i ] << "; Read y = " << y << " Expected y = " << m_Y[ j ] << " \n";
281 std::cout << " Absolute differences are " << fabs(x - m_X[i]) << " and " << fabs(y - m_Y[j]) << "\n";
282 std::string problem;
283 problem = " The TwoD_Node_Mesh.read method is trying to read a \n";
284 problem += " file whose nodal points are in a different position. \n";
285 throw ExceptionRuntime(problem);
286 }
287 } else {
288 m_X[ i ] = x;
289 m_Y[ j ] = y;
290 }
291 }
292 }
293 }

◆ read() [3/3]

template<typename _Type >
void CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::read ( std::string  filename,
const bool  reset = false 
)

A simple method for reading data from a file.

Parameters
filenameThe filename to write the data to (will overwrite)
resetWill reset the nodal positions using those from the file

Referenced by CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::TwoD_Mapped_Node_Mesh().

◆ scale()

template<typename _Type >
void CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::scale ( const _Type &  value)
inline

Rescale all values stored in the mapped mesh by a scalar.

Parameters
valueThe scalar that is to multiply all mesh content

Definition at line 226 of file TwoD_Mapped_Node_Mesh.h.

226 {
227 m_vars.scale(value);
228 }

References CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_vars.

◆ set_nodes_vars()

template<typename _Type >
void CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::set_nodes_vars ( const std::size_t  nodex,
const std::size_t  nodey,
const DenseVector< _Type > &  U 
)

Set the variables stored at A SPECIFIED node.

Parameters
nodexThe x nodal index to be set
nodeyThe y nodal index to be set
UThe vector of VARIABLES to be written to this nodal point

Definition at line 20 of file TwoD_Mapped_Node_Mesh.cpp.

20 {
21#ifdef PARANOID
22 if(U.size() > m_nv) {
23 std::string problem;
24 problem = " The TwoD_Mapped_Node_Mesh.set_nodes_vars method is trying to use a \n";
25 problem += " vector that has more entries than variables stored in the mesh. \n";
26 throw ExceptionRuntime(problem);
27 }
28#endif
29 // assign contents of U to the member data
30 std::size_t offset((nodex * m_ny + nodey) * m_nv);
31 for(std::size_t var = 0; var < m_nv; ++var) {
32 m_vars[ offset++ ] = U[ var ];
33 }
34 }
@ U
Definition: BVPKarman.cpp:20

References U.

◆ xnodes()

template<typename _Type >
DenseVector< double > & CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::xnodes

Access the vector of x-nodal positions.

Returns
A vector of the nodal positions for this mesh

Definition at line 176 of file TwoD_Mapped_Node_Mesh.cpp.

176 {
177 return m_X;
178 }

◆ ynodes()

template<typename _Type >
DenseVector< double > & CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::ynodes

Access the vector of y-nodal positions.

Returns
A vector of the nodal positions for this mesh

Definition at line 181 of file TwoD_Mapped_Node_Mesh.cpp.

181 {
182 return m_Y;
183 }

Member Data Documentation

◆ m_bottom

template<typename _Type >
double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_bottom
protected

◆ m_compX

template<typename _Type >
DenseVector<double> CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compX
protected

◆ m_compY

template<typename _Type >
DenseVector<double> CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_compY
protected

◆ m_left

template<typename _Type >
double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_left
protected

◆ m_nv

template<typename _Type >
std::size_t CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nv
protected

◆ m_nx

template<typename _Type >
std::size_t CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_nx
protected

◆ m_ny

template<typename _Type >
std::size_t CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_ny
protected

◆ m_right

template<typename _Type >
double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_right
protected

◆ m_top

template<typename _Type >
double CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_top
protected

◆ m_vars

template<typename _Type >
DenseVector<_Type> CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_vars
protected

◆ m_X

template<typename _Type >
DenseVector<double> CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_X
protected

◆ m_Y

template<typename _Type >
DenseVector<double> CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::m_Y
protected

The documentation for this class was generated from the following files:

© 2012

R.E. Hewitt