46 #include <Teuchos_RCP.hpp> 48 #include <Shards_CellTopology.hpp> 50 #include <Intrepid_Basis.hpp> 51 #include <Intrepid_RealSpaceTools.hpp> 52 #include <Intrepid_Types.hpp> 64 Intrepid::Basis<Scalar, Intrepid::FieldContainer<double>>> &face_basis,
65 const Intrepid::FieldContainer<double> &point,
66 const Intrepid::FieldContainer<double> &face_nodes,
67 const Intrepid::FieldContainer<double> &face_node_normals )
68 : d_face_basis( face_basis )
70 , d_face_nodes( face_nodes )
71 , d_face_node_normals( face_node_normals )
73 d_space_dim = d_point.dimension( 0 );
74 d_topo_dim = d_space_dim - 1;
75 d_cardinality = d_face_basis->getCardinality();
78 Intrepid::FieldContainer<double>( d_cardinality, num_points );
79 d_grad_evals = Intrepid::FieldContainer<double>( d_cardinality, num_points,
81 d_eval_points = Intrepid::FieldContainer<double>( num_points, d_topo_dim );
89 const Intrepid::FieldContainer<double> &u )
92 for (
int i = 0; i < d_topo_dim; ++i )
94 d_eval_points( 0, i ) = u( 0, i );
98 d_face_basis->getValues( d_basis_evals, d_eval_points,
99 Intrepid::OPERATOR_VALUE );
102 d_face_basis->getValues( d_grad_evals, d_eval_points,
103 Intrepid::OPERATOR_GRAD );
111 const Intrepid::FieldContainer<double> &u,
112 Intrepid::FieldContainer<double> &F )
const 114 DTK_REQUIRE( 2 == u.rank() );
115 DTK_REQUIRE( 2 == F.rank() );
116 DTK_REQUIRE( F.dimension( 0 ) == u.dimension( 0 ) );
117 DTK_REQUIRE( F.dimension( 1 ) == u.dimension( 1 ) );
120 for (
int i = 0; i < d_space_dim; ++i )
122 F( 0, i ) = -d_point( i );
123 for (
int n = 0; n < d_cardinality; ++n )
125 F( 0, i ) += d_basis_evals( n, 0 ) *
126 ( d_face_nodes( n, i ) -
127 u( 0, d_topo_dim ) * d_face_node_normals( n, i ) );
137 const Intrepid::FieldContainer<double> &u,
138 Intrepid::FieldContainer<double> &J )
const 140 DTK_REQUIRE( 2 == u.rank() );
141 DTK_REQUIRE( 3 == J.rank() );
142 DTK_REQUIRE( J.dimension( 0 ) == u.dimension( 0 ) );
143 DTK_REQUIRE( J.dimension( 1 ) == u.dimension( 1 ) );
144 DTK_REQUIRE( J.dimension( 2 ) == u.dimension( 1 ) );
147 for (
int i = 0; i < d_space_dim; ++i )
150 for (
int j = 0; j < d_topo_dim; ++j )
153 for (
int n = 0; n < d_cardinality; ++n )
156 d_grad_evals( n, 0, j ) *
157 ( d_face_nodes( n, i ) -
158 u( 0, d_topo_dim ) * d_face_node_normals( n, i ) );
163 J( 0, i, d_space_dim - 1 ) = 0.0;
164 for (
int n = 0; n < d_cardinality; ++n )
166 J( 0, i, d_space_dim - 1 ) -=
167 d_basis_evals( n, 0 ) * d_face_node_normals( n, i );
180 const Intrepid::FieldContainer<double> &point,
181 const Intrepid::FieldContainer<double> &face_edge_nodes,
182 const Intrepid::FieldContainer<double> &face_edge_node_normals,
185 , d_face_edge_nodes( face_edge_nodes )
186 , d_face_edge_node_normals( face_edge_node_normals )
189 DTK_CHECK( d_point.dimension( 0 ) == d_face_edge_nodes.dimension( 1 ) );
190 DTK_CHECK( d_point.dimension( 0 ) ==
191 d_face_edge_node_normals.dimension( 1 ) );
192 DTK_CHECK( d_face_edge_nodes.dimension( 0 ) ==
193 d_face_edge_node_normals.dimension( 0 ) );
194 d_space_dim = point.dimension( 0 );
195 d_topo_dim = d_space_dim - 1;
199 Intrepid::FieldContainer<double>( d_cardinality, num_points );
200 d_grad_evals = Intrepid::FieldContainer<double>( d_cardinality, num_points,
202 d_eval_points = Intrepid::FieldContainer<double>( num_points, d_topo_dim );
203 d_proj_normal = Intrepid::FieldContainer<double>( d_space_dim );
207 d_face_normal_nodes = Intrepid::FieldContainer<double>( 2, d_space_dim );
208 for (
int j = 0; j < d_space_dim; ++j )
210 d_face_normal_nodes( 0, j ) =
211 face_edge_nodes( 0, j ) + d_c * d_face_edge_node_normals( 0, j );
212 d_face_normal_nodes( 1, j ) =
213 face_edge_nodes( 1, j ) + d_c * d_face_edge_node_normals( 1, j );
222 const Intrepid::FieldContainer<double> &u )
225 d_basis_evals( 0, 0 ) = ( 1.0 - u( 0, 0 ) ) * ( 1.0 - u( 0, 1 ) );
226 d_basis_evals( 1, 0 ) = u( 0, 0 ) * ( 1.0 - u( 0, 1 ) );
227 d_basis_evals( 2, 0 ) = u( 0, 0 ) * u( 0, 1 );
228 d_basis_evals( 3, 0 ) = ( 1.0 - u( 0, 0 ) ) * u( 0, 1 );
231 d_grad_evals( 0, 0, 0 ) = u( 0, 1 ) - 1.0;
232 d_grad_evals( 1, 0, 0 ) = 1.0 - u( 0, 1 );
233 d_grad_evals( 2, 0, 0 ) = u( 0, 1 );
234 d_grad_evals( 3, 0, 0 ) = -u( 0, 1 );
235 d_grad_evals( 0, 0, 1 ) = u( 0, 0 ) - 1.0;
236 d_grad_evals( 1, 0, 1 ) = -u( 0, 0 );
237 d_grad_evals( 2, 0, 1 ) = u( 0, 0 );
238 d_grad_evals( 3, 0, 1 ) = 1.0 - u( 0, 0 );
241 Intrepid::FieldContainer<double> v1( d_space_dim );
242 Intrepid::FieldContainer<double> v2( d_space_dim );
243 for (
int i = 0; i < d_space_dim; ++i )
245 v1( i ) = d_face_edge_nodes( 1, i ) - d_face_edge_nodes( 0, i ) +
246 d_c * u( 0, 1 ) * ( d_face_edge_node_normals( 1, i ) -
247 d_face_edge_node_normals( 0, i ) );
248 v2( i ) = u( 0, 0 ) * d_face_edge_node_normals( 1, i ) +
249 ( 1 - u( 0, 0 ) ) * d_face_edge_node_normals( 0, i );
251 Intrepid::RealSpaceTools<Scalar>::vecprod( d_proj_normal, v1, v2 );
259 const Intrepid::FieldContainer<double> &u,
260 Intrepid::FieldContainer<double> &F )
const 262 DTK_REQUIRE( 2 == u.rank() );
263 DTK_REQUIRE( 2 == F.rank() );
264 DTK_REQUIRE( F.dimension( 0 ) == u.dimension( 0 ) );
265 DTK_REQUIRE( F.dimension( 1 ) == u.dimension( 1 ) );
268 for (
int i = 0; i < d_space_dim; ++i )
270 F( 0, i ) = d_basis_evals( 0, 0 ) * d_face_edge_nodes( 0, i ) +
271 d_basis_evals( 1, 0 ) * d_face_edge_nodes( 1, i ) +
272 d_basis_evals( 2, 0 ) * d_face_normal_nodes( 1, i ) +
273 d_basis_evals( 3, 0 ) * d_face_normal_nodes( 0, i ) +
274 u( 0, 2 ) * d_proj_normal( i ) - d_point( i );
283 const Intrepid::FieldContainer<double> &u,
284 Intrepid::FieldContainer<double> &J )
const 286 DTK_REQUIRE( 2 == u.rank() );
287 DTK_REQUIRE( 3 == J.rank() );
288 DTK_REQUIRE( J.dimension( 0 ) == u.dimension( 0 ) );
289 DTK_REQUIRE( J.dimension( 1 ) == u.dimension( 1 ) );
290 DTK_REQUIRE( J.dimension( 2 ) == u.dimension( 1 ) );
293 for (
int i = 0; i < d_space_dim; ++i )
296 for (
int j = 0; j < d_topo_dim; ++j )
299 d_grad_evals( 0, 0, j ) * d_face_edge_nodes( 0, i ) +
300 d_grad_evals( 1, 0, j ) * d_face_edge_nodes( 1, i ) +
301 d_grad_evals( 2, 0, j ) * d_face_normal_nodes( 1, i ) +
302 d_grad_evals( 3, 0, j ) * d_face_normal_nodes( 0, i );
306 J( 0, i, d_space_dim - 1 ) = d_proj_normal( i );
318 const Intrepid::FieldContainer<double> &point,
319 const Intrepid::FieldContainer<double> &edge_nodes,
320 const Intrepid::FieldContainer<double> &edge_node_normals,
321 const Intrepid::FieldContainer<double> &edge_node_binormals )
323 , d_edge_nodes( edge_nodes )
324 , d_edge_node_normals( edge_node_normals )
325 , d_edge_node_binormals( edge_node_binormals )
327 DTK_CHECK( 1 == d_point.rank() );
328 DTK_CHECK( 2 == d_edge_nodes.rank() );
329 DTK_CHECK( 2 == d_edge_node_normals.rank() );
330 DTK_CHECK( 2 == d_edge_node_binormals.rank() );
331 DTK_CHECK( d_point.dimension( 0 ) == d_edge_nodes.dimension( 1 ) );
332 DTK_CHECK( d_point.dimension( 0 ) == d_edge_node_normals.dimension( 1 ) );
333 DTK_CHECK( d_edge_nodes.dimension( 0 ) ==
334 d_edge_node_normals.dimension( 0 ) );
335 DTK_CHECK( d_point.dimension( 0 ) == d_edge_node_binormals.dimension( 1 ) );
336 DTK_CHECK( d_edge_nodes.dimension( 0 ) ==
337 d_edge_node_binormals.dimension( 0 ) );
338 d_space_dim = d_point.dimension( 0 );
346 const Intrepid::FieldContainer<double> &u )
355 const Intrepid::FieldContainer<double> &u,
356 Intrepid::FieldContainer<double> &F )
const 358 DTK_REQUIRE( 2 == u.rank() );
359 DTK_REQUIRE( 2 == F.rank() );
360 DTK_REQUIRE( F.dimension( 0 ) == u.dimension( 0 ) );
361 DTK_REQUIRE( F.dimension( 1 ) == u.dimension( 1 ) );
364 for (
int i = 0; i < d_space_dim; ++i )
367 d_edge_nodes( 0, i ) +
368 u( 0, 0 ) * ( d_edge_nodes( 1, i ) - d_edge_nodes( 0, i ) ) +
369 u( 0, 1 ) * ( d_edge_node_normals( 0, i ) +
370 u( 0, 0 ) * ( d_edge_node_normals( 1, i ) -
371 d_edge_node_normals( 0, i ) ) ) +
372 u( 0, 2 ) * ( d_edge_node_binormals( 0, i ) +
373 u( 0, 0 ) * ( d_edge_node_binormals( 1, i ) -
374 d_edge_node_binormals( 0, i ) ) ) -
384 const Intrepid::FieldContainer<double> &u,
385 Intrepid::FieldContainer<double> &J )
const 387 DTK_REQUIRE( 2 == u.rank() );
388 DTK_REQUIRE( 3 == J.rank() );
389 DTK_REQUIRE( J.dimension( 0 ) == u.dimension( 0 ) );
390 DTK_REQUIRE( J.dimension( 1 ) == u.dimension( 1 ) );
391 DTK_REQUIRE( J.dimension( 2 ) == u.dimension( 1 ) );
394 for (
int i = 0; i < d_space_dim; ++i )
396 J( 0, i, 0 ) = ( d_edge_nodes( 1, i ) - d_edge_nodes( 0, i ) ) +
397 u( 0, 1 ) * ( d_edge_node_normals( 1, i ) -
398 d_edge_node_normals( 0, i ) ) +
399 u( 0, 2 ) * ( d_edge_node_binormals( 1, i ) -
400 d_edge_node_binormals( 0, i ) );
402 J( 0, i, 1 ) = d_edge_node_normals( 0, i ) +
403 u( 0, 0 ) * ( d_edge_node_normals( 1, i ) -
404 d_edge_node_normals( 0, i ) );
406 J( 0, i, 2 ) = d_edge_node_binormals( 0, i ) +
407 u( 0, 0 ) * ( d_edge_node_binormals( 1, i ) -
408 d_edge_node_binormals( 0, i ) );
ProjectPointToFaceNonlinearProblem(const Teuchos::RCP< Intrepid::Basis< Scalar, Intrepid::FieldContainer< double >>> &face_basis, const Intrepid::FieldContainer< double > &point, const Intrepid::FieldContainer< double > &face_nodes, const Intrepid::FieldContainer< double > &face_node_normals)
Constructor.
void evaluateResidual(const Intrepid::FieldContainer< double > &u, Intrepid::FieldContainer< double > &F) const
Evaluate the nonlinear residual.
void updateState(const Intrepid::FieldContainer< double > &u)
Update the state of the problem given the new solution vector.
void updateState(const Intrepid::FieldContainer< double > &u)
Update the state of the problem given the new solution vector.
void evaluateJacobian(const Intrepid::FieldContainer< double > &u, Intrepid::FieldContainer< double > &J) const
Evaluate the jacobian.
void evaluateResidual(const Intrepid::FieldContainer< double > &u, Intrepid::FieldContainer< double > &F) const
Evaluate the nonlinear residual.
ProjectPointFeatureToEdgeFeatureNonlinearProblem(const Intrepid::FieldContainer< double > &point, const Intrepid::FieldContainer< double > &edge_nodes, const Intrepid::FieldContainer< double > &edge_node_normals, const Intrepid::FieldContainer< double > &edge_node_binormals)
Constructor.
void evaluateJacobian(const Intrepid::FieldContainer< double > &u, Intrepid::FieldContainer< double > &J) const
Evaluate the jacobian.
void evaluateResidual(const Intrepid::FieldContainer< double > &u, Intrepid::FieldContainer< double > &F) const
Evaluate the nonlinear residual.
Assertions and Design-by-Contract for error handling.
void evaluateJacobian(const Intrepid::FieldContainer< double > &u, Intrepid::FieldContainer< double > &J) const
Evaluate the jacobian.
void updateState(const Intrepid::FieldContainer< double > &u)
Update the state of the problem given the new solution vector.
PointInFaceVolumeOfInfluenceNonlinearProblem(const Intrepid::FieldContainer< double > &point, const Intrepid::FieldContainer< double > &face_edge_nodes, const Intrepid::FieldContainer< double > &face_edge_node_normals, const double c)
Constructor.
Nonlinear problem definitions for projection primitives.