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

A one dimensional mesh utility object. More...

#include <OneD_Node_Mesh.h>

Public Member Functions

 OneD_Node_Mesh ()
 Default constructor. More...
 
 OneD_Node_Mesh (const DenseVector< _Xtype > &nodes, const std::size_t nvars)
 ctor for a given nodal distribution More...
 
template<typename _sourceType >
 OneD_Node_Mesh (const OneD_Node_Mesh< _sourceType > &source)
 implicit conversion ctor for D_complex from double data More...
 
 OneD_Node_Mesh (std::string filename, const std::size_t nnodes, const std::size_t nvars)
 ctor from an existing file More...
 
virtual ~OneD_Node_Mesh ()
 Destructor. More...
 
_Type & operator() (const std::size_t i, const std::size_t var)
 Access a variable at a node. More...
 
const _Type & operator() (const std::size_t i, const std::size_t var) const
 Access a variable at a node. More...
 
const _Xtype & coord (const std::size_t &node) const
 Access a nodal position. More...
 
_Xtype & coord (const std::size_t &node)
 Access a nodal position. More...
 
void set_nodes_vars (const std::size_t node, const DenseVector< _Type > &u)
 Set the variables stored at A SPECIFIED node. More...
 
DenseVector< _Type > get_nodes_vars (const std::size_t &node) const
 Get the variables stored at A SPECIFIED node. More...
 
DenseVector< _Type > get_interpolated_vars (const _Xtype &pos) const
 Get the variable data at an interpolated position using a first order scheme. More...
 
std::size_t get_nnodes () const
 
std::size_t get_nvars () const
 
const DenseVector< _Xtype > & nodes () const
 
DenseVector< double > find_roots1 (const std::size_t &var, double value=0.0) const
 Find a list of approximate locations at which a specified variable attains a given value. More...
 
_Type integral2 (std::size_t var=0) const
 Integrate over the domain. More...
 
_Xtype squared_integral2 (std::size_t var=0) const
 Compute the integral of the absolute variable squared: |variable|^2. More...
 
_Type integral4 (std::size_t var=0) const
 Integrate over the domain with a Simpson rule. More...
 
const DenseVector< _Type > & vars_as_vector () const
 For each nodal point we push each variable into a vector in sequence. More...
 
void set_vars_from_vector (const DenseVector< _Type > &vec)
 Set the variables of this mesh from a vector. More...
 
const std::vector< DenseVector< _Type > > & get_vars () const
 
void dump () const
 A simple method for dumping data to std::cout. More...
 
void dump_gnu (std::string filename, int precision=10) const
 A simple method for dumping data to a file for gnuplot. More...
 
void remesh1 (const DenseVector< _Xtype > &z)
 Interpolate this mesh data (linearly) into a new mesh with nodal points defined in the argument list. More...
 
void scale (_Type x)
 Scale the whole contents of the mesh. More...
 
void normalise (const std::size_t &var)
 Normalise all data in the mesh based on one variable. More...
 
double max_abs (unsigned var)
 Find the maximum stored absolute value in the mesh for a given variable in a range of the domain. More...
 
void read (std::string filename, const bool reset=false)
 Assign mesh contents using a filename. More...
 
DenseVector< double > find_roots1 (const std::size_t &var, double value) const
 
void remesh1 (const DenseVector< double > &newX)
 
void remesh1 (const DenseVector< double > &newX)
 
void remesh1 (const DenseVector< std::complex< double > > &z)
 
DenseVector< double > get_interpolated_vars (const double &x_pos) const
 
DenseVector< std::complex< double > > get_interpolated_vars (const double &x_pos) const
 
DenseVector< std::complex< double > > get_interpolated_vars (const std::complex< double > &pos) const
 
DenseVector< double > find_roots1 (const std::size_t &var, double value) const
 
void read (std::string filename, bool reset)
 
void read (std::string filename, bool reset)
 
void read (std::string filename, bool reset)
 
void dump_gnu (std::string filename, int precision) const
 
void dump_gnu (std::string filename, int precision) const
 
void dump_gnu (std::string filename, int precision) const
 

Protected Attributes

std::size_t m_nv
 
DenseVector< _Xtype > m_X
 
DenseVector< _Type > m_vars
 

Detailed Description

template<typename _Type, typename _Xtype = double>
class CppNoddy::OneD_Node_Mesh< _Type, _Xtype >

A one dimensional mesh utility object.

Data can be placed into the mesh for interpolation to a new mesh or computation of integrals of the data. The default typing assumes that the mesh is along the real line. You can store (complex) data across a set of points in the complex plane using the second typename.

Definition at line 24 of file OneD_Node_Mesh.h.

Constructor & Destructor Documentation

◆ OneD_Node_Mesh() [1/4]

template<typename _Type , typename _Xtype = double>
CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::OneD_Node_Mesh ( )
inline

Default constructor.

Definition at line 28 of file OneD_Node_Mesh.h.

28 {
29 // default to zero variables
30 m_nv = 0;
31 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_nv.

◆ OneD_Node_Mesh() [2/4]

template<typename _Type , typename _Xtype = double>
CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::OneD_Node_Mesh ( const DenseVector< _Xtype > &  nodes,
const std::size_t  nvars 
)
inline

ctor for a given nodal distribution

Parameters
nodesThe positions of the nodal points
nvarsThe number of variables to store in the mesh

Definition at line 36 of file OneD_Node_Mesh.h.

36 :
37 m_nv(nvars), m_X(nodes) {
38#ifdef PARANOID
39 for(unsigned i = 1; i < m_X.size(); ++i) {
40 /*
41 if ( m_X[i+1] < m_X[i] )
42 {
43 std::string problem;
44 problem = "The OneD_Node_Mesh has been passed a vector of nodes that are\n";
45 problem += "not in INCREASING order. This will screw up some of the methods\n";
46 problem += "in the class. We should fix this .....\n";
47 throw ExceptionRuntime( problem );
48 }
49 */
50 }
51#endif
52 // set the contents to zero
53 m_vars = DenseVector<_Type>(m_nv * m_X.size(), _Type(0.0));
54 }
const DenseVector< _Xtype > & nodes() const
DenseVector< _Type > m_vars
DenseVector< _Xtype > m_X

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_nv, CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_vars, and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_X.

◆ OneD_Node_Mesh() [3/4]

template<typename _Type , typename _Xtype >
template<typename _sourceType >
CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::OneD_Node_Mesh ( const OneD_Node_Mesh< _sourceType > &  source)
inline

implicit conversion ctor for D_complex from double data

Definition at line 259 of file OneD_Node_Mesh.h.

259 {
260 // there is implicit conversion of DenseVector so this will
261 // allow copy of OneD_Node_Mesh<double> to OneD_Node_Mesh<D_complex>.
262 m_nv = source.get_nvars();
263 m_X = source.nodes();
264 m_vars = source.vars_as_vector();
265 }
double source(const double &x, const double &y, const double &t)
Definition: IBVPLinear.cpp:26

◆ OneD_Node_Mesh() [4/4]

template<typename _Type , typename _Xtype = double>
CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::OneD_Node_Mesh ( std::string  filename,
const std::size_t  nnodes,
const std::size_t  nvars 
)
inline

ctor from an existing file

Parameters
filenameFilename of the data file
nnodesNumber of nodes
nvarsNumber of variables stored at each node

Definition at line 64 of file OneD_Node_Mesh.h.

64 :
65 m_nv(nvars) {
66 m_X = DenseVector<_Xtype>(nnodes, _Xtype(0.0) ); //coordinates, currently empty
67 m_vars = DenseVector<_Type>(nvars*nnodes, _Type(0.0) ); // nodal values
68 read(filename, true); // true => reset the nodal coordinates using file
69 }
void read(std::string filename, const bool reset=false)
Assign mesh contents using a filename.

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_vars, CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_X, and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::read().

◆ ~OneD_Node_Mesh()

template<typename _Type , typename _Xtype = double>
virtual CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::~OneD_Node_Mesh ( )
inlinevirtual

Destructor.

Definition at line 72 of file OneD_Node_Mesh.h.

73 {}

Member Function Documentation

◆ coord() [1/2]

template<typename _Type , typename _Xtype >
_Xtype & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord ( const std::size_t &  node)
inline

Access a nodal position.

Parameters
nodeThe nodal position to return
Returns
The spatial position of this node

Definition at line 283 of file OneD_Node_Mesh.h.

283 {
284 return m_X[ node ];
285 }

◆ coord() [2/2]

template<typename _Type , typename _Xtype >
const _Xtype & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord ( const std::size_t &  node) const
inline

Access a nodal position.

Parameters
nodeThe nodal position to return
Returns
The spatial position of this node

Definition at line 278 of file OneD_Node_Mesh.h.

278 {
279 return m_X[ node ];
280 }

Referenced by CppNoddy::FT::dft(), CppNoddy::FT::dft_with_shift(), CppNoddy::HST::Orr_Sommerfeld::global_evp(), CppNoddy::FT::idft(), CppNoddy::FT::idft_with_ishift(), CppNoddy::FT::ishift(), main(), CppNoddy::Utility::max_abs_location(), CppNoddy::Utility::max_abs_location_range(), and CppNoddy::FT::shift().

◆ dump()

template<typename _Type , typename _Xtype >
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::dump

A simple method for dumping data to std::cout.

Definition at line 499 of file OneD_Node_Mesh.cpp.

499 {
500 std::cout << "Number of nodes = " << m_X.size() << "\n";
501 std::cout << "Nodal positions :\n";
502 m_X.dump();
503 std::cout << "\n";
504 std::cout << "Number of vars = " << m_nv << "\n";
505 std::cout << "Interleaved mesh data : \n";
506 m_vars.dump();
507 std::cout << "Mesh dump complete\n";
508 }

Referenced by main().

◆ dump_gnu() [1/4]

void CppNoddy::OneD_Node_Mesh< double, double >::dump_gnu ( std::string  filename,
int  precision 
) const

Definition at line 511 of file OneD_Node_Mesh.cpp.

511 {
512 std::ofstream dump;
513 dump.open(filename.c_str());
514 dump.precision(precision);
515 dump.setf(std::ios::showpoint);
516 dump.setf(std::ios::showpos);
517 dump.setf(std::ios::scientific);
518 for(std::size_t i = 0; i < m_X.size(); ++i) {
519 dump << m_X[ i ] << " ";
520 for(std::size_t var = 0; var < m_nv; ++var) {
521 dump << m_vars[ i * m_nv + var ] << " ";
522 }
523 dump << "\n";
524 }
525 }
void dump() const
A simple method for dumping data to std::cout.

◆ dump_gnu() [2/4]

void CppNoddy::OneD_Node_Mesh< std::complex< double >, double >::dump_gnu ( std::string  filename,
int  precision 
) const

Definition at line 528 of file OneD_Node_Mesh.cpp.

528 {
529 std::ofstream dump;
530 dump.open(filename.c_str());
531 dump.precision(precision);
532 dump.setf(std::ios::showpoint);
533 dump.setf(std::ios::showpos);
534 dump.setf(std::ios::scientific);
535 for(std::size_t i = 0; i < m_X.size(); ++i) {
536 dump << m_X[ i ] << " ";
537 for(std::size_t var = 0; var < m_nv; ++var) {
538 dump << real(m_vars[ i * m_nv + var ]) << " " << imag(m_vars[ i * m_nv + var ]) << " ";
539 }
540 dump << "\n";
541 }
542 }
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

◆ dump_gnu() [3/4]

void CppNoddy::OneD_Node_Mesh< std::complex< double >, std::complex< double > >::dump_gnu ( std::string  filename,
int  precision 
) const

Definition at line 545 of file OneD_Node_Mesh.cpp.

545 {
546 std::ofstream dump;
547 dump.open(filename.c_str());
548 dump.precision(precision);
549 dump.setf(std::ios::showpoint);
550 dump.setf(std::ios::showpos);
551 dump.setf(std::ios::scientific);
552 for(std::size_t i = 0; i < m_X.size(); ++i) {
553 dump << real(m_X[ i ]) << " " << imag(m_X[ i ]) << " ";
554 for(std::size_t var = 0; var < m_nv; ++var) {
555 dump << real(m_vars[ i * m_nv + var ]) << " " << imag(m_vars[ i * m_nv + var ]) << " ";
556 }
557 dump << "\n";
558 }
559 }

◆ dump_gnu() [4/4]

template<typename _Type , typename _Xtype = double>
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::dump_gnu ( std::string  filename,
int  precision = 10 
) const

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

Parameters
filenameThe filename to write the data to (will overwrite)
precisionPrecision of the output strings

Referenced by main().

◆ find_roots1() [1/3]

DenseVector< double > CppNoddy::OneD_Node_Mesh< double, double >::find_roots1 ( const std::size_t &  var,
double  value 
) const

Definition at line 65 of file OneD_Node_Mesh.cpp.

65 {
66 DenseVector<double> roots;
67 for(std::size_t node = 0; node < m_X.size() - 1; ++node) {
68 std::size_t offset(node * m_nv + var);
69 // find bracket nodes
70 if((m_vars[ offset ] - value) * (m_vars[ offset + m_nv ] - value) < 0.0) {
71 double deriv = (m_vars[ offset + m_nv ] - m_vars[ offset ]) / (m_X[ node + 1 ] - m_X[ node ]);
72 double x = m_X[ node ] + (value - m_vars[ offset ]) / deriv;
73 // add the left hand node to the roots vector
74 roots.push_back(x);
75 }
76 }
77 return roots;
78 }

References CppNoddy::DenseVector< _Type >::push_back().

◆ find_roots1() [2/3]

DenseVector< double > CppNoddy::OneD_Node_Mesh< std::complex< double >, double >::find_roots1 ( const std::size_t &  var,
double  value 
) const

Definition at line 306 of file OneD_Node_Mesh.cpp.

306 {
307 std::string problem;
308 problem = " The OneD_Node_Mesh.find_roots1 method has been called with \n";
309 problem += " a mesh containing complex data.\n";
310 throw ExceptionRuntime(problem);
311 }

◆ find_roots1() [3/3]

template<typename _Type , typename _Xtype = double>
DenseVector< double > CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::find_roots1 ( const std::size_t &  var,
double  value = 0.0 
) const

Find a list of approximate locations at which a specified variable attains a given value.

First order only.

Parameters
varThe variable to be examined for zeros
valueThe value to find

Referenced by main().

◆ get_interpolated_vars() [1/4]

template<typename _Type , typename _Xtype = double>
DenseVector< _Type > CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_interpolated_vars ( const _Xtype &  pos) const

Get the variable data at an interpolated position using a first order scheme.

Doesn't really make sense unless the data is along the real line.

Parameters
posThe position to interpolate at.
Returns
A vector of interpolated variables.

Referenced by CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_interpolated_vars(), CppNoddy::TwoD_Node_Mesh< _Type >::get_interpolated_vars(), and main().

◆ get_interpolated_vars() [2/4]

DenseVector< double > CppNoddy::OneD_Node_Mesh< double, double >::get_interpolated_vars ( const double &  x_pos) const

Definition at line 209 of file OneD_Node_Mesh.cpp.

209 {
210 for(unsigned node = 0; node < m_X.size() - 1; ++node) {
211 // find bracketing nodes - incl shameless hack for evaluations at the boundary
212 //if ( ( m_X[ node ] < x_pos || std::abs( m_X[ node ] - x_pos ) < 1.e-7 ) &&
213 //( m_X[ node + 1 ] > x_pos || std::abs( m_X[ node + 1 ] - x_pos ) < 1.e-7 ) )
214 if(((m_X[ node ] < x_pos) && (m_X[ node + 1 ] > x_pos))
215 || (std::abs(m_X[ node ] - x_pos) < 1.e-7) || (std::abs(m_X[ node + 1 ] - x_pos) < 1.e-7)) {
216 // distance from left node
217 double delta_x(x_pos - m_X[ node ]);
218 // empty data to return
219 DenseVector<double> left;
220 DenseVector<double> right;
221 DenseVector<double> deriv;
222 // interpolate data linearly
223 left = get_nodes_vars(node);
224 right = get_nodes_vars(node + 1);
225 deriv = (right - left) / (m_X[ node + 1 ] - m_X[ node ]);
226 // overwrite right
227 right = left + deriv * delta_x;
228 return right;
229 }
230 }
231 std::cout << "You asked for a position of " << x_pos << " in a range " << m_X[ 0 ] << " to " << m_X[ m_X.size() - 1 ] << "\n";
232 std::string problem;
233 problem = "You have asked the OneD_Node_Mesh class to interpolate data at\n";
234 problem += "a point that is outside the range covered by the mesh object.\n";
235 throw ExceptionRuntime(problem);
236 }
DenseVector< _Type > get_nodes_vars(const std::size_t &node) const
Get the variables stored at A SPECIFIED node.

◆ get_interpolated_vars() [3/4]

DenseVector< std::complex< double > > CppNoddy::OneD_Node_Mesh< std::complex< double >, double >::get_interpolated_vars ( const double &  x_pos) const

Definition at line 240 of file OneD_Node_Mesh.cpp.

240 {
241 for(unsigned node = 0; node < m_X.size() - 1; ++node) {
242 // find bracketing nodes - incl shameless hack for evaluations at the boundary
243 if((m_X[ node ] < x_pos || std::abs(m_X[ node ] - x_pos) < 1.e-7) &&
244 (m_X[ node + 1 ] > x_pos || std::abs(m_X[ node + 1 ] - x_pos) < 1.e-7)) {
245 // distance from left node
246 double delta_x(x_pos - m_X[ node ]);
247 // empty data to return
248 DenseVector<std::complex<double> > left;
249 DenseVector<std::complex<double> > right;
250 DenseVector<std::complex<double> > deriv;
251 // interpolate data linearly
252 left = get_nodes_vars(node);
253 right = get_nodes_vars(node + 1);
254 deriv = (right - left) / (m_X[ node + 1 ] - m_X[ node ]);
255 // overwrite right
256 right = left + deriv * delta_x;
257 return right;
258 }
259 }
260 std::cout << "You asked for a position of " << x_pos << " in a range " << m_X[ 0 ] << " to " << m_X[ m_X.size() - 1 ] << "\n";
261 std::string problem;
262 problem = "You have asked the OneD_Node_Mesh class to interpolate data at\n";
263 problem += "a point that is outside the range covered by the mesh object.\n";
264 throw ExceptionRuntime(problem);
265 }

◆ get_interpolated_vars() [4/4]

DenseVector< std::complex< double > > CppNoddy::OneD_Node_Mesh< std::complex< double >, std::complex< double > >::get_interpolated_vars ( const std::complex< double > &  pos) const

Definition at line 269 of file OneD_Node_Mesh.cpp.

269 {
270 double x_pos(pos.real());
271#ifdef PARANOID
272 std::cout << "WARNING: You are interpolating complex data on a complex mesh with 'get_interpolated_vars'.\n";
273 std::cout << " This does a simple piecewise linear interpolating assuming a single valued path. \n";
274#endif
275 for(unsigned node = 0; node < m_X.size() - 1; ++node) {
276 // find bracketing nodes - incl shameless hack for evaluations at the boundary
277 if((m_X[ node ].real() < x_pos || std::abs(m_X[ node ].real() - x_pos) < 1.e-7) &&
278 (m_X[ node + 1 ].real() > x_pos || std::abs(m_X[ node + 1 ].real() - x_pos) < 1.e-7)) {
279 // distance from left node -- real coordinate is given. We also need to
280 // interpolate between the two complex nodes -- hence imaginary coordinate is implict from
281 // bracketing (complex) nodes
282 std::complex<double> delta_z = (m_X[ node + 1 ] - m_X[ node ]) * (x_pos - m_X[ node ].real()) / (m_X[ node + 1 ].real() - m_X[ node ].real());
283 // empty data to return
284 DenseVector<std::complex<double> > left;
285 DenseVector<std::complex<double> > right;
286 DenseVector<std::complex<double> > deriv;
287 // interpolate data linearly
288 left = get_nodes_vars(node);
289 right = get_nodes_vars(node + 1);
290 // derivative of the data
291 deriv = (right - left) / (m_X[ node + 1 ] - m_X[ node ]);
292 // overwrite right
293 right = left + deriv * delta_z;
294 return right;
295 }
296 }
297 std::cout << "You asked for a position of " << x_pos << " in a range " << m_X[ 0 ] << " to " << m_X[ m_X.size() - 1 ] << "\n";
298 std::string problem;
299 problem = "You have asked the OneD_Node_Mesh class to interpolate data at\n";
300 problem += "a point that is outside the range covered by the mesh object.\n";
301 problem += "Even for complex nodes we assume the path is single valued.\n";
302 throw ExceptionRuntime(problem);
303 }

◆ get_nnodes()

template<typename _Type , typename _Xtype >
std::size_t CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nnodes

◆ get_nodes_vars()

template<typename _Type , typename _Xtype >
DenseVector< _Type > CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nodes_vars ( const std::size_t &  node) const

Get the variables stored at A SPECIFIED node.

Parameters
nodeThe nodal index to be returned
Returns
The vector of VARIABLES stored at this nodal point

Definition at line 33 of file OneD_Node_Mesh.cpp.

33 {
34#ifdef PARANOID
35 if((node >= m_X.size()) || (node < 0)) {
36 std::string problem;
37 problem = " The OneD_Node_Mesh.get_nodes_vars method is trying to access \n";
38 problem += " a nodal point outside of the range stored. \n";
39 throw ExceptionRange(problem, m_X.size(), node);
40 }
41#endif
42 DenseVector<_Type> nodes_vars;
43 for(std::size_t var = 0; var < m_nv; ++var) {
44 nodes_vars.push_back(m_vars[ node * m_nv + var ]);
45 }
46 return nodes_vars;
47 }

References CppNoddy::DenseVector< _Type >::push_back().

Referenced by CppNoddy::FT::dft_with_shift(), CppNoddy::FT::idft(), and CppNoddy::FT::idft_with_ishift().

◆ get_nvars()

template<typename _Type , typename _Xtype >
std::size_t CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nvars
Returns
The number of variables that have data stored at each nodal point

Definition at line 55 of file OneD_Node_Mesh.cpp.

55 {
56 return m_nv;
57 }

Referenced by CppNoddy::FT::idft(), CppNoddy::FT::ishift(), and CppNoddy::FT::shift().

◆ get_vars()

template<typename _Type , typename _Xtype = double>
const std::vector< DenseVector< _Type > > & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_vars ( ) const
Returns
An STL vector of dense vectors of the variables in the mesh

◆ integral2()

template<typename _Type , typename _Xtype >
_Type CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::integral2 ( std::size_t  var = 0) const

Integrate over the domain.

Typically useful for finite volume methods.

Parameters
varThe variable-index to be integrated over the mesh using a trapezium rule.
Returns
The integral value.

Definition at line 314 of file OneD_Node_Mesh.cpp.

314 {
315 _Type sum = 0.0;
316 _Xtype dx = 0.0;
317 // sum interior segments
318 for(std::size_t node = 0; node < m_X.size() - 1; ++node) {
319 dx = (m_X[ node + 1 ] - m_X[ node ]);
320 sum += 0.5 * dx * (m_vars[ node * m_nv + var ] + m_vars[(node+1) * m_nv + var ]);
321 }
322 // return the value
323 return sum;
324 }

Referenced by main().

◆ integral4()

template<typename _Type , typename _Xtype >
_Type CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::integral4 ( std::size_t  var = 0) const

Integrate over the domain with a Simpson rule.

Parameters
varThe variable-index to be integrated over the mesh using a trapezium rule.
Returns
The integral value.

Definition at line 341 of file OneD_Node_Mesh.cpp.

341 {
342 if((m_X.size()) % 2 == 0) {
343 std::string problem;
344 problem = " The OneD_Node_Mesh.Simpson_integral method is trying to run \n";
345 problem += " on a mesh with an even number of points. \n";
346 throw ExceptionRuntime(problem);
347 }
348
349 _Type f0, f1, f2;
350 _Xtype x0, x1, x2;
351
352 _Type sum = 0.0;
353
354 // sum interior segments
355 for(std::size_t node = 0; node < m_X.size() - 2; node += 2) {
356 x0 = m_X[ node ];
357 x1 = m_X[ node + 1 ];
358 x2 = m_X[ node + 2 ];
359 f0 = m_vars[ node * m_nv + var ];
360 f1 = m_vars[(node+1) * m_nv + var ];
361 f2 = m_vars[(node+2) * m_nv + var ];
362 sum += (x2 - x0)
363 * (
364 f1 * pow(x0 - x2, 2) + f0 * (x1 - x2) * (2. * x0 - 3. * x1 + x2)
365 - f2 * (x0 - x1) * (x0 - 3. * x1 + 2. * x2)
366 )
367 / (6. * (x0 - x1) * (x1 - x2));
368 // sum += (x1-x0)*( f0 + 4*f1 + f2 ) / 3.0 for equal spacing
369 }
370 // return the value
371 return sum;
372 }

Referenced by main().

◆ max_abs()

template<typename _Type , typename _Xtype = double>
double CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::max_abs ( unsigned  var)
inline

Find the maximum stored absolute value in the mesh for a given variable in a range of the domain.

Parameters
varThe variable index whose maximum is being asked for
leftOnly examine the sub-range x>left
rightOnly examine the sub-range x<right
Returns
The value of the maximum (abs value) 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 228 of file OneD_Node_Mesh.h.

228 {
229 double max(0.0);
230 // step through the nodes
231 for(unsigned node = 0; node < m_X.size(); ++node) {
232 if(std::abs(m_vars[ node * m_nv + var ]) > max) {
233 max = std::abs(m_vars[ node * m_nv + var ]);
234 }
235 }
236 return max;
237 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_nv, CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_vars, and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_X.

Referenced by CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::normalise().

◆ nodes()

template<typename _Type , typename _Xtype >
const DenseVector< _Xtype > & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::nodes
Returns
A vector of the nodal positions for this mesh

Definition at line 60 of file OneD_Node_Mesh.cpp.

60 {
61 return m_X;
62 }

Referenced by CppNoddy::HST::Rayleigh< _Type >::eigenvector().

◆ normalise()

template<typename _Type , typename _Xtype = double>
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::normalise ( const std::size_t &  var)
inline

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.

Definition at line 184 of file OneD_Node_Mesh.h.

184 {
185 double maxval(max_abs(var));
186 m_vars.scale(1./maxval);
187 }
double max_abs(unsigned var)
Find the maximum stored absolute value in the mesh for a given variable in a range of the domain.

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_vars, and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::max_abs().

Referenced by main().

◆ operator()() [1/2]

template<typename _Type , typename _Xtype >
_Type & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::operator() ( const std::size_t  i,
const std::size_t  var 
)
inline

Access a variable at a node.

Parameters
iThe index of the node to be accessed
varThe variable to return the data for
Returns
The variable stored at the node

Definition at line 268 of file OneD_Node_Mesh.h.

268 {
269 return m_vars[ i * m_nv + var ];
270 }

◆ operator()() [2/2]

template<typename _Type , typename _Xtype >
const _Type & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::operator() ( const std::size_t  i,
const std::size_t  var 
) const
inline

Access a variable at a node.

Parameters
iThe index of the node to be accessed
varThe variable to return the data for
Returns
The variable stored at the node

Definition at line 273 of file OneD_Node_Mesh.h.

273 {
274 return m_vars[ i * m_nv + var ];
275 }

◆ read() [1/4]

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

Definition at line 377 of file OneD_Node_Mesh.cpp.

377 {
378 std::ifstream dump;
379 std::cout << "[INFO] Reading OneD_Node_Mesh<double,double> data from " + filename + "\n";
380 dump.open(filename.c_str());
381 if (dump.good() != true) {
382 std::string problem;
383 problem = " The OneD_Node_Mesh.read method is trying to read a \n";
384 problem += " file (" + filename + ") that doesn't exist.\n";
385 throw ExceptionRuntime(problem);
386 }
387 dump.precision(15);
388 dump.setf(std::ios::showpoint);
389 dump.setf(std::ios::showpos);
390 dump.setf(std::ios::scientific);
391 std::size_t N = get_nnodes();
392 for (std::size_t i = 0; i < N; ++i) {
393 double x;
394 dump >> x;
395 for (std::size_t var = 0; var < m_nv; ++var) {
396 double value;
397 dump >> value;
398 m_vars[ i * m_nv + var ] = value;
399 }
400 if (reset != true) {
401 // if not reseting the mesh we should check the node positions
402 if (std::fabs(x - m_X[ i ]) > 1.e-6) {
403 std::cout << " Read x = " << x << " Expected x = " << m_X[ i ] << "\n";
404 std::cout << " Absolute differences are " << fabs(x - m_X[i]) << "\n";
405 std::string problem;
406 problem = " The TwoD_Node_Mesh.read method is trying to read a \n";
407 problem += " file whose nodal points are in a different position. \n";
408 throw ExceptionRuntime(problem);
409 }
410 } else {
411 m_X[ i ] = x;
412 }
413 }
414 }
std::size_t get_nnodes() const

◆ read() [2/4]

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

Definition at line 417 of file OneD_Node_Mesh.cpp.

417 {
418 std::ifstream dump;
419 std::cout << "[INFO] Reading OneD_Node_Mesh<D_complex,double> data from " + filename + "\n";
420 dump.open(filename.c_str());
421 if (dump.good() != true) {
422 std::string problem;
423 problem = " The OneD_Node_Mesh.read method is trying to read a \n";
424 problem += " file (" + filename + ") that doesn't exist.\n";
425 throw ExceptionRuntime(problem);
426 }
427 dump.precision(15);
428 dump.setf(std::ios::showpoint);
429 dump.setf(std::ios::showpos);
430 dump.setf(std::ios::scientific);
431 std::size_t N = get_nnodes();
432 for (std::size_t i = 0; i < N; ++i) {
433 double x;
434 dump >> x;
435 for (std::size_t var = 0; var < m_nv; ++var) {
436 double rValue, iValue;
437 dump >> rValue; dump >> iValue;
438 m_vars[ i * m_nv + var ] = D_complex(rValue,iValue);
439 }
440 if (reset != true) {
441 // if not reseting the mesh we should check the node positions
442 if (std::fabs(x - m_X[ i ]) > 1.e-6) {
443 std::cout << " Read x = " << x << " Expected x = " << m_X[ i ] << "\n";
444 std::cout << " Absolute differences are " << fabs(x - m_X[i]) << "\n";
445 std::string problem;
446 problem = " The TwoD_Node_Mesh.read method is trying to read a \n";
447 problem += " file whose nodal points are in a different position. \n";
448 throw ExceptionRuntime(problem);
449 }
450 } else {
451 m_X[ i ] = x;
452 }
453 }
454 }
std::complex< double > D_complex
A complex double precision number using std::complex.
Definition: Types.h:98

◆ read() [3/4]

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

Definition at line 457 of file OneD_Node_Mesh.cpp.

457 {
458 std::ifstream dump;
459 std::cout << "[INFO] Reading OneD_Node_Mesh<D_complex, D_complex> data from " + filename + "\n";
460 dump.open(filename.c_str());
461 if (dump.good() != true) {
462 std::string problem;
463 problem = " The OneD_Node_Mesh.read method is trying to read a \n";
464 problem += " file (" + filename + ") that doesn't exist.\n";
465 throw ExceptionRuntime(problem);
466 }
467 dump.precision(15);
468 dump.setf(std::ios::showpoint);
469 dump.setf(std::ios::showpos);
470 dump.setf(std::ios::scientific);
471 std::size_t N = get_nnodes();
472 for (std::size_t i = 0; i < N; ++i) {
473 double rX, iX;
474 dump >> rX; dump >> iX;
475 D_complex X(rX,iX);
476 for (std::size_t var = 0; var < m_nv; ++var) {
477 double rValue, iValue;
478 dump >> rValue; dump >> iValue;
479 m_vars[ i * m_nv + var ] = D_complex(rValue,iValue);
480 }
481 if (reset != true) {
482 // if not reseting the mesh we should check the node positions
483 if (std::fabs( X - m_X[ i ]) > 1.e-6) {
484 std::cout << " Read x = " << X << " Expected x = " << m_X[ i ] << "\n";
485 std::cout << " Absolute differences are " << fabs( X - m_X[i]) << "\n";
486 std::string problem;
487 problem = " The TwoD_Node_Mesh.read method is trying to read a \n";
488 problem += " file whose nodal points are in a different position. \n";
489 throw ExceptionRuntime(problem);
490 }
491 } else {
492 m_X[ i ] = X;
493 }
494 }
495 }

◆ read() [4/4]

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

Assign mesh contents using a filename.

Parameters
filenameFilename that contains the data
resetA boolean, if true then coordinate data is overwritten using the file data. Default is false.

Referenced by CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::OneD_Node_Mesh().

◆ remesh1() [1/4]

template<typename _Type , typename _Xtype = double>
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::remesh1 ( const DenseVector< _Xtype > &  z)

Interpolate this mesh data (linearly) into a new mesh with nodal points defined in the argument list.

Parameters
zThe nodal coordinates to be used in the new mesh.

Referenced by main(), and CppNoddy::HST::Rayleigh< _Type >::remesh1().

◆ remesh1() [2/4]

void CppNoddy::OneD_Node_Mesh< double, double >::remesh1 ( const DenseVector< double > &  newX)

Definition at line 99 of file OneD_Node_Mesh.cpp.

99 {
100#ifdef PARANOID
101 if(std::abs(m_X[ 0 ] - newX[ 0 ]) > 1.e-10 ||
102 std::abs(m_X[ m_X.size() - 1 ] - newX[ newX.size() - 1 ]) > 1.e-10) {
103 std::string problem;
104 problem = " The OneD_Node_Mesh.remesh method has been called with \n";
105 problem += " a passed coordinate vector that has different start and/or \n";
106 problem += " end points from the instantiated object. \n";
107 throw ExceptionRuntime(problem);
108 }
109
110 for(std::size_t i = 0; i < newX.size() - 1; ++i) {
111 if(newX[ i ] >= newX[ i + 1 ]) {
112 std::string problem;
113 problem = " The OneD_Node_Mesh.remesh method has been passed \n";
114 problem += " a non-monotonic coordinate vector. \n";
115 throw ExceptionRuntime(problem);
116 }
117 }
118#endif
119 // copy current state of this mesh
120 DenseVector<double> copy_of_vars(m_vars);
121 // resize the local storage
122 m_vars.resize(newX.size() * m_nv);
123
124 // first nodal values are assumed to be untouched
125 // loop thru destination mesh node at a time
126 for(std::size_t node = 1; node < newX.size() - 1; ++node) {
127 // loop through the source mesh and find the bracket-nodes
128 for(std::size_t i = 0; i < m_X.size(); ++i) {
129 if((m_X[ i ] <= newX[ node ]) && (newX[ node ] < m_X[ i + 1 ])) {
130 // linearly interpolate each variable in the mesh
131 for(std::size_t var = 0; var < m_nv; ++var) {
132 double dX = newX[ node ] - m_X[ i ];
133 double dvarsdX = (copy_of_vars[(i+1)*m_nv + var ] - copy_of_vars[ i*m_nv + var ]) / (m_X[ i + 1 ] - m_X[ i ]);
134 m_vars[ node * m_nv + var ] = copy_of_vars[ i * m_nv + var ] + dX * dvarsdX;
135 }
136 }
137 }
138 }
139
140 // add the last nodal values to the resized vector
141 for(std::size_t var = 0; var < m_nv; ++var) {
142 m_vars[(newX.size() - 1) * m_nv + var ] = copy_of_vars[(m_X.size() - 1) * m_nv + var ];
143 }
144 // replace the old nodes with the new ones
145 m_X = newX;
146 }
std::size_t size() const
A pass-thru definition to get the size of the vector.
Definition: DenseVector.h:330

◆ remesh1() [3/4]

void CppNoddy::OneD_Node_Mesh< std::complex< double >, double >::remesh1 ( const DenseVector< double > &  newX)

Definition at line 149 of file OneD_Node_Mesh.cpp.

149 {
150#ifdef PARANOID
151 if(std::abs(m_X[ 0 ] - newX[ 0 ]) > 1.e-10 ||
152 std::abs(m_X[ m_X.size() - 1 ] - newX[ newX.size() - 1 ]) > 1.e-10) {
153 std::string problem;
154 problem = " The OneD_Node_Mesh.remesh method has been called with \n";
155 problem += " a passed coordinate vector that has different start and/or \n";
156 problem += " end points from the instantiated object. \n";
157 throw ExceptionRuntime(problem);
158 }
159
160 for(std::size_t i = 0; i < newX.size() - 1; ++i) {
161 if(newX[ i ] >= newX[ i + 1 ]) {
162 std::string problem;
163 problem = " The OneD_Node_Mesh.remesh method has been passed \n";
164 problem += " a non-monotonic coordinate vector. \n";
165 throw ExceptionRuntime(problem);
166 }
167 }
168#endif
169 // copy current state of this mesh
170 DenseVector<std::complex<double> > copy_of_vars(m_vars);
171 // resize the local storage
172 m_vars.resize(newX.size() * m_nv);
173
174 // first nodal values are assumed to be untouched
175 // loop thru destination mesh node at a time
176 for(std::size_t node = 1; node < newX.size() - 1; ++node) {
177 // loop through the source mesh and find the bracket-nodes
178 for(std::size_t i = 0; i < m_X.size(); ++i) {
179 if((m_X[ i ] <= newX[ node ]) && (newX[ node ] < m_X[ i + 1 ])) {
180 // linearly interpolate each variable in the mesh
181 for(std::size_t var = 0; var < m_nv; ++var) {
182 double dX = newX[ node ] - m_X[ i ];
183 // if the paranoid checks above are satisfied, then the m_X[ i + 1 ] should still be in bounds
184 std::complex<double> dvarsdX = (copy_of_vars[(i+1)*m_nv + var ] - copy_of_vars[ i*m_nv + var ]) / (m_X[ i + 1 ] - m_X[ i ]);
185 m_vars[ node * m_nv + var ] = copy_of_vars[ i * m_nv + var ] + dX * dvarsdX;
186 }
187 }
188 }
189 }
190
191 // add the last nodal values to the resized vector
192 for(std::size_t var = 0; var < m_nv; ++var) {
193 m_vars[(newX.size() - 1) * m_nv + var ] = copy_of_vars[(m_X.size() - 1) * m_nv + var ];
194 }
195 // replace the old nodes with the new ones
196 m_X = newX;
197
198 }

◆ remesh1() [4/4]

void CppNoddy::OneD_Node_Mesh< std::complex< double >, std::complex< double > >::remesh1 ( const DenseVector< std::complex< double > > &  z)

Definition at line 201 of file OneD_Node_Mesh.cpp.

201 {
202 std::string problem;
203 problem = " The OneD_Node_Mesh.remesh method has been called with \n";
204 problem += " a complex data set on a complex mesh.\n";
205 throw ExceptionRuntime(problem);
206 }

◆ scale()

template<typename _Type , typename _Xtype = double>
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::scale ( _Type  x)
inline

Scale the whole contents of the mesh.

Parameters
xThe value to multiply the contents of the mesh by

Definition at line 176 of file OneD_Node_Mesh.h.

176 {
177 m_vars.scale(x);
178 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_vars.

◆ set_nodes_vars()

template<typename _Type , typename _Xtype >
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_nodes_vars ( const std::size_t  node,
const DenseVector< _Type > &  u 
)

Set the variables stored at A SPECIFIED node.

Parameters
nodeThe nodal index to be set
uThe vector of VARIABLES to be written to this nodal point

Definition at line 17 of file OneD_Node_Mesh.cpp.

17 {
18#ifdef PARANOID
19 if(U.size() != m_nv) {
20 std::string problem;
21 problem = " The OneD_Node_Mesh.set_node method is trying to add \n";
22 problem += " an DenseVector of variables of a different size to that \n";
23 problem += " stored at other nodal points. \n";
24 throw ExceptionGeom(problem, m_nv, U.size());
25 }
26#endif
27 for(std::size_t var = 0; var < U.size(); ++var) {
28 m_vars[ node * m_nv + var ] = U[ var ];
29 }
30 }
@ U
Definition: BVPKarman.cpp:20

References U.

Referenced by CppNoddy::ODE_EVP< _Type >::add_tagged_to_mesh(), CppNoddy::FT::dft(), CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_interpolated_vars(), CppNoddy::OneD_TVDLF_Mesh::get_slope(), CppNoddy::OneD_TVDLF_Mesh::get_soln(), CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_xsection_at_xnode(), CppNoddy::TwoD_Node_Mesh< _Type >::get_xsection_at_xnode(), CppNoddy::TwoD_Mapped_Node_Mesh< _Type >::get_xsection_at_ynode(), and CppNoddy::TwoD_Node_Mesh< _Type >::get_xsection_at_ynode().

◆ set_vars_from_vector()

template<typename _Type , typename _Xtype >
void CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_vars_from_vector ( const DenseVector< _Type > &  vec)

Set the variables of this mesh from a vector.

Parameters
vecThe vector to be used.

Definition at line 86 of file OneD_Node_Mesh.cpp.

86 {
87#ifdef PARANOID
88 if(vec.size() != m_nv * m_X.size()) {
89 std::string problem;
90 problem = "The set_vars_from_vector method has been passed a vector\n";
91 problem += "of a length that is of an incompatible size for this mesh object\n";
92 throw ExceptionRuntime(problem);
93 }
94#endif
95 m_vars = vec;
96 }

References CppNoddy::DenseVector< _Type >::size().

Referenced by main().

◆ squared_integral2()

template<typename _Type , typename _Xtype >
_Xtype CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::squared_integral2 ( std::size_t  var = 0) const

Compute the integral of the absolute variable squared: |variable|^2.

Parameters
varThe variable-index to be integrated
Returns
The integral of the square of the absolute value.

Definition at line 327 of file OneD_Node_Mesh.cpp.

327 {
328 _Xtype sum = 0.0;
329 double dx = 0.0;
330 // sum interior segments
331 for(std::size_t node = 0; node < m_X.size() - 1; ++node) {
332 dx = std::abs(m_X[ node + 1 ] - m_X[ node ]);
333 sum += 0.5 * dx * (std::pow(std::abs(m_vars[ node * m_nv + var ]), 2)
334 + std::pow(std::abs(m_vars[(node + 1) * m_nv + var ]), 2));
335 }
336 // return the value
337 return std::sqrt(sum);
338 }

◆ vars_as_vector()

template<typename _Type , typename _Xtype >
const DenseVector< _Type > & CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::vars_as_vector

For each nodal point we push each variable into a vector in sequence.

So the returned vector has the data v_00,,...,v_0n,v_10,...,v1n,....v_mn, where the first index denotes the nodal point 0-m and the second the variable 0-n.

Returns
A vector of the variables stored in the mesh

Definition at line 81 of file OneD_Node_Mesh.cpp.

81 {
82 return m_vars;
83 }

Member Data Documentation

◆ m_nv

template<typename _Type , typename _Xtype = double>
std::size_t CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_nv
protected

◆ m_vars

template<typename _Type , typename _Xtype = double>
DenseVector<_Type> CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_vars
protected

◆ m_X

template<typename _Type , typename _Xtype = double>
DenseVector<_Xtype> CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::m_X
protected

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

© 2012

R.E. Hewitt