https://mooseframework.inl.gov
Public Member Functions | Private Types | Private Attributes | List of all members
NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm > Class Template Reference

Adaptor functor for MOOSE RealVectorValue/RankTwoTensor residual and Jacobian types. More...

#include <NestedSolve.h>

Public Member Functions

 RealVectorEigenAdaptorFunctor (R &residual_lambda, J &jacobian_lambda)
 
int operator() (V &guess, V &residual)
 
int df (V &guess, typename NestedSolveTempl< is_ad >::template CorrespondingJacobian< V > &jacobian)
 
const RealgetResidualNorm ()
 

Private Types

using V = typename NestedSolveTempl< is_ad >::DynamicVector
 

Private Attributes

R & _residual_lambda
 
J & _jacobian_lambda
 
Real _residual_norm
 

Detailed Description

template<bool is_ad, typename R, typename J, bool store_residual_norm>
class NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >

Adaptor functor for MOOSE RealVectorValue/RankTwoTensor residual and Jacobian types.

Performs type conversion.

Definition at line 704 of file NestedSolve.h.

Member Typedef Documentation

◆ V

template<bool is_ad, typename R, typename J, bool store_residual_norm>
using NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::V = typename NestedSolveTempl<is_ad>::DynamicVector
private

Definition at line 706 of file NestedSolve.h.

Constructor & Destructor Documentation

◆ RealVectorEigenAdaptorFunctor()

template<bool is_ad, typename R, typename J, bool store_residual_norm>
NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::RealVectorEigenAdaptorFunctor ( R &  residual_lambda,
J &  jacobian_lambda 
)
inline

Definition at line 709 of file NestedSolve.h.

Member Function Documentation

◆ df()

template<bool is_ad, typename R, typename J, bool store_residual_norm>
int NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::df ( V guess,
typename NestedSolveTempl< is_ad >::template CorrespondingJacobian< V > &  jacobian 
)
inline

Definition at line 725 of file NestedSolve.h.

726  {
727  typename NestedSolveTempl<is_ad>::NSRealVectorValue guess_moose(guess(0), guess(1), guess(2));
728  typename NestedSolveTempl<is_ad>::NSRankTwoTensor jacobian_moose;
729  _jacobian_lambda(guess_moose, jacobian_moose);
730  jacobian =
731  Eigen::Map<typename NestedSolveTempl<is_ad>::DynamicMatrix>(&jacobian_moose(0, 0), 3, 3);
732  return 0;
733  }
Moose::GenericType< RankTwoTensor, is_ad > NSRankTwoTensor
Definition: NestedSolve.h:53
Moose::GenericType< RealVectorValue, is_ad > NSRealVectorValue
Definition: NestedSolve.h:52

◆ getResidualNorm()

template<bool is_ad, typename R, typename J, bool store_residual_norm>
const Real& NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::getResidualNorm ( )
inline

Definition at line 735 of file NestedSolve.h.

736  {
737  if constexpr (store_residual_norm)
738  return _residual_norm;
739  }

◆ operator()()

template<bool is_ad, typename R, typename J, bool store_residual_norm>
int NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::operator() ( V guess,
V residual 
)
inline

Definition at line 713 of file NestedSolve.h.

714  {
715  typename NestedSolveTempl<is_ad>::NSRealVectorValue guess_moose(guess(0), guess(1), guess(2));
716  typename NestedSolveTempl<is_ad>::NSRealVectorValue residual_moose;
717  _residual_lambda(guess_moose, residual_moose);
718  residual(0) = residual_moose(0);
719  residual(1) = residual_moose(1);
720  residual(2) = residual_moose(2);
721  if constexpr (store_residual_norm)
723  return 0;
724  }
static Real normSquare(const V &v)
Compute squared norm of v (dropping derivatives as this is only for convergence checking) ...
Definition: NestedSolve.h:555
Moose::GenericType< RealVectorValue, is_ad > NSRealVectorValue
Definition: NestedSolve.h:52

Member Data Documentation

◆ _jacobian_lambda

template<bool is_ad, typename R, typename J, bool store_residual_norm>
J& NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::_jacobian_lambda
private

◆ _residual_lambda

template<bool is_ad, typename R, typename J, bool store_residual_norm>
R& NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::_residual_lambda
private

◆ _residual_norm

template<bool is_ad, typename R, typename J, bool store_residual_norm>
Real NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::_residual_norm
private

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