27 GISMO_ASSERT( 1 == m_nodes.rows(),
"Inconsistent quadrature mapping");
29 const T h = (endVal-startVal) / (T)(2);
32 nodes = (h * (m_nodes.array()+1)) + startVal;
35 weights.noalias() = (0==h?T(0.5):h) * m_weights;
43 GISMO_ASSERT( ab.rows() == m_nodes.rows(),
"Inconsistent quadrature mapping");
44 nodes.resize( m_nodes.rows(), m_nodes.cols() );
46 const gsVector<T> h = (ab.col(1)-ab.col(0)) / (T)(2) ;
47 nodes.noalias() = ( h.asDiagonal() * (m_nodes.array()+1).matrix() ).colwise() + ab.col(0);
55 GISMO_ASSERT( 1 == m_nodes.rows(),
"Inconsistent quadrature mapping.");
56 GISMO_ASSERT( breaks.size()>1,
"At least 2 breaks are needed.");
58 const size_t nint = breaks.size() - 1;
59 const index_t nnodes = numNodes();
61 nodes .resize(1, nint*nnodes);
62 weights.resize( nint*nnodes );
64 for (
size_t i = 0; i!=nint; ++i)
66 const T startVal = breaks[i ];
67 const T endVal = breaks[i+1];
68 const T h = (endVal-startVal) / (T)(2);
71 nodes.middleCols(i*nnodes,nnodes) = (h * (m_nodes.array()+1)) + startVal;
74 weights.segment(i*nnodes,nnodes) = (0==h?T(0.5):h) * m_weights;
83 this->computeTensorProductRule_into(nodes,weights,m_nodes,m_weights);
94 "Nodes and weights do not agree." );
101 numNodes[i] = weights[i].rows();
104 "Inconsistent sizes in nodes and weights.");
107 targetWeights.resize( targetNodes.cols() );
112 targetWeights[r] = weights[0][curr[0]];
114 targetWeights[r] *= weights[i][curr[i]];
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class representing a reference quadrature rule.
Definition gsQuadRule.h:29
void computeTensorProductRule(const std::vector< gsVector< T > > &nodes, const std::vector< gsVector< T > > &weights)
Computes the tensor product rule from coordinate-wise 1D nodes and weights.
Definition gsQuadRule.hpp:80
void mapToAll(const std::vector< T > &breaks, gsMatrix< T > &nodes, gsVector< T > &weights) const
Maps a univariate quadrature rule (i.e., points and weights) from the reference interval to a number ...
Definition gsQuadRule.hpp:52
virtual void mapTo(const gsVector< T > &lower, const gsVector< T > &upper, gsMatrix< T > &nodes, gsVector< T > &weights) const
Maps quadrature rule (i.e., points and weights) from the reference domain to an element.
Definition gsQuadRule.h:177
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
gsMatrix< T > gsPointGrid(gsVector< T > const &a, gsVector< T > const &b, gsVector< unsigned > const &np)
Construct a Cartesian grid of uniform points in a hypercube, using np[i] points in direction i.
Definition gsPointGrid.hpp:82
bool nextLexicographic(Vec &cur, const Vec &size)
Iterates through a tensor lattice with the given size. Updates cur and returns true if another entry ...
Definition gsCombinatorics.h:196
Provides combinatorial unitilies.
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides functions to generate structured point data.
The G+Smo namespace, containing all definitions for the library.