DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_ProjectionPrimitiveNonlinearProblems.cpp
1 //---------------------------------------------------------------------------//
2 /*
3  Copyright (c) 2012, Stuart R. Slattery
4  All rights reserved.
5 
6  Redistribution and use in source and binary forms, with or without
7  modification, are permitted provided that the following conditions are
8  met:
9 
10  *: Redistributions of source code must retain the above copyright
11  notice, this list of conditions and the following disclaimer.
12 
13  *: Redistributions in binary form must reproduce the above copyright
14  notice, this list of conditions and the following disclaimer in the
15  documentation and/or other materials provided with the distribution.
16 
17  *: Neither the name of the University of Wisconsin - Madison nor the
18  names of its contributors may be used to endorse or promote products
19  derived from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
24  A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25  HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26  SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
27  LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
28  DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
29  THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
31  OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 */
33 //---------------------------------------------------------------------------//
39 //---------------------------------------------------------------------------//
40 
41 #include <cmath>
42 
43 #include "DTK_DBC.hpp"
45 
46 #include <Teuchos_RCP.hpp>
47 
48 #include <Shards_CellTopology.hpp>
49 
50 #include <Intrepid_Basis.hpp>
51 #include <Intrepid_RealSpaceTools.hpp>
52 #include <Intrepid_Types.hpp>
53 
54 namespace DataTransferKit
55 {
56 //---------------------------------------------------------------------------//
57 // ProjectPointToFaceNonlinearProblem Implementation.
58 //---------------------------------------------------------------------------//
63  const Teuchos::RCP<
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 )
69  , d_point( point )
70  , d_face_nodes( face_nodes )
71  , d_face_node_normals( face_node_normals )
72 {
73  d_space_dim = d_point.dimension( 0 );
74  d_topo_dim = d_space_dim - 1;
75  d_cardinality = d_face_basis->getCardinality();
76  int num_points = 1;
77  d_basis_evals =
78  Intrepid::FieldContainer<double>( d_cardinality, num_points );
79  d_grad_evals = Intrepid::FieldContainer<double>( d_cardinality, num_points,
80  d_topo_dim );
81  d_eval_points = Intrepid::FieldContainer<double>( num_points, d_topo_dim );
82 }
83 
84 //---------------------------------------------------------------------------//
89  const Intrepid::FieldContainer<double> &u )
90 {
91  // Extract the current natural coordinates from the solution vector.
92  for ( int i = 0; i < d_topo_dim; ++i )
93  {
94  d_eval_points( 0, i ) = u( 0, i );
95  }
96 
97  // Evaluate the basis at the current natural coordinates.
98  d_face_basis->getValues( d_basis_evals, d_eval_points,
99  Intrepid::OPERATOR_VALUE );
100 
101  // Evaluate the basis gradient at the current natural coordinates.
102  d_face_basis->getValues( d_grad_evals, d_eval_points,
103  Intrepid::OPERATOR_GRAD );
104 }
105 
106 //---------------------------------------------------------------------------//
111  const Intrepid::FieldContainer<double> &u,
112  Intrepid::FieldContainer<double> &F ) const
113 {
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 ) );
118 
119  // Build the nonlinear residual.
120  for ( int i = 0; i < d_space_dim; ++i )
121  {
122  F( 0, i ) = -d_point( i );
123  for ( int n = 0; n < d_cardinality; ++n )
124  {
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 ) );
128  }
129  }
130 }
131 
132 //---------------------------------------------------------------------------//
137  const Intrepid::FieldContainer<double> &u,
138  Intrepid::FieldContainer<double> &J ) const
139 {
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 ) );
145 
146  // Build the Jacobian.
147  for ( int i = 0; i < d_space_dim; ++i )
148  {
149  // Start with the basis gradient contributions.
150  for ( int j = 0; j < d_topo_dim; ++j )
151  {
152  J( 0, i, j ) = 0.0;
153  for ( int n = 0; n < d_cardinality; ++n )
154  {
155  J( 0, i, j ) +=
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 ) );
159  }
160  }
161 
162  // Then add the point normal contribution.
163  J( 0, i, d_space_dim - 1 ) = 0.0;
164  for ( int n = 0; n < d_cardinality; ++n )
165  {
166  J( 0, i, d_space_dim - 1 ) -=
167  d_basis_evals( n, 0 ) * d_face_node_normals( n, i );
168  }
169  }
170 }
171 
172 //---------------------------------------------------------------------------//
173 // PointInFaceVolumeOfInfluenceNonlinearProblem Implementation.
174 //---------------------------------------------------------------------------//
180  const Intrepid::FieldContainer<double> &point,
181  const Intrepid::FieldContainer<double> &face_edge_nodes,
182  const Intrepid::FieldContainer<double> &face_edge_node_normals,
183  const double c )
184  : d_point( point )
185  , d_face_edge_nodes( face_edge_nodes )
186  , d_face_edge_node_normals( face_edge_node_normals )
187  , d_c( c )
188 {
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;
196  d_cardinality = 4;
197  int num_points = 1;
198  d_basis_evals =
199  Intrepid::FieldContainer<double>( d_cardinality, num_points );
200  d_grad_evals = Intrepid::FieldContainer<double>( d_cardinality, num_points,
201  d_topo_dim );
202  d_eval_points = Intrepid::FieldContainer<double>( num_points, d_topo_dim );
203  d_proj_normal = Intrepid::FieldContainer<double>( d_space_dim );
204 
205  // Compute the two extra psuedo-nodes along the node normal
206  // directions.
207  d_face_normal_nodes = Intrepid::FieldContainer<double>( 2, d_space_dim );
208  for ( int j = 0; j < d_space_dim; ++j )
209  {
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 );
214  }
215 }
216 
217 //---------------------------------------------------------------------------//
222  const Intrepid::FieldContainer<double> &u )
223 {
224  // Evaluate the basis at the current natural coordinates.
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 );
229 
230  // Evaluate the basis gradient at the current natural coordinates.
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 );
239 
240  // Compute the current projection normal direction.
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 )
244  {
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 );
250  }
251  Intrepid::RealSpaceTools<Scalar>::vecprod( d_proj_normal, v1, v2 );
252 }
253 
254 //---------------------------------------------------------------------------//
259  const Intrepid::FieldContainer<double> &u,
260  Intrepid::FieldContainer<double> &F ) const
261 {
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 ) );
266 
267  // Build the nonlinear residual.
268  for ( int i = 0; i < d_space_dim; ++i )
269  {
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 );
275  }
276 }
277 
278 //---------------------------------------------------------------------------//
283  const Intrepid::FieldContainer<double> &u,
284  Intrepid::FieldContainer<double> &J ) const
285 {
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 ) );
291 
292  // Build the Jacobian.
293  for ( int i = 0; i < d_space_dim; ++i )
294  {
295  // Start with the basis gradient contributions.
296  for ( int j = 0; j < d_topo_dim; ++j )
297  {
298  J( 0, i, 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 );
303  }
304 
305  // Then add the point normal contribution.
306  J( 0, i, d_space_dim - 1 ) = d_proj_normal( i );
307  }
308 }
309 
310 //---------------------------------------------------------------------------//
311 // ProjectPointFeatureToEdgeFeatureNonlinearProblem Implementation.
312 //---------------------------------------------------------------------------//
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 )
322  : d_point( point )
323  , d_edge_nodes( edge_nodes )
324  , d_edge_node_normals( edge_node_normals )
325  , d_edge_node_binormals( edge_node_binormals )
326 {
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 );
339 }
340 
341 //---------------------------------------------------------------------------//
346  const Intrepid::FieldContainer<double> &u )
347 { /* ... */
348 }
349 
350 //---------------------------------------------------------------------------//
355  const Intrepid::FieldContainer<double> &u,
356  Intrepid::FieldContainer<double> &F ) const
357 {
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 ) );
362 
363  // Build the nonlinear residual.
364  for ( int i = 0; i < d_space_dim; ++i )
365  {
366  F( 0, 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 ) ) ) -
375  d_point( i );
376  }
377 }
378 
379 //---------------------------------------------------------------------------//
384  const Intrepid::FieldContainer<double> &u,
385  Intrepid::FieldContainer<double> &J ) const
386 {
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 ) );
392 
393  // Build the Jacobian.
394  for ( int i = 0; i < d_space_dim; ++i )
395  {
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 ) );
401 
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 ) );
405 
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 ) );
409  }
410 }
411 
412 //---------------------------------------------------------------------------//
413 
414 } // end namespace DataTransferKit
415 
416 //---------------------------------------------------------------------------//
417 // end DTK_ProjectionPrimitiveNonlinearProblems.cpp
418 //---------------------------------------------------------------------------//
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.
DTK_BasicEntitySet.cpp.