41 #include "DTK_IntrepidShapeFunction.hpp" 51 const shards::CellTopology &topology,
52 const Teuchos::ArrayView<const double> &reference_point,
53 Teuchos::Array<double> &values )
const 56 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
57 basis = getIntrepidBasis( topology );
60 Teuchos::Array<int> point_dims( 2 );
62 point_dims[1] = reference_point.size();
63 Intrepid::FieldContainer<double> point_container(
64 point_dims, const_cast<double *>( reference_point.getRawPtr() ) );
67 values.resize( basis->getCardinality() );
68 Teuchos::Array<int> value_dims( 2 );
69 value_dims[0] = basis->getCardinality();
71 Intrepid::FieldContainer<double> value_container( value_dims,
75 basis->getValues( value_container, point_container,
76 Intrepid::OPERATOR_VALUE );
83 const shards::CellTopology &topology,
84 const Teuchos::ArrayView<const double> &reference_point,
85 Teuchos::Array<Teuchos::Array<double>> &gradients )
const 88 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
89 basis = getIntrepidBasis( topology );
92 int space_dim = reference_point.size();
93 Teuchos::Array<int> point_dims( 2 );
95 point_dims[1] = space_dim;
96 Intrepid::FieldContainer<double> point_container(
97 point_dims, const_cast<double *>( reference_point.getRawPtr() ) );
100 int cardinality = basis->getCardinality();
101 Intrepid::FieldContainer<double> grad_container( cardinality, 1,
103 basis->getValues( grad_container, point_container,
104 Intrepid::OPERATOR_GRAD );
107 gradients.resize( cardinality );
108 for (
int n = 0; n < cardinality; ++n )
110 gradients[n].resize( space_dim );
111 for (
int d = 0; d < space_dim; ++d )
113 gradients[n][d] = grad_container( n, 0, d );
120 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
121 IntrepidShapeFunction::getIntrepidBasis(
122 const shards::CellTopology &topology )
const 125 unsigned basis_key = topology.getKey();
126 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
128 if ( d_basis.count( basis_key ) )
130 basis = d_basis.find( basis_key )->second;
134 basis = IntrepidBasisFactory::create( topology );
135 d_basis.emplace( basis_key, basis );
Assertions and Design-by-Contract for error handling.
void evaluateGradient(const shards::CellTopology &topology, const Teuchos::ArrayView< const double > &reference_point, Teuchos::Array< Teuchos::Array< double >> &gradients) const
Given an topology and a reference point, evaluate the gradient of the shape function of the topology ...
Factory for Intrepid basis functions.
void evaluateValue(const shards::CellTopology &topology, const Teuchos::ArrayView< const double > &reference_point, Teuchos::Array< double > &values) const
Given an topology and a reference point, evaluate the shape function of the topology at that point...