Trilinos::Sacado example
Public Member Functions | Public Attributes | Static Public Attributes | List of all members
Sacado_Wrapper::SymTensor< dim > Class Template Reference

#include <Sacado_Wrapper.h>

Inheritance diagram for Sacado_Wrapper::SymTensor< dim >:
Inheritance graph
[legend]
Collaboration diagram for Sacado_Wrapper::SymTensor< dim >:
Collaboration graph
[legend]

Public Member Functions

 SymTensor ()
 
void init (SymmetricTensor< 2, dim > &tensor_double)
 
void set_dofs (unsigned int nbr_total_dofs=n_dofs)
 
void get_tangent (SymmetricTensor< 4, dim > &Tangent, SymmetricTensor< 2, dim, fad_double > &sigma)
 
void get_tangent (SymmetricTensor< 2, dim > &Tangent, fad_double &argument)
 
void get_values (SymmetricTensor< 2, dim > &tensor_double)
 
SymmetricTensor< 2, dim > get_value ()
 

Public Attributes

std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
 
unsigned int start_index = 0
 

Static Public Attributes

static const unsigned int n_dofs = ((dim==2)?3:6)
 

Constructor & Destructor Documentation

◆ SymTensor()

template<int dim>
Sacado_Wrapper::SymTensor< dim >::SymTensor ( )
inline
70  {
71  get_index_map<dim>( std_map_indicies );
72  }
std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
Definition: Sacado_Wrapper.h:75

Member Function Documentation

◆ get_tangent() [1/2]

template<int dim>
void Sacado_Wrapper::SymTensor< dim >::get_tangent ( SymmetricTensor< 4, dim > &  Tangent,
SymmetricTensor< 2, dim, fad_double > &  sigma 
)

Referenced by sacado_test_10(), sacado_test_3B(), and sacado_test_4().

130  {
131  for ( unsigned int i=0; i<dim; ++i)
132  for ( unsigned int j=0; j<dim; ++j )
133  {
134  double *derivs = &sigma[i][j].fastAccessDx(0); // Access the derivatives of the (i,j)-th component of \a sigma
135 
136  // We loop over all the dofs. To be able to use this independent of the chosen dimension \a dim, we use a ternary operator
137  // to decide whether we have to loop over 6 derivatives or just 3.
138  for(unsigned int x=start_index;x<(start_index+n_dofs);++x)
139  {
140  unsigned int k=std_map_indicies[x].first;
141  unsigned int l=std_map_indicies[x].second;
142 
143  if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
144  {
145  Tangent[i][j][k][l] = 0.5*derivs[x];
146  Tangent[i][j][l][k] = 0.5*derivs[x];
147  }
148  else
149  Tangent[i][j][k][l] = derivs[x];
150  }
151  }
152  }
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:78
unsigned int start_index
Definition: Sacado_Wrapper.h:77
std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
Definition: Sacado_Wrapper.h:75

◆ get_tangent() [2/2]

template<int dim>
void Sacado_Wrapper::SymTensor< dim >::get_tangent ( SymmetricTensor< 2, dim > &  Tangent,
fad_double argument 
)
156  {
157  double *derivs = &argument.fastAccessDx(0); // Access derivatives
158  for(unsigned int x=start_index;x<(start_index+n_dofs);++x)
159  {
160  unsigned int k=std_map_indicies[x].first;
161  unsigned int l=std_map_indicies[x].second;
162 
163  Tangent[k][l] = derivs[x]; // ToDo: check whether the 0.5* is necessary here too
164 
165  // Correct the off-diagonal terms by the factor of 0.5
166  if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
167  Tangent[k][l] *= 0.5; // ToDo: check whether the 0.5* is necessary here too
168  }
169  }
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:78
unsigned int start_index
Definition: Sacado_Wrapper.h:77
std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
Definition: Sacado_Wrapper.h:75

◆ get_value()

template<int dim>
SymmetricTensor< 2, dim > Sacado_Wrapper::SymTensor< dim >::get_value ( )

References Sacado_Wrapper::SymTensor< dim >::get_values().

182  {
183  SymmetricTensor<2,dim> tmp;
184  (*this).get_values(tmp);
185  return tmp;
186  }

◆ get_values()

template<int dim>
void Sacado_Wrapper::SymTensor< dim >::get_values ( SymmetricTensor< 2, dim > &  tensor_double)

Referenced by Sacado_Wrapper::SymTensor< dim >::get_value().

174  {
175  for ( unsigned int i=0; i<dim; ++i)
176  for ( unsigned int j=0; j<dim; ++j )
177  tensor_double[i][j] = ((*this)[i][j]).val();
178  }

◆ init()

template<int dim>
void Sacado_Wrapper::SymTensor< dim >::init ( SymmetricTensor< 2, dim > &  tensor_double)

Referenced by sacado_test_10(), sacado_test_3B(), and sacado_test_4().

99  {
100  for ( unsigned int i=0; i<dim; ++i)
101  for ( unsigned int j=0; j<dim; ++j )
102  (*this)[i][j] = tensor_double[i][j];
103  }

◆ set_dofs()

template<int dim>
void Sacado_Wrapper::SymTensor< dim >::set_dofs ( unsigned int  nbr_total_dofs = n_dofs)

Referenced by sacado_test_10(), sacado_test_3B(), and Sacado_Wrapper::DoFs_summary< dim >::set_dofs().

111  {
112  // Instead of calling the *.diff(*) on the components one-by-one we could also use the following for-loop, so
113  // we also use the map to set the dofs
114  for ( unsigned int x=start_index; x<(start_index+n_dofs); ++x )
115  {
116  unsigned int i=std_map_indicies[x].first;
117  unsigned int j=std_map_indicies[x].second;
118  (*this)[i][j].diff(x,nbr_total_dofs);
119  }
120  }
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:78
unsigned int start_index
Definition: Sacado_Wrapper.h:77
std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
Definition: Sacado_Wrapper.h:75

Member Data Documentation

◆ n_dofs

template<int dim>
const unsigned int Sacado_Wrapper::SymTensor< dim >::n_dofs = ((dim==2)?3:6)
static

◆ start_index

template<int dim>
unsigned int Sacado_Wrapper::SymTensor< dim >::start_index = 0

◆ std_map_indicies

template<int dim>
std::map<unsigned int,std::pair<unsigned int,unsigned int> > Sacado_Wrapper::SymTensor< dim >::std_map_indicies

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