Trilinos::Sacado example
Sacado_Wrapper.h
Go to the documentation of this file.
1 #ifndef Sacado_Wrapper_H
2 #define Sacado_Wrapper_H
3 
4 // @section includes Include Files
5 // The data type SymmetricTensor and some related operations, such as trace, symmetrize, deviator, ... for tensor calculus
6 #include <deal.II/base/symmetric_tensor.h>
7 
8 // @todo Check whether the following three headers are needed at all
9 #include <iostream>
10 #include <fstream>
11 #include <cmath>
12 
13 // Sacado (from Trilinos, data types, operations, ...)
14 #include <Sacado.hpp>
15 
16 using namespace dealii;
17 
18 using fad_double = Sacado::Fad::DFad<double>; // this data type now represents a double, but also contains the derivative of this variable with respect to the defined dofs (set via command *.diff(*))
19 typedef Sacado::Fad::DFad<double> DFadType;
20 /*
21  * A function to set up the map that relates the indices (e.g. first entry in vector is the (0,0) component, whereas the second is the (0,1) component
22  */
23 template<int dim>
24 void get_index_map( std::map<unsigned int,std::pair<unsigned int,unsigned int>> &std_map_indicies )
25 {
26  std::pair<unsigned int, unsigned int> tmp_pair;
27 
28  switch ( dim )
29  {
30  case 2:
31  tmp_pair.first=0; tmp_pair.second=0;
32  std_map_indicies[0] = tmp_pair;
33 
34  tmp_pair.first=0; tmp_pair.second=1;
35  std_map_indicies[1] = tmp_pair;
36 
37  tmp_pair.first=1; tmp_pair.second=1;
38  std_map_indicies[2] = tmp_pair;
39  break;
40  case 3:
41  tmp_pair.first=0; tmp_pair.second=0;
42  std_map_indicies[0] = tmp_pair;
43 
44  tmp_pair.first=0; tmp_pair.second=1;
45  std_map_indicies[1] = tmp_pair;
46 
47  tmp_pair.first=0; tmp_pair.second=2;
48  std_map_indicies[2] = tmp_pair;
49 
50  tmp_pair.first=1; tmp_pair.second=1;
51  std_map_indicies[3] = tmp_pair;
52 
53  tmp_pair.first=1; tmp_pair.second=2;
54  std_map_indicies[4] = tmp_pair;
55 
56  tmp_pair.first=2; tmp_pair.second=2;
57  std_map_indicies[5] = tmp_pair;
58  break;
59  }
60 }
61 
62 
63 namespace Sacado_Wrapper
64 {
65  template <int dim>
66  class SymTensor: public SymmetricTensor<2,dim, fad_double>
67  {
68  public:
70  {
71  get_index_map<dim>( std_map_indicies );
72  }
73 
74  // To simplify the access to the dofs we define a map that relate the components of our strain tensor to the dof-nbr
75  std::map<unsigned int,std::pair<unsigned int,unsigned int>> std_map_indicies;
76 
77  unsigned int start_index = 0;
78  static const unsigned int n_dofs = ((dim==2)?3:6); // ToDo: Find a way to use eps.n_independent_components(); that gives ((dim==2)?3:6) for a SymmetricTensor
79 
80  void init ( SymmetricTensor<2,dim> &tensor_double );
81 
82  void set_dofs( unsigned int nbr_total_dofs=n_dofs );
83 
84  void get_tangent (SymmetricTensor<4,dim> &Tangent, SymmetricTensor<2,dim, fad_double> &sigma);
85 
86  void get_tangent( SymmetricTensor<2,dim> &Tangent, fad_double &argument );
87 
88  void get_values ( SymmetricTensor<2,dim> &tensor_double );
89 
90  SymmetricTensor<2,dim> get_value ();
91  };
92 
93 
94  /*
95  * Initialization of the SymTensor data type with the values from a normal double SymmetricTensor
96  */
97  template<int dim>
98  void SymTensor<dim>::init( SymmetricTensor<2,dim> &tensor_double )
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  }
104 
105 
106  /*
107  * Set the dofs as the components of the SymTensor data type
108  */
109  template<int dim>
110  void SymTensor<dim>::set_dofs( unsigned int nbr_total_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  }
121 
122 
123  /*
124  * Assemble the tangent based on the derivatives saved in the argument sigma
125  * @param Tangent The desired tangent that will be computed as d_sigma/d_this, where "this" could be your strain tensor
126  * @param sigma The input argument used to extract the derivatives
127  */
128  template<int dim>
129  void SymTensor<dim>::get_tangent( SymmetricTensor<4,dim> &Tangent, SymmetricTensor<2,dim, fad_double> &sigma )
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  }
153 
154  template<int dim>
155  void 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  }
170 
171 
172  template<int dim>
173  void SymTensor<dim>::get_values ( SymmetricTensor<2,dim> &tensor_double )
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  }
179 
180  template<int dim>
181  SymmetricTensor<2,dim> SymTensor<dim>::get_value ()
182  {
183  SymmetricTensor<2,dim> tmp;
184  (*this).get_values(tmp);
185  return tmp;
186  }
187 
188 
189  //###########################################################################################################//
190 
191 
192  template <int dim>
193  class SymTensor2: public SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> >
194  {
195  public:
197  {
198  get_index_map<dim>( std_map_indicies );
199  }
200 
201  // To simplify the access to the dofs we define a map that relate the components of our strain tensor to the dof-nbr
202  std::map<unsigned int,std::pair<unsigned int,unsigned int>> std_map_indicies;
203 
204  unsigned int start_index = 0;
205  static const unsigned int n_dofs = ((dim==2)?3:6); // ToDo: Find a way to use eps.n_independent_components(); that gives ((dim==2)?3:6) for a SymmetricTensor
206 
207  void init_set_dofs ( SymmetricTensor<2,dim> &tensor_double, const unsigned int nbr_total_dofs=n_dofs );
208 
209  void get_tangent (SymmetricTensor<4,dim> &Tangent, SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &sigma);
210 
211  void get_tangent( SymmetricTensor<2,dim> &Tangent, Sacado::Fad::DFad<DFadType> &argument );
212 
213  void get_tangent( SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &Tangent, Sacado::Fad::DFad<DFadType> &argument );
214 
215  void get_curvature( SymmetricTensor<4,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument );
216 
217 //
218 // void get_values ( SymmetricTensor<2,dim> &tensor_double );
219  };
220 
221  /*
222  * Initialization of the \a SymTensor2 data type with the values from a normal double \a SymmetricTensor
223  */
224  template<int dim>
225  void SymTensor2<dim>::init_set_dofs( SymmetricTensor<2,dim> &tensor_double, const unsigned int nbr_total_dofs )
226  {
227  for ( unsigned int x=0; x<n_dofs; ++x )
228  {
229  unsigned int i=std_map_indicies[x].first;
230  unsigned int j=std_map_indicies[x].second;
231  ((*this)[i][j]).diff( x, nbr_total_dofs); // set up the "inner" derivatives
232  ((*this)[i][j]).val() = fad_double(nbr_total_dofs, x, tensor_double[i][j]); // set up the "outer" derivatives
233  }
234  }
235 
236 
237  template<int dim>
238  void SymTensor2<dim>::get_tangent( SymmetricTensor<2,dim> &Tangent, Sacado::Fad::DFad<DFadType> &argument )
239  {
240  for ( unsigned int x=0; x<n_dofs; ++x )
241  {
242  unsigned int i=std_map_indicies[x].first;
243  unsigned int j=std_map_indicies[x].second;
244  if ( i!=j )
245  Tangent[i][j] = 0.5 * argument.dx(x).val();
246  else
247  Tangent[i][j] = argument.dx(x).val();
248  }
249  }
250 
251  template<int dim>
252  void SymTensor2<dim>::get_tangent( SymmetricTensor<4,dim> &Tangent, SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &argument )
253  {
254  for(unsigned int x=0;x<n_dofs;++x)
255  for(unsigned int y=0;y<n_dofs;++y)
256  {
257  const unsigned int i=std_map_indicies[y].first;
258  const unsigned int j=std_map_indicies[y].second;
259  const unsigned int k=std_map_indicies[x].first;
260  const unsigned int l=std_map_indicies[x].second;
261 
262  if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
263  Tangent[i][j][k][l] = 0.5 * argument[k][l].dx(y).val(); // ToDo: check this on paper (was more like a gut feeling)
264  else
265  Tangent[i][j][k][l] = argument[k][l].dx(y).val();
266  }
267  }
268 
269 
270  /*
271  * Compute the tangent still containing all the second derivatives
272  */
273  template<int dim>
274  void SymTensor2<dim>::get_tangent( SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &Tangent, Sacado::Fad::DFad<DFadType> &argument )
275  {
276  for ( unsigned int x=0; x<n_dofs; ++x )
277  {
278  unsigned int i=std_map_indicies[x].first;
279  unsigned int j=std_map_indicies[x].second;
280  if ( i!=j )
281  Tangent[i][j] = 0.5 * argument.dx(x);
282  else
283  Tangent[i][j] = argument.dx(x);
284  }
285  }
286 
287 
288  template<int dim>
289  void SymTensor2<dim>::get_curvature( SymmetricTensor<4,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument )
290  {
291  for(unsigned int x=0;x<n_dofs;++x)
292  for(unsigned int y=0;y<n_dofs;++y)
293  {
294  const unsigned int i=std_map_indicies[y].first;
295  const unsigned int j=std_map_indicies[y].second;
296  const unsigned int k=std_map_indicies[x].first;
297  const unsigned int l=std_map_indicies[x].second;
298 
299  double deriv = argument.dx(x).dx(y); // Access the derivatives of the (i,j)-th component of \a sigma
300 
301  if ( k!=l && i!=j )
302  Curvature[i][j][k][l] = 0.25* deriv;
303  else if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
304  {
305  Curvature[i][j][k][l] = 0.5*deriv;
306  Curvature[i][j][l][k] = 0.5*deriv;
307  }
308  else
309  Curvature[i][j][k][l] = deriv;
310  }
311  }
312 
313  //###########################################################################################################//
314 
315 
316 
317  /*
318  * The same as the class above, just for data type fad_double
319  */
320  template<int dim>
321  class SW_double: public fad_double
322  {
323  public:
324 // SW_double( double &sdf );
325 // :
326 // (*this)(sdf)
327 // {
328 // }
329 
330  // Assignment operator: \n
331  // According to https://stackoverflow.com/questions/31029007/c-assignment-operator-in-derived-class
332  // the operator= is not derived from the base class fad_double and needs to be set explicitly:
333  SW_double & operator=(double double_init) { fad_double::operator =( double_init ) ;return *this;}
334  SW_double & operator=(fad_double fad_assignment) { fad_double::operator =( fad_assignment ) ;return *this;}
335 
336  static const unsigned int n_dofs = 1; // a double is just a single number, so it represents a single dof
337  unsigned int start_index = 0;
338 
339  void init ( const double &double_init );
340 
341  void set_dofs ( unsigned int nbr_total_dofs=n_dofs );
342 
343  void get_tangent (SymmetricTensor<2,dim> &Tangent, SymmetricTensor<2,dim, fad_double> &sigma);
344 
345  void get_tangent ( double &Tangent, fad_double &argument );
346 
347  void get_values ( double &return_double );
348  };
349 
350  template<int dim>
351  void SW_double<dim>::init ( const double &double_init )
352  {
353  (*this) = double_init;
354  }
355 
356  template<int dim>
357  void SW_double<dim>::set_dofs ( unsigned int nbr_total_dofs )
358  {
359  (*this).diff( this->start_index, nbr_total_dofs );
360  }
361 
362  template<int dim>
363  void SW_double<dim>::get_tangent (SymmetricTensor<2,dim> &Tangent, SymmetricTensor<2,dim, fad_double> &sigma)
364  {
365  // reassemble the tangent as a SECOND order tensor
366  for ( unsigned int i=0; i<dim; ++i)
367  for ( unsigned int j=0; j<dim; ++j )
368  {
369  double *derivs = &sigma[i][j].fastAccessDx(0); // Access derivatives
370  Tangent[i][j] = derivs[ this->start_index ]; // ToDo: check whether the 0.5* is necessary here too
371  }
372  }
373 
374  template<int dim>
375  void SW_double<dim>::get_tangent ( double &Tangent, fad_double &argument )
376  {
377  double *derivs = &argument.fastAccessDx(0);
378  Tangent = derivs[ this->start_index ];
379  }
380 
381  template<int dim>
382  void SW_double<dim>::get_values ( double &return_double )
383  {
384  return_double = (*this).val();
385  }
386 
387  //###########################################################################################################//
388 
389 
390 
391  /*
392  * The same as the class above, just for data type fad_double
393  */
394  template<int dim>
395  class SW_double2: public Sacado::Fad::DFad<DFadType>
396  {
397  public:
398 // SW_double2( double &sdf );
399 // :
400 // (*this)(sdf)
401 // {
402 // }
403 
404  // Assignment operator: \n
405  // According to https://stackoverflow.com/questions/31029007/c-assignment-operator-in-derived-class
406  // the operator= is not derived from the base class fad_double and needs to be set explicitly:
407  SW_double2 & operator=(double double_init) { Sacado::Fad::DFad<DFadType>::operator =( double_init ) ;return *this;}
408  SW_double2 & operator=(Sacado::Fad::DFad<DFadType> fad_assignment) { Sacado::Fad::DFad<DFadType>::operator =( fad_assignment ) ;return *this;}
409 
410  static const unsigned int n_dofs = 1; // a double is just a single number, so it represents a single dof
411  unsigned int start_index = 0;
412 
413  void init_set_dofs ( const double &double_init, unsigned int nbr_total_dofs=n_dofs );
414 
415  void get_tangent (double &Tangent, Sacado::Fad::DFad<DFadType> &argument);
416  void get_tangent (SymmetricTensor<2,dim> &Tangent, SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &argument, SymTensor2<dim> &eps);
417 
418 
419  void get_curvature (double &Curvature, Sacado::Fad::DFad<DFadType> &argument);
420 
421  void get_curvature (SymmetricTensor<2,dim> &Curvature, SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &argument, SymTensor2<dim> &eps );
422  void get_curvature (SymmetricTensor<2,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument, SymTensor2<dim> &eps );
423 
424 // void get_values ( double &return_double );
425  };
426 
427  template<int dim>
428  void SW_double2<dim>::init_set_dofs ( const double &double_init, unsigned int nbr_total_dofs )
429  {
430  (*this).diff( start_index , nbr_total_dofs );
431  (*this).val() = fad_double(nbr_total_dofs, start_index, double_init);
432  }
433 
434 
435  template<int dim>
436  void SW_double2<dim>::get_tangent (double &Tangent, Sacado::Fad::DFad<DFadType> &argument)
437  {
438  Tangent = argument.dx(this->start_index).val();
439  }
440 
441 
442  template<int dim>
443  void SW_double2<dim>::get_tangent (SymmetricTensor<2,dim> &Tangent, SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &argument, SymTensor2<dim> &eps)
444  {
445  for(unsigned int x=0;x<eps.n_dofs;++x)
446  {
447  const unsigned int i=eps.std_map_indicies[x].first;
448  const unsigned int j=eps.std_map_indicies[x].second;
449  // ToDo: find a better way to loop over the indices than using eps as input argument
450  Tangent[i][j] = argument[i][j].dx(start_index).val();
451  }
452  }
453 
454 
455  template<int dim>
456  void SW_double2<dim>::get_curvature (double &Curvature, Sacado::Fad::DFad<DFadType> &argument)
457  {
458  Curvature = argument.dx(this->start_index).dx(this->start_index);
459  }
460 
461 
462  /*
463  * Compute Curvature d_argument / d_phi, where argument is already the derivative with respect to d_eps
464  */
465  template<int dim>
466  void SW_double2<dim>::get_curvature (SymmetricTensor<2,dim> &Curvature, SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > &argument, SymTensor2<dim> &eps )
467  {
468  for(unsigned int x=0;x<eps.n_dofs;++x)
469  {
470  const unsigned int i=eps.std_map_indicies[x].first;
471  const unsigned int j=eps.std_map_indicies[x].second;
472 
473 // if ( i!=j )
474 // Curvature[i][j] = 0.5 * argument[i][j].val().dx(start_index); // ToDo: check the factor 0.5
475 // else
476  Curvature[i][j] = argument[i][j].val().dx(start_index);
477  }
478  }
479 
480  /*
481  * Compute Curvature d2_argument / d_phi_d_eps
482  */
483  template<int dim>
484  void SW_double2<dim>::get_curvature (SymmetricTensor<2,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument, SymTensor2<dim> &eps )
485  {
486  for(unsigned int x=0;x<eps.n_dofs;++x)
487  {
488  const unsigned int i=eps.std_map_indicies[x].first;
489  const unsigned int j=eps.std_map_indicies[x].second;
490 
491  if ( i!=j )
492  Curvature[i][j] = 0.5 * argument.dx(start_index).dx(x); // ToDo: check the factor 0.5
493  else
494  Curvature[i][j] = argument.dx(start_index).dx(x);
495  }
496  }
497 
498 
499 
500  //###########################################################################################################//
501 
502 
503 
504  template<int dim>
506  {
507  public:
508  void set_dofs( SymTensor<dim> &eps, SW_double<dim> &double_arg );
509  void init_set_dofs( SymTensor2<dim> &eps, SymmetricTensor<2,dim> &eps_init, SW_double2<dim> &double_arg, double &double_init );
510 
511  void get_curvature( SymmetricTensor<2,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument, SymTensor2<dim> &eps, SW_double2<dim> &double_arg );
512  void get_curvature( SymmetricTensor<2,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument, SW_double2<dim> &double_arg, SymTensor2<dim> &eps );
513 
514 
515  void set_dofs( SymTensor<dim> &eps, SW_double<dim> &double_arg1, SW_double<dim> &double_arg2, SW_double<dim> &double_arg3 );
516  // e.g. for strain, gamma_p and gamma_d
517  void set_dofs( SymTensor<dim> &eps, SW_double<dim> &double_arg1, SW_double<dim> &double_arg2 );
518  };
519 
520  template<int dim>
522  {
523  const unsigned int nbr_total_dofs = eps.n_dofs + double_arg.n_dofs;
524 
525  eps.start_index = 0;
526  double_arg.start_index = eps.n_dofs;
527 
528  eps.set_dofs( nbr_total_dofs );
529  double_arg.set_dofs( nbr_total_dofs );
530  }
531 
532  template<int dim>
533  void DoFs_summary<dim>::init_set_dofs(SymTensor2<dim> &eps, SymmetricTensor<2,dim> &eps_init, SW_double2<dim> &double_arg, double &double_init )
534 
535  {
536  const unsigned int nbr_total_dofs = eps.n_dofs + double_arg.n_dofs;
537 
538  eps.start_index = 0;
539  double_arg.start_index = eps.n_dofs;
540 
541  eps.init_set_dofs( eps_init, nbr_total_dofs);
542  double_arg.init_set_dofs( double_init, nbr_total_dofs );
543  }
544 
545  template<int dim>
546  void DoFs_summary<dim>::set_dofs(SymTensor<dim> &eps, SW_double<dim> &double_arg1, SW_double<dim> &double_arg2, SW_double<dim> &double_arg3 )
547  {
548  const unsigned int nbr_total_dofs = eps.n_dofs + double_arg1.n_dofs + double_arg2.n_dofs + double_arg3.n_dofs ;
549 
550  eps.start_index = 0;
551  double_arg1.start_index = eps.n_dofs;
552  double_arg2.start_index = eps.n_dofs + double_arg1.n_dofs;
553  double_arg3.start_index = eps.n_dofs + double_arg1.n_dofs + double_arg2.n_dofs;
554 
555  eps.set_dofs( nbr_total_dofs );
556  double_arg1.set_dofs( nbr_total_dofs );
557  double_arg2.set_dofs( nbr_total_dofs );
558  double_arg3.set_dofs( nbr_total_dofs );
559  }
560 
561  template<int dim>
563  {
564  const unsigned int nbr_total_dofs = eps.n_dofs + double_arg1.n_dofs + double_arg2.n_dofs ;
565 
566  eps.start_index = 0;
567  double_arg1.start_index = eps.n_dofs;
568  double_arg2.start_index = eps.n_dofs + double_arg1.n_dofs;
569 
570  eps.set_dofs( nbr_total_dofs );
571  double_arg1.set_dofs( nbr_total_dofs );
572  double_arg2.set_dofs( nbr_total_dofs );
573  }
574 
575  /*
576  *
577  * example call: DoFs_summary.get_curvature(d2_energy_d_eps_d_phi, energy, eps_fad, phi_fad)
578  */
579  template<int dim>
580  void DoFs_summary<dim>::get_curvature( SymmetricTensor<2,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument, SymTensor2<dim> &eps, SW_double2<dim> &double_arg )
581  {
582  SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > d_arg_d_eps;
583  eps.get_tangent(d_arg_d_eps, argument);
584  double_arg.get_curvature(Curvature, d_arg_d_eps, eps);
585  }
586 
587 
588  /*
589  *
590  * example call: DoFs_summary.get_curvature(d2_energy_d_phi_d_eps, energy, phi_fad, eps_fad)
591  */
592  template<int dim>
593  void DoFs_summary<dim>::get_curvature( SymmetricTensor<2,dim> &Curvature, Sacado::Fad::DFad<DFadType> &argument, SW_double2<dim> &double_arg, SymTensor2<dim> &eps )
594  {
595  double_arg.get_curvature(Curvature, argument, eps);
596  }
597 }
598 
599 
600 
601 
602 
603 
604 #endif // Sacado_Wrapper_H
605 
Definition: Sacado_Wrapper.h:66
std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
Definition: Sacado_Wrapper.h:202
void get_index_map(std::map< unsigned int, std::pair< unsigned int, unsigned int >> &std_map_indicies)
Definition: Sacado_Wrapper.h:24
unsigned int start_index
Definition: Sacado_Wrapper.h:337
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:78
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:205
Sacado::Fad::DFad< double > DFadType
Definition: Sacado_Wrapper.h:19
SW_double & operator=(fad_double fad_assignment)
Definition: Sacado_Wrapper.h:334
SW_double2 & operator=(double double_init)
Definition: Sacado_Wrapper.h:407
void set_dofs(unsigned int nbr_total_dofs=n_dofs)
Definition: Sacado_Wrapper.h:357
void init_set_dofs(SymmetricTensor< 2, dim > &tensor_double, const unsigned int nbr_total_dofs=n_dofs)
Definition: Sacado_Wrapper.h:225
Definition: Sacado_Wrapper.h:395
void get_tangent(SymmetricTensor< 4, dim > &Tangent, SymmetricTensor< 2, dim, Sacado::Fad::DFad< DFadType > > &sigma)
Definition: Sacado_Wrapper.h:252
void get_curvature(double &Curvature, Sacado::Fad::DFad< DFadType > &argument)
Definition: Sacado_Wrapper.h:456
void init_set_dofs(const double &double_init, unsigned int nbr_total_dofs=n_dofs)
Definition: Sacado_Wrapper.h:428
Definition: Sacado_Wrapper.h:505
unsigned int start_index
Definition: Sacado_Wrapper.h:77
SW_double & operator=(double double_init)
Definition: Sacado_Wrapper.h:333
std::map< unsigned int, std::pair< unsigned int, unsigned int > > std_map_indicies
Definition: Sacado_Wrapper.h:75
SymTensor()
Definition: Sacado_Wrapper.h:69
Sacado::Fad::DFad< double > fad_double
Definition: Sacado_Wrapper.h:18
unsigned int start_index
Definition: Sacado_Wrapper.h:204
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:410
SymTensor2()
Definition: Sacado_Wrapper.h:196
Sacado::Fad::DFad< double > fad_double
Definition: Sacado-auxiliary_functions.h:7
Definition: Sacado_Wrapper.h:321
Definition: Sacado_Wrapper.h:63
unsigned int start_index
Definition: Sacado_Wrapper.h:411
Definition: Sacado_Wrapper.h:193
void set_dofs(unsigned int nbr_total_dofs=n_dofs)
Definition: Sacado_Wrapper.h:110
void get_values(SymmetricTensor< 2, dim > &tensor_double)
Definition: Sacado_Wrapper.h:173
static const unsigned int n_dofs
Definition: Sacado_Wrapper.h:336
SW_double2 & operator=(Sacado::Fad::DFad< DFadType > fad_assignment)
Definition: Sacado_Wrapper.h:408