Trilinos::Sacado example
mainpage.h
Go to the documentation of this file.
1 /**
2 
3 @{
4 
5 \brief Example usage of Sacado for tensor calculus
6 
7 \author jfriedlein
8 
9 @tableofcontents
10 
11 @mainpage Trilinos::Sacado example documentation
12 
13 @section intro Introduction
14 The way we see and use Sacado here is as follows: \n
15 If you usually compute the following equation \n
16 \f[ c = a + b \f]
17 for instance with data types double as \n
18 \f[ 1.0 + 1.0 \rightarrow 2.0 \f]
19 your results is just a double number \f$ c \f$ that contains the value \f$ 2 \f$. \n
20 
21 Using Sacado, on the other hand, the variable \f$ c_{fad} \f$ is now of, for example, data type
22 @code
23  Sacado::Fad::DFad<double> c_fad;
24 @endcode
25 As a result, \f$ c_{fad} \f$ now contains not just the number \f$ 2 \f$, but also all the derivatives of \f$ c_{fad} \f$
26 with respect to the previously defined degrees of freedom (set via command *.diff(*)). \n
27 The following figure tries to visualize this for first derivatives:
28 \image html Sacado_data-type.png
29 and for second derivatives (requires different data type):
30 \image html Sacado_data-type_1_2_derivatives.png
31 Every variable that was declared as such a data type now contains besides the actual value \f$ c \f$, also its derivatives.
32 
33 @subsection subsec_overview Overview
34 @todo Add the example 10 (templated functions) to the overview and some text on why to use it
35 @todo add the example 9 (factor 0.5) to the background group
36 
37 This overview shall give you a first impression what to expect from each of the examples.
38 The background/basics group gives you the promised look under the hood.
39 Whereas the application group shows you how you can use the Sacado_Wrapper to quickly compute tangents.
40 \image html overview_of_examples.png
41 
42 If you right away want to use Sacado, then you might skip the first examples and jump to \ref Ex3B "example 3B".
43 There we show how to use the "Sacado_Wrapper" that does everything from \ref Ex2 "example 2" and \ref Ex3 "example 3" in just a view lines of code. This does not mean that
44 the here shown approach is the fastest or most efficient, it is just simple and easy to use.
45 
46 Furthermore, if you, for instance, compute problems with two-fields (e.g. displacement and scalar damage) and you need
47 tangents with respect to both a tensor (e.g. strain tensor) and a scalar (e.g. damage variable), you can use the Sacado_Wrapper as shown in \ref Ex4 "example 4".
48 
49 @subsection subsec_more_basics Some more basics
50 One can access the double value of \f$ c_{fad} \f$ with the Sacado command *.val():
51 @code
52  double c_value = c_fad.val();
53 @endcode
54 The derivatives of \f$ c_{fad} \f$ can be accessed with the command *.dx():
55 @code
56  double d_c_d_a = c_fad.dx(0);
57  double c_c_d_b = c_fad.dx(1);
58 @endcode
59 The arguments of \a dx, namely 0 and 1 are the numbers corresponding to the dof that belong to \a a and \a b. More details on how to set this up
60 and use it, are given in \ref Ex1 "example 1".
61 
62 @subsection subsec_resources Some resources/links
63 You can use Sacado to compute general derivatives of functions (with or without tensors) with respect to variables (double, Tensors, ...).
64 @todo Link is broken
65 
66 - Basics on Sacado, automatic differentation and coding examples (important remarks, e.g. on point of non-differentiability): \n
67  https://software.sandia.gov/SESS/past_seminars/111307_Phipps.pdf
68 - Basics and available autodiff libraries provided by deal.ii: \n
69  https://www.dealii.org/current/doxygen/deal.II/group__auto__symb__diff.html
70 - Template-based generic programming: \n
71  downloads.hindawi.com/journals/sp/2012/202071.pdf
72 - Tutorial 33 by deal.ii: introduction to Sacado and implementation to assemble the residuum and compute its derivative \n
73  https://www.dealii.org/current/doxygen/deal.II/step_33.html
74 
75 The here shown examples shall solely show how Sacado can be applied and give some background and a look under the hood.
76 The code is neither elegant nor efficient, but it works. A more user-friendly version is provided by means of the "Sacado_Wrapper". \n
77 @todo Check if factor 0.5 is also necessary for d_sigma / d_phi
78 
79 @note This documentation and code only protocol my first steps with Sacado. They are not guaranteed to be correct neither are they verified.
80 Any comments, criticism, corrections, feedback, improvements, ... are appreciated and very well welcomed.
81 
82 @section code The commented program
83 
84 \code
85 
86 /*
87  * Author: jfriedlein, 2019
88  * dsoldner, 2019
89  */
90 
91 \endcode
92 @section includes Include Files
93 The data type SymmetricTensor and some related operations, such as trace, symmetrize, deviator, ... for tensor calculus
94 \code
95 #include <deal.II/base/symmetric_tensor.h>
96 
97 \endcode
98 C++ headers (some basics, standard stuff)
99 \code
100 #include <iostream>
101 #include <fstream>
102 #include <cmath>
103 
104 \endcode
105 TimerOutput: to measure the time spent in certain parts of the code
106 \code
107 #include <deal.II/base/timer.h>
108 
109 \endcode
110 Sacado (from Trilinos, data types, operations, ...)
111 \code
112 #include <Sacado.hpp>
113 
114 #include "Sacado_Wrapper.h"
116 
117 
118 
119 \endcode
120 Those headers are related to data types and autodiff, but don't seem to be needed (nor make a difference)
121 \code
122 //# include <deal.II/base/numbers.h>
123 //# include <deal.II/differentiation/ad/ad_number_traits.h>
124 //# include <deal.II/differentiation/ad/sacado_number_types.h>
125 
126 \endcode
127 According to the basics of deal.ii-programming (see dealii.org and https://www.dealii.org/current/doxygen/deal.II/step_1.html for a start)
128 \code
129 using namespace dealii;
130 
131 \endcode
132 Defining a data type for the Sacado variables (here we simply used the standard types from the deal.ii step-33 tutorial's introduction)
133 \code
134 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(*))
135 
136 
137 \endcode
138 @section Ex1 1. example: simple scalar equation
139 1. example: simple scalar equation from deal.ii-tutorial step-33 (see the introduction there to get a first impression, https://www.dealii.org/current/doxygen/deal.II/step_33.html)
140 @todo clean up the documentation of the classes
141 \code
142 //
143 void sacado_test_scalar ()
144 {
145  std::cout << "Scalar Test:" << std::endl;
146 \endcode
147 Define the variables used in the computation (inputs/independent variables: a, b; output/result: c; auxiliaries/passive variables: *) as the Sacado-data type
148 \code
149  fad_double a,b,c;
150 \endcode
151 Initialize the input variables a and b; This (a,b) = (1,2) will be the point where the derivatives are computed.
152 Compare: y=x² -> (dy/dx)(\@x=1) = 2. We can only compute the derivative numerically at a certain point.
153 \code
154  a = 1;
155  b = 2;
156 
157  a.diff(0,2); // Set a to be dof 0, in a 2-dof system.
158  b.diff(1,2); // Set b to be dof 1, in a 2-dof system.
159 \endcode
160 Our equation here is very simply. But you can use nested equations and many standard mathematical operations, such as sqrt, pow, sin, ...
161 \code
162  c = 2*a + std::cos(a*b);
163  double *derivs = &c.fastAccessDx(0); // Access the derivatives of
164 \endcode
165 Output the derivatives of c with respect to the two above defined degrees of freedom (dof)
166 \code
167  std::cout << "Derivatives at the point (" << a << "," << b << ")" << std::endl;
168  std::cout << "dc/da = " << derivs[0] << ", dc/db=" << derivs[1] << std::endl;
169 }
170 
171 
172 \endcode
173 @section Ex2 2. example: Preparation for the use of Sacado with tensors
174 Here we want to introduce tensors for the first time. Hence, we limit ourselves to a trivial equation relating the strain tensor \a eps
175 with dim x dim components with the stress tensor \a sigma. Both here used tensors are symmetric, hence we use the SymmetricTensor class and
176 have to keep some details in mind (see below factor 0.5 related to Voigt-Notation). Don't be scared by the enormous number of repetitive
177 lines of code, everything shown in this example and the following will be handled by the Sacado_Wrapper with roughly four lines of code.
178 \code
179 /*
180  * 2. example: use of tensors
181  */
182 void sacado_test_2 ()
183 {
184  std::cout << "Test 2:" << std::endl;
185 
186 \endcode
187 First we set the dimension \a dim: 2D->dim=2; 3D->dim=3 \n This defines the "size" of the tensors and the number of dofs. \ref Ex2 "Example 2" only works in 3D, whereas the following Ex3 is set up dimension-independent.
188 \code
189  const unsigned int dim = 3;
190 
191 \endcode
192 Declare our input, auxiliary and output variables as SymmetricTensors consisting of fad_doubles (instead of the standard SymmetricTensor out of doubles)
193 \code
194  SymmetricTensor<2,dim, fad_double> sigma, eps;
195 
196 \endcode
197 Init the strain tensor (the point at which the derivative shall be computed)
198 \code
199  eps[0][0] = 1;
200  eps[1][1] = 2;
201  eps[2][2] = 3;
202  eps[0][1] = 4;
203  eps[0][2] = 5;
204  eps[1][2] = 6;
205 
206 \endcode
207 Now we declare the dofs. The derivative to a tensor requires all components, therefore we set the components of the strain tensor here one by one as the dofs.
208 Because our tensors are symmetric, we only need 6 components in 3D instead of 9 for a full second order tensor
209 \code
210  eps[0][0].diff(0,6);
211  eps[1][1].diff(1,6);
212  eps[2][2].diff(2,6);
213  eps[0][1].diff(3,6);
214  eps[0][2].diff(4,6);
215  eps[1][2].diff(5,6);
216 
217 \endcode
218 The equation describing the stresses (here just a simple test case)
219 \code
220  sigma = eps;
221 
222 \endcode
223 Let's output the computed stress tensor.
224 \code
225  std::cout << sigma << std::endl;
226 \endcode
227 The resulting values of \a sigma are fairly boring, due to our simple equation. It is the additional output generated by
228 this, that is interesting here: \n
229 output: \n
230 1 [ 1 0 0 0 0 0 ] 4 [ 0 0 0 1 0 0 ] 5 [ 0 0 0 0 1 0 ] 4 [ 0 0 0 1 0 0 ] 2 [ 0 1 0 0 0 0 ] 6 [ 0 0 0 0 0 1 ] 5 [ 0 0 0 0 1 0 ] 6 [ 0 0 0 0 0 1 ] 3 [ 0 0 1 0 0 0 ] \n
231 The numbers 1, 4, 5, 4, ... are the entries in the stress tensor \a sigma. In square brackets we see the derivatives of sigma with respect to all the dofs set previously
232 given in the order we defined them above. Meaning: The first entry in the square brackets corresponds to the 0-th dof set by
233 @code eps[0][0].diff(0,6); @endcode referring to the component (0,0) in the strain tensor \a eps.
234 
235 Computing the derivatives for certain components of the resulting tangent modulus: \n
236 We now access these lists of derivatives (output above in square brackets) for one component of the stress tensor \a sigma at a time.
237 \code
238  {
239 \endcode
240 Access the derivatives corresponding to the component (0,0) of the stress tensor \a sigma
241 \code
242  double *derivs = &sigma[0][0].fastAccessDx(0);
243 \endcode
244 The following output will show us the same derivatives that we already saw above, just formatted differently \n
245 output: d_sigma[0][0]/d_eps = 1 , 0 , 0 , 0 , 0 , 0 ,
246 \code
247  std::cout << "d_sigma[0][0]/d_eps = ";
248  for ( unsigned int i=0; i<6; ++i)
249  std::cout << derivs[i] << " , ";
250  std::cout << std::endl;
251  }
252  {
253 \endcode
254 Access the derivatives corresponding to the component (1,2) of the stress tensor \a sigma
255 \code
256  double *derivs = &sigma[1][2].fastAccessDx(0);
257 \endcode
258 output: d_sigma[1][2]/d_eps = 0 , 0 , 0 , 0 , 0 , 1 ,
259 \code
260  std::cout << "d_sigma[1][2]/d_eps = ";
261  for ( unsigned int i=0; i<6; ++i)
262  std::cout << derivs[i] << " , ";
263  std::cout << std::endl;
264  }
265 }
266 
267 
268 \endcode
269 @section Ex3 3. example: Using a slightly more complicated stress equation
270 \code
271 void sacado_test_3 ()
272 {
273  std::cout << "Test 3:" << std::endl;
274 
275  const unsigned int dim = 3;
276 
277 \endcode
278 Here we also define some constant, for instance the bulk modulus \a kappa and the second Lamè parameter \a mu.
279 We now also define one of our constants as fad_double. By doing this we can use the normal multiplication (see below).
280 \code
281  double kappa_param = 5;
282  fad_double kappa (kappa_param);
283 \endcode
284 The second constant remains as a double just to show the difference.
285 \code
286  double mu = 2;
287 
288  SymmetricTensor<2,dim, fad_double> sigma, eps;
289 
290 \endcode
291 To simplify the access to the dofs we define a map that relate the components of our strain tensor to the dof-nbr
292 \code
293  std::map<unsigned int,std::pair<unsigned int,unsigned int>> std_map_indicies;
294 
295 \endcode
296 The point at which the derivative shall be computed: \n
297 As mentioned previously, we will implement this example for 2D and 3D, hence we once have to set up a strain tensor
298 and the derivatives for 3D with 6 independent components ...
299 \code
300  if(dim==3)
301  {
302  eps[0][0] = 1;
303  eps[1][1] = 2;
304  eps[2][2] = 3;
305 
306  eps[0][1] = 4;
307  eps[0][2] = 5;
308  eps[1][2] = 6;
309 
310 
311  eps[0][0].diff(0,6);
312  eps[0][1].diff(1,6);
313  eps[0][2].diff(2,6);
314  eps[1][1].diff(3,6);
315  eps[1][2].diff(4,6);
316  eps[2][2].diff(5,6);
317 
318 \endcode
319 By using the map and the following pairs, we have to set up the relation between strain components and dofs only once
320 and can use the map to access the entries of the list later, without possibly mixing up indices and creating errors.
321 Please don't be confused, but the dofs in the Wrapper are set up
322 in a different order that we showed earlier. Earlier: (0,0)-(1,1)-(2,2)-...; Now: (0,0)-(0,1)-(0,2)-...
323 \code
324  std::pair<unsigned int, unsigned int> tmp_pair;
325  tmp_pair.first=0; tmp_pair.second=0;
326  std_map_indicies[0] = tmp_pair;
327 
328  tmp_pair.first=0; tmp_pair.second=1;
329  std_map_indicies[1] = tmp_pair;
330 
331  tmp_pair.first=0; tmp_pair.second=2;
332  std_map_indicies[2] = tmp_pair;
333 
334  tmp_pair.first=1; tmp_pair.second=1;
335  std_map_indicies[3] = tmp_pair;
336 
337  tmp_pair.first=1; tmp_pair.second=2;
338  std_map_indicies[4] = tmp_pair;
339 
340  tmp_pair.first=2; tmp_pair.second=2;
341  std_map_indicies[5] = tmp_pair;
342  }
343 \endcode
344 ... and once for 2D with just 3 independent components.
345 \code
346  else if(dim==2)
347  {
348  eps[0][0] = 1;
349  eps[1][1] = 2;
350 
351  eps[0][1] = 4;
352 
353 
354  eps[0][0].diff(0,3);
355  eps[0][1].diff(1,3);
356  eps[1][1].diff(2,3);
357 
358  std::pair<unsigned int, unsigned int> tmp_pair;
359  tmp_pair.first=0; tmp_pair.second=0;
360  std_map_indicies[0] = tmp_pair;
361 
362  tmp_pair.first=0; tmp_pair.second=1;
363  std_map_indicies[1] = tmp_pair;
364 
365  tmp_pair.first=1; tmp_pair.second=1;
366  std_map_indicies[2] = tmp_pair;
367  }
368  else
369  {
370  throw std::runtime_error("only dim==2 or dim==3 allowed");
371  }
372 
373 \endcode
374 Instead of calling the *.diff(*) on the components one-by-one we could also use the following for-loop, so
375 we also use the map to set the dofs (as we will do in the Wrapper later).
376 @code
377 for ( unsigned int x=0; x<((dim==2)?3:6); ++x )
378 {
379  unsigned int i=std_map_indicies[x].first;
380  unsigned int j=std_map_indicies[x].second;
381  eps[i][j].diff(x,((dim==2)?3:6));
382 }
383 @endcode
384 
385 For our slightly more complicated stress equation we need the unit and deviatoric tensors.
386 We can simply define them by writing the values of the already existing deal.ii functions into newly
387 defined SymmetricTensors build from fad_doubles.
388 \code
389  SymmetricTensor<2,dim, fad_double> stdTensor_I (( unit_symmetric_tensor<dim,fad_double>()) );
390  SymmetricTensor<4,dim, fad_double> stdTensor_Idev ( (deviator_tensor<dim,fad_double>()) );
391 
392 \endcode
393 With everything set and defined, we can compute our stress \a sigma according to:
394 \f[ \sigma = \kappa \cdot trace(\varepsilon) \cdot \boldsymbol{I} + 2 \cdot \mu \cdot \varepsilon^{dev} \f]
395 Here you can see that we can directly multiply the constant and the tensors when kappa is also declared as fad_double
396 \code
397  sigma = kappa * (trace(eps) * stdTensor_I);
398 \endcode
399 We didn't do the same for mu to once again emphasize the difference between constants as double and as fad_double. \n
400 The remaining code uses a normal double constant.
401 \code
402  SymmetricTensor<2,dim,fad_double> tmp = deviator<dim,fad_double>(symmetrize<dim,fad_double>(eps)); tmp*=(mu*2);
403  sigma += tmp;
404 \endcode
405 The fairly cumbersome computation is caused by the way the operators are set up for tensors out of fad_doubles.
406 
407 \code
408  std::cout << "sigma=" << sigma << std::endl;
409 
410 \endcode
411 Now we want to actually build our tangent modulus called \a C_Sacado that contains all the derivatives and relates
412 the stress tensor with the strain tensor. \n
413 The fourth-order tensor \a C_Sacado is our final goal, we don't have to compute anything that is related to Sacado with
414 this tensor, so we can finally return to our standard SymmetricTensor out of doubles. The latter is necessary to use
415 the tangent in the actual FE code.
416 \code
417  SymmetricTensor<4,dim> C_Sacado;
418 
419 \endcode
420 As in \ref Ex2 "example 2" we access the components of the stress tensor one by one. In order to capture all of them we loop over the
421 components i and j of the stress tensor.
422 \code
423  for ( unsigned int i=0; i<dim; ++i)
424  for ( unsigned int j=0; j<dim; ++j )
425  {
426  double *derivs = &sigma[i][j].fastAccessDx(0); // Access the derivatives of the (i,j)-th component of \a sigma
427 
428 \endcode
429 To visually ensure that every stress component has in fact all 6 derivatives for 3D or 3 for 2D, we output the size:
430 \code
431  std::cout<<"size: "<<sigma[i][j].size()<<std::endl;
432 
433 \endcode
434 We loop over all the dofs. To be able to use this independent of the chosen dimension \a dim, we use a ternary operator
435 to decide whether we have to loop over 6 derivatives or just 3.
436 \code
437  for(unsigned int x=0;x<((dim==2)?3:6);++x)
438  {
439  unsigned int k=std_map_indicies[x].first;
440  unsigned int l=std_map_indicies[x].second;
441 
442  if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
443  {
444  C_Sacado[i][j][k][l] = 0.5*derivs[x];
445  C_Sacado[i][j][l][k] = 0.5*derivs[x];
446  }
447  else
448  C_Sacado[i][j][k][l] = derivs[x];
449  }
450 
451  }
452 
453 \endcode
454 After resembling the fourth-order tensor, we now have got our tangent saved in \a C_Sacado ready to be used
455 
456 To ensure that Sacado works properly, we can compute the analytical tangent for comparison
457 \code
458  double kappa_d = 5;
459  double mu_d = 2;
460 \endcode
461 Our stress equation in this example is still simple enough to derive the tangent analytically by hand:
462 \f[ \overset{4}{C_{analy}} = \kappa \cdot \boldsymbol{I} \otimes \boldsymbol{I} + 2 \cdot \mu \cdot \overset{4}{I^{dev}} \f]
463 \code
464  SymmetricTensor<4,dim> C_analy = kappa_d * outer_product(unit_symmetric_tensor<dim>(), unit_symmetric_tensor<dim>()) + 2* mu_d * deviator_tensor<dim>();
465 
466 
467 \endcode
468 We again define our strain tensor \a eps_d (*_d for standard double in contrast to fad_double)
469 \code
470  SymmetricTensor<2,dim> eps_d;
471 
472  if(dim==3)
473  {
474  eps_d[0][0] = 1;
475  eps_d[1][1] = 2;
476  eps_d[2][2] = 3;
477 
478  eps_d[0][1] = 4;
479  eps_d[0][2] = 5;
480  eps_d[2][1] = 6;
481 
482  }
483  else if(dim==2)
484  {
485  eps_d[0][0] = 1;
486  eps_d[1][1] = 2;
487 
488  eps_d[1][0] = 4;
489 
490  }
491  else
492  {
493  throw std::runtime_error("only dim==2 or dim==3 allowed");
494  }
495 \endcode
496 @todo use boldsymbol for tensors
497 
498 To output the stress tensor we first have to compute it. We do this here via
499 \f[ \sigma = \overset{4}{C_{analy}} : \varepsilon \f]
500 The output exactly matched the result obtained with Sacado.
501 @note Checking the Sacado stress tensor against an analytically computed or otherwise determined stress tensor is absolutely no way to check whether
502 the tangent computed via Sacado is correct. When we compute the stress tensor with Sacado and for example mix up a + and - sign, this might not matter
503 at all if the number that is added or subtracted is small. However, for the tangent this nasty sign can be very critical. Just keep in mind: the
504 tangent has 81 components and the stress tensor just 9, so how does one want to verify 81 variables by comparing 9?
505 
506 \code
507  std::cout << "sigma_analy: " << (C_analy*eps_d) << std::endl;
508 
509 \endcode
510 That's the reason we compare all the entries in the Sacado and the analytical tensor one by one
511 \code
512  for (unsigned int i=0; i<dim; ++i)
513  for ( unsigned int j=0; j<dim; ++j)
514  for ( unsigned int k=0; k<dim; ++k)
515  for ( unsigned int l=0; l<dim; ++l)
516  std::cout << "C_analy["<<i<<"]["<<j<<"]["<<k<<"]["<<l<<"] = " << C_analy[i][j][k][l] << " vs C_Sacado: " << C_Sacado[i][j][k][l] << std::endl;
517 
518 
519 \endcode
520 To simplify the comparison we compute a scalar error as the sum of the absolute differences of each component
521 \code
522  double error_Sacado_vs_analy=0;
523  for (unsigned int i=0; i<dim; ++i)
524  for ( unsigned int j=0; j<dim; ++j)
525  for ( unsigned int k=0; k<dim; ++k)
526  for ( unsigned int l=0; l<dim; ++l)
527  error_Sacado_vs_analy += std::fabs(C_Sacado[i][j][k][l] - C_analy[i][j][k][l]);
528 
529 
530 \endcode
531 As desired: The numerical error is zero (0 in double precision) and the tensor components are equal
532 \code
533  std::cout << "numerical error: " << error_Sacado_vs_analy << std::endl;
534 }
535 
536 
537 \endcode
538 @section Ex3B 3B. Example: Using the wrapper for Ex3
539 \code
540 void sacado_test_3B ()
541 {
542  std::cout << "Test 3B:" << std::endl;
543  const unsigned int dim=3;
544 
545 \endcode
546 The following declarations are usually input arguments. So you receive the strain tensor and the constants out of doubles.
547 \code
548  SymmetricTensor<2,dim> eps_d;
549  eps_d[0][0] = 1;
550  eps_d[1][1] = 2;
551  eps_d[2][2] = 3;
552 
553  eps_d[0][1] = 4;
554  eps_d[0][2] = 5;
555  eps_d[1][2] = 6;
556 
557  double kappa = 5;
558  double mu = 2;
559 
560 \endcode
561 Now we start working with Sacado: \n
562 When we use the index notation to compute e.g. our stress we do not need to declare our constants (here kappa, mu) as
563 fad_double.
564 
565 We declare our strain tensor as the special data type Sacado_Wrapper::SymTensor from the file "Sacado_Wrapper.h"
566 where this data type was derived from the SymmetricTensor<2,dim,fad_double>.
567 \code
568  Sacado_Wrapper::SymTensor<dim> eps;
569 
570 \endcode
571 Next we initialize our Sacado strain tensor with the values of the inputed double strain tensor:
572 \code
573  eps.init(eps_d);
574 
575 \endcode
576 We define all the entries in the symmetric tensor \a eps as the dofs. So we can later derive any variable
577 with respect to the strain tensor \a eps.
578 @note Keep the order first declare \a eps, then initialise it and afterwards set it as dofs
579 
580 \code
581  eps.set_dofs();
582 
583 \endcode
584 Now we declare our output and auxiliary variables as Sacado-Tensors.
585 \code
586  SymmetricTensor<2,dim,fad_double> sigma;
587 
588  SymmetricTensor<2,dim, fad_double> stdTensor_I (( unit_symmetric_tensor<dim,fad_double>()) );
589 
590 \endcode
591 Our stress equation is now computed in index notation to simplify the use of the constants and
592 especially the use of the \a deviator.
593 \code
594  for ( unsigned int i=0; i<dim; ++i)
595  for ( unsigned int j=0; j<dim; ++j )
596  sigma[i][j] = kappa * trace(eps) * stdTensor_I[i][j] + 2. * mu * deviator(eps)[i][j];
597 
598 \endcode
599 Finally we declare our desired tangent as the fourth order tensor \a C_Sacado and compute the tangent via
600 the command \a get_tangent.
601 \code
602  SymmetricTensor<4,dim> C_Sacado;
603  eps.get_tangent(C_Sacado, sigma);
604 
605 \endcode
606 We could again compare the herein computed tangent with the analytical tangent from Ex2, but as before
607 the results are fairly boring, because Sacado hits the analytical tangent exactly --- no surprise for such
608 simple equations.
609 
610 And that's it. By using the Sacado_wrapper we can achieve everything from Ex2 (besides the equations)
611 with just four lines of code namely:
612 - eps.init(eps_d); // To initialize the Sacado strain tensor
613 - eps.set_dofs(); // To declare the components of eps as the dofs
614 - eps.get_tangent(*); // To get the tangent
615 \code
616 }
617 
618 
619 \endcode
620 @section Ex4 4. Example: Computing derivatives with respect to a tensor and a scalar
621 \code
622 void sacado_test_4 ()
623 {
624  std::cout << "Test 4:" << std::endl;
625  const unsigned int dim=3;
626 
627 \endcode
628 The following declarations are usually input arguments. So you receive the strain tensor \q eps_d,
629 the damage variable \a phi and the constants \a kappa and \a mu out of doubles.
630 \code
631  SymmetricTensor<2,dim> eps_d;
632  eps_d[0][0] = 1;
633  eps_d[1][1] = 2;
634  eps_d[2][2] = 3;
635 
636  eps_d[0][1] = 4;
637  eps_d[0][2] = 5;
638  eps_d[1][2] = 6;
639 
640  double phi_d = 0.3;
641 
642 \endcode
643 We don't need these constants in the current example.
644 double kappa = 5;
645 double mu = 2;
646 
647 
648 We set up our strain tensor as in Ex3B.
649 \code
650  Sacado_Wrapper::SymTensor<dim> eps;
651  Sacado_Wrapper::SW_double<dim> phi;
652 
653 \endcode
654 Initialize the strain tensor and the damage variable
655 \code
656  eps.init(eps_d);
657  phi.init(phi_d);
658 
659 \endcode
660 Set the dofs, where the argument sets the total nbr of dofs (3 or 6 for the sym. tensor and 1 for the double)
661 \code
662 // eps.set_dofs(eps.n_independent_components+1/*an additional dof for phi*/);
663 //
664 \endcode
665 In order to also compute derivatives with respect to the scalar \a phi, we add this scalar to our list
666 of derivatives. Because we have already defined 3 or 6 dofs our additional dof will be placed at the end
667 of this list. We set this up with the member variable start_index ...
668 \code
669 // phi.start_index=eps.n_independent_components;
670 \endcode
671 and again using the input argument representing the total number of dofs
672 \code
673 // phi.set_dofs(eps.n_independent_components+1);
674 
675 \endcode
676 All of the above 3 lines of code are automatically done by the DoFs_summary class. So, to
677 set our dofs we just create an instance and call set_dofs with our variables containing the desired dofs.
678 \code
679  Sacado_Wrapper::DoFs_summary<dim> DoFs_summary;
680  DoFs_summary.set_dofs(eps, phi);
681 
682 
683 \endcode
684 Compute the stress tensor and damage variable \a d (here we just use some arbitrary equations for testing): \n
685 Let us first declare our output (and auxiliary) variables as Sacado data types.
686 \code
687  SymmetricTensor<2,dim,fad_double> sigma;
688  fad_double d;
689 \endcode
690 @todo It would be nice to use the data types from the Sacado_Wrapper for all the Sacado variables. But
691 somehow the operators (multiply*, ...) seem to cause conflicts again.
692 
693 The actual computation in the following scope uses the exact same equation as your normal computation e. g. via the data type double.
694 Hence, you could either directly compute your stress, etc. via the Sacado variables or you define
695 template functions that contain your equations and are either called templated with double or fad_double.
696 When using the first option, please consider the computation time that is generally higher for a computation
697 with fad_double than with normal doubles (own experience in a special case: slower by factor 30).
698 The second option with templates does not suffer these issues.
699 \code
700  {
701  d = phi*phi + 25 + trace(eps) + eps.norm();
702  std::cout << "d=" << d << std::endl;
703 
704  for ( unsigned int i=0; i<dim; ++i)
705  for ( unsigned int j=0; j<dim; ++j )
706  sigma[i][j] = phi * d * eps[i][j];
707 \endcode
708 ToDo: strangely when phi is a fad_double then the multiplication phi * eps works directly without
709 having to use the index notation
710 \code
711  std::cout << "sigma=" << sigma << std::endl << std::endl;
712  }
713 
714 
715 \endcode
716 Get the tangents
717 d_sigma / d_eps: SymmetricTensor with respect to SymmetricTensor
718 \code
719  SymmetricTensor<4,dim> C_Sacado;
720  eps.get_tangent(C_Sacado, sigma);
721  std::cout << "C_Sacado=" << C_Sacado << std::endl;
722 
723 \endcode
724 Compute the analytical tangent:
725 \code
726  SymmetricTensor<4,dim> C_analy;
727  C_analy = ( std::pow(phi_d, 3) + 25*phi_d + phi_d*trace(eps_d) + phi_d*eps_d.norm() ) * identity_tensor<dim>()
728  + phi_d * outer_product( eps_d, unit_symmetric_tensor<dim>())
729  + phi_d * outer_product( eps_d, eps_d ) * 1./eps_d.norm();
730 \endcode
731 @note Be aware of the difference between \f[ eps_d \otimes \boldsymbol{1} \text{ and } \boldsymbol{1} \otimes eps_d \f]
732 
733 \code
734  std::cout << "C_analy =" << C_analy << std::endl;
735 
736 \endcode
737 To simplify the comparison we compute a scalar error as the sum of the absolute differences of each component
738 \code
739  double error_Sacado_vs_analy=0;
740  for (unsigned int i=0; i<dim; ++i)
741  for ( unsigned int j=0; j<dim; ++j)
742  for ( unsigned int k=0; k<dim; ++k)
743  for ( unsigned int l=0; l<dim; ++l)
744  error_Sacado_vs_analy += std::fabs(C_Sacado[i][j][k][l] - C_analy[i][j][k][l]);
745  std::cout << "numerical error: " << error_Sacado_vs_analy << std::endl << std::endl;
746 
747 \endcode
748 d_d / d_eps: double with respect to SymmetricTensor
749 \code
750  SymmetricTensor<2,dim> d_d_d_eps;
751  eps.get_tangent(d_d_d_eps, d);
752  std::cout << "d_d_d_eps =" << d_d_d_eps << std::endl;
753  SymmetricTensor<2,dim> d_d_d_eps_analy;
754  d_d_d_eps_analy = unit_symmetric_tensor<dim>() + eps_d / eps_d.norm();
755  std::cout << "d_d_d_eps_analy=" << d_d_d_eps_analy << std::endl << std::endl;
756 
757 \endcode
758 d_sigma / d_phi: SymmetricTensor with respect to double
759 \code
760  SymmetricTensor<2,dim> d_sigma_d_phi;
761  phi.get_tangent(d_sigma_d_phi, sigma);
762  std::cout << "d_sigma_d_phi =" << d_sigma_d_phi << std::endl;
763  SymmetricTensor<2,dim> d_sigma_d_phi_analy;
764  d_sigma_d_phi_analy = ( phi_d*phi_d + 25 + trace(eps_d) + eps_d.norm() + 2 * phi_d*phi_d ) * eps_d;
765  std::cout << "d_sigma_d_phi_analy=" << d_sigma_d_phi_analy << std::endl << std::endl;
766 
767 \endcode
768 Retrieve the values stored in \a sigma:
769 \code
770  SymmetricTensor<2,dim> sigma_d;
771  for ( unsigned int i=0; i<dim; ++i)
772  for ( unsigned int j=0; j<dim; ++j )
773  sigma_d[i][j] = sigma[i][j].val();
774  std::cout << "sigma_d = " << sigma_d << std::endl;
775 
776 \endcode
777 d_d / d_phi: double with respect to double
778 \code
779  double d_d_d_phi;
780  phi.get_tangent(d_d_d_phi, d);
781  std::cout << "d_d_d_phi=" << d_d_d_phi << std::endl;
782 
783 \endcode
784 Retrieve the value stored in d
785 \code
786  double d_double = d.val();
787 
788 \endcode
789 Taylor-series for point x
790 @todo Maybe add some more text on the linearization
791 \code
792  double x = phi_d + 0.05;
793  std::cout << "d_lin = d_d_d_phi * (phi_d - x) + d@(phi_d) = " << d_d_d_phi * (x-phi_d) + d_double << std::endl;
794  std::cout << "d@x = " << x*x + 25 + trace(eps_d) + eps_d.norm() << std::endl;
795 
796 \endcode
797 And that's it. By using the Sacado_wrapper we can compute derivatives with respect to
798 a tensor and a scalar at the same time (besides the equations)
799 in essence with just the following lines of code namely:
800 - eps.init(eps_d); phi.init(phi_d); // To initialize the Sacado strain tensor and scalar damage variable
801 - DoFs_summary.set_dofs(eps, phi); // To declare the components of eps and phi as the dofs
802 - eps.get_tangent(*); // To get tangents with respect to eps
803 - phi.get_tangent(*); // To get tangents with respect to phi
804 \code
805 }
806 
807 
808 
809 \endcode
810 @section Ex5 5. Example: Using a vector-valued equation
811 \code
812 void sacado_test_5 ()
813 {
814  const unsigned int dim=3;
815  std::cout << "Test 5:" << std::endl;
816  Tensor<1,dim,fad_double> c;
817  fad_double a,b;
818  unsigned int n_dofs=2;
819  a = 1; b = 2; // at the point (a,b) = (1,2)
820  a.diff(0,2); // Set a to be dof 0, in a 2-dof system.
821  b.diff(1,2); // Set b to be dof 1, in a 2-dof system.
822 \endcode
823 c is now a vector with three components
824 \code
825  c[0] = 2*a+3*b;
826  c[1] = 4*a+5*b;
827  c[2] = 6*a+7*b;
828 
829 \endcode
830 Access to the derivatives works as before.
831 \code
832  for(unsigned int i=0;i<dim;++i)
833  {
834  const fad_double &derivs = c[i]; // Access derivatives
835  for(unsigned int j=0;j<n_dofs;++j)
836  {
837  std::cout << "Derivatives at the point (" << a << "," << b << ") for "
838  <<i<<"th component wrt "<<j<<"th direction "<< std::endl;
839  std::cout << "dc_i/dxj = " << derivs.fastAccessDx(j) << std::endl;
840  }
841  }
842 }
843 
844 
845 
846 \endcode
847 @section Ex6 6. Example: First and second derivatives - Scalar equation
848 The here shown example was copied from https://github.com/trilinos/Trilinos/blob/master/packages/sacado/example/dfad_dfad_example.cpp
849 and modified to get a first impression on how we can work with first and second derivatives
850 \code
851 void sacado_test_6 ()
852 {
853  std::cout << "Test 6:" << std::endl;
854 \endcode
855 Define the variables used in the computation (inputs: a, b; output: c; auxiliaries: *) as doubles
856 \code
857  double a=1;
858  double b=2;
859 
860 \endcode
861 Number of independent variables (scalar a and b)
862 \code
863  int num_dofs = 2;
864 
865 \endcode
866 Define another data type containing even more Sacado data types
867 @todo try to merge the fad_double data type with this templated data type
868 \code
869  typedef Sacado::Fad::DFad<double> DFadType;
870  Sacado::Fad::DFad<DFadType> afad(num_dofs, 0, a);
871  Sacado::Fad::DFad<DFadType> bfad(num_dofs, 1, b);
872  Sacado::Fad::DFad<DFadType> cfad;
873 
874 \endcode
875 Output the variables: We se that the values of \a a and \a b are set but the derivatives have not yet been fully declared
876 \code
877  std::cout << "afad=" << afad << std::endl;
878  std::cout << "bfad=" << bfad << std::endl;
879  std::cout << "cfad=" << cfad << std::endl;
880 
881 \endcode
882 Now we set the "inner" derivatives.
883 \code
884  afad.val() = fad_double(num_dofs, 0, a); // set afad.val() as the first dof and init it with the double a
885  bfad.val() = fad_double(num_dofs, 1, b);
886 
887 \endcode
888 Compute function and derivative with AD
889 \code
890  cfad = 2*afad + std::cos(afad*bfad);
891 
892 \endcode
893 After this, we output the variables again and see that some additional derivatives have been declared. Furthermore,
894 \a cfad is filled with the values and derivatives
895 \code
896  std::cout << "afad=" << afad << std::endl;
897  std::cout << "bfad=" << bfad << std::endl;
898  std::cout << "cfad=" << cfad << std::endl;
899 
900 \endcode
901 Extract value and derivatives
902 \code
903  double c_ad = cfad.val().val(); // r
904  double dcda_ad = cfad.dx(0).val(); // dr/da
905  double dcdb_ad = cfad.dx(1).val(); // dr/db
906  double d2cda2_ad = cfad.dx(0).dx(0); // d^2r/da^2
907  double d2cdadb_ad = cfad.dx(0).dx(1); // d^2r/dadb
908  double d2cdbda_ad = cfad.dx(1).dx(0); // d^2r/dbda
909  double d2cdb2_ad = cfad.dx(1).dx(1); // d^2/db^2
910 
911 \endcode
912 Now we can print the actual double value of c and some of the derivatives:
913 \code
914  std::cout << "c_ad=" << c_ad << std::endl;
915  std::cout << "Derivatives at the point (" << a << "," << b << ")" << std::endl;
916  std::cout << "dc/da = " << dcda_ad << ", dc/db=" << dcdb_ad << std::endl;
917  std::cout << "d²c/da² = " << d2cda2_ad << ", d²c/db²=" << d2cdb2_ad << std::endl;
918  std::cout << "d²c/dadb = " << d2cdadb_ad << ", d²c/dbda=" << d2cdbda_ad << std::endl;
919 }
920 
921 
922 \endcode
923 @section Ex7 7. Example: First and second derivatives - Using tensors (The full story)
924 \code
925 void sacado_test_7 ()
926 {
927  const unsigned int dim=3;
928 
929  std::cout << "Test 7:" << std::endl;
930 
931 \endcode
932 Defining the inputs (material parameters, strain tensor)
933 \code
934  double lambda=1;
935  double mu=2;
936  SymmetricTensor<2,dim, double> eps;
937 
938  eps[0][0] = 1.;
939  eps[1][1] = 2.;
940  eps[2][2] = 3.;
941 
942  eps[0][1] = 4.;
943  eps[0][2] = 5.;
944  eps[1][2] = 6.;
945 
946 \endcode
947 Here we skip the one-field example and right away show the equations for a two-field problem
948 with \a eps and \a phi.
949 \code
950  double phi=0.3;
951 
952 \endcode
953 Setup of the map relating the indices (as before)
954 \code
955  std::map<unsigned int,std::pair<unsigned int,unsigned int>> std_map_indicies;
956 
957  std::pair<unsigned int, unsigned int> tmp_pair;
958  tmp_pair.first=0; tmp_pair.second=0;
959  std_map_indicies[0] = tmp_pair;
960 
961  tmp_pair.first=0; tmp_pair.second=1;
962  std_map_indicies[1] = tmp_pair;
963 
964  tmp_pair.first=0; tmp_pair.second=2;
965  std_map_indicies[2] = tmp_pair;
966 
967  tmp_pair.first=1; tmp_pair.second=1;
968  std_map_indicies[3] = tmp_pair;
969 
970  tmp_pair.first=1; tmp_pair.second=2;
971  std_map_indicies[4] = tmp_pair;
972 
973  tmp_pair.first=2; tmp_pair.second=2;
974  std_map_indicies[5] = tmp_pair;
975 
976 \endcode
977 Number of independent variables (6 for the tensor and 1 for the scalar phi)
978 \code
979  const unsigned int nbr_dofs = 6+1;
980 
981 \endcode
982 Declaring the special data types containing all derivatives
983 \code
984  typedef Sacado::Fad::DFad<double> DFadType;
985  SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > eps_fad, eps_fad_squared;
986  Sacado::Fad::DFad<DFadType> phi_fad;
987 
988 
989 \endcode
990 Setting the dofs
991 \code
992  for ( unsigned int x=0; x<6; ++x )
993  {
994  unsigned int i=std_map_indicies[x].first;
995  unsigned int j=std_map_indicies[x].second;
996  (eps_fad[i][j]).diff( x, nbr_dofs); // set up the "inner" derivatives
997  (eps_fad[i][j]).val() = fad_double(nbr_dofs, x, eps[i][j]); // set up the "outer" derivatives
998  }
999 
1000  phi_fad.diff( 6, nbr_dofs );
1001  phi_fad.val() = fad_double(nbr_dofs, 6, phi); // set up the "outer" derivatives
1002 
1003  std::cout << "eps_fad=" << eps_fad << std::endl;
1004  std::cout << "phi_fad=" << phi_fad << std::endl;
1005 
1006 \endcode
1007 Compute eps² = eps_ij * eps_jk in index notation
1008 \code
1009  for ( unsigned int i=0; i<dim; ++i)
1010  for ( unsigned int k=0; k<dim; ++k )
1011  for ( unsigned int j=0; j<dim; ++j )
1012  if ( i>=k )
1013  eps_fad_squared[i][k] += eps_fad[i][j] * eps_fad[j][k];
1014 
1015 \endcode
1016 Compute the strain energy density
1017 \code
1018  Sacado::Fad::DFad<DFadType> energy;
1019  energy = lambda/2. * trace(eps_fad)*trace(eps_fad) + mu * trace(eps_fad_squared) + 25 * phi_fad * trace(eps_fad);
1020 
1021 \endcode
1022 Give some insight into the storage of the values and derivatives
1023 \code
1024  std::cout << "energy=" << energy << std::endl;
1025 
1026 \endcode
1027 Compute sigma as \f[ \frac{\partial \Psi}{\partial \boldsymbol{\varepsilon}} \f]
1028 \code
1029  SymmetricTensor<2,dim> sigma_Sac;
1030  for ( unsigned int x=0; x<6; ++x )
1031  {
1032  unsigned int i=std_map_indicies[x].first;
1033  unsigned int j=std_map_indicies[x].second;
1034  if ( i!=j )
1035  sigma_Sac[i][j] = 0.5 * energy.dx(x).val();
1036  else
1037  sigma_Sac[i][j] = energy.dx(x).val();
1038  }
1039  std::cout << "sigma_Sacado=" << sigma_Sac << std::endl;
1040 
1041  double d_energy_d_phi = energy.dx(6).val();
1042  std::cout << "d_energy_d_phi=" << d_energy_d_phi << std::endl;
1043 
1044  double d2_energy_d_phi_2 = energy.dx(6).dx(6);
1045  std::cout << "d2_energy_d_phi_2=" << d2_energy_d_phi_2 << std::endl;
1046 
1047 \endcode
1048 Analytical stress tensor:
1049 \code
1050  SymmetricTensor<2,dim> sigma;
1051  sigma = lambda*trace(eps)*unit_symmetric_tensor<dim>() + 2. * mu * eps + 25 * phi * unit_symmetric_tensor<dim>();
1052  std::cout << "analy. sigma=" << sigma << std::endl;
1053 
1054 
1055 \endcode
1056 Sacado-Tangent
1057 \code
1058  SymmetricTensor<4,dim> C_Sac;
1059  for(unsigned int x=0;x<6;++x)
1060  for(unsigned int y=0;y<6;++y)
1061  {
1062  const unsigned int i=std_map_indicies[y].first;
1063  const unsigned int j=std_map_indicies[y].second;
1064  const unsigned int k=std_map_indicies[x].first;
1065  const unsigned int l=std_map_indicies[x].second;
1066 
1067  double deriv = energy.dx(x).dx(y); // Access the derivatives of the (i,j)-th component of \a sigma
1068 
1069  if ( k!=l && i!=j )
1070  C_Sac[i][j][k][l] = 0.25* deriv;
1071  else if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
1072  {
1073  C_Sac[i][j][k][l] = 0.5*deriv;
1074  C_Sac[i][j][l][k] = 0.5*deriv;
1075  }
1076  else
1077  C_Sac[i][j][k][l] = deriv;
1078  }
1079 
1080 \endcode
1081 Analytical tangent
1082 \code
1083  SymmetricTensor<4,dim> C_analy;
1084  C_analy = lambda * outer_product(unit_symmetric_tensor<dim>(), unit_symmetric_tensor<dim>()) + 2. * mu * identity_tensor<dim>();
1085 
1086  double error_Sacado_vs_analy=0;
1087  for (unsigned int i=0; i<dim; ++i)
1088  for ( unsigned int j=0; j<dim; ++j)
1089  for ( unsigned int k=0; k<dim; ++k)
1090  for ( unsigned int l=0; l<dim; ++l)
1091  error_Sacado_vs_analy += std::fabs(C_Sac[i][j][k][l] - C_analy[i][j][k][l]);
1092 
1093  std::cout << "Numerical error=" << error_Sacado_vs_analy << std::endl;
1094 }
1095 
1096 
1097 
1098 \endcode
1099 @section Ex8 8. Example: First and second derivatives - Using the Wrapper
1100 \code
1101 void sacado_test_8 ()
1102 
1103 {
1104  const unsigned int dim=3;
1105 
1106  std::cout << "Test 8:" << std::endl;
1107 
1108 \endcode
1109 Defining the inputs (material parameters, strain tensor)
1110 \code
1111  double lambda=1;
1112  double mu=2;
1113  SymmetricTensor<2,dim, double> eps;
1114  double phi = 0.3;
1115 
1116  eps[0][0] = 1.;
1117  eps[1][1] = 2.;
1118  eps[2][2] = 3.;
1119 
1120  eps[0][1] = 4.;
1121  eps[0][2] = 5.;
1122  eps[1][2] = 6.;
1123 
1124 \endcode
1125 Declaring the special data types containing all derivatives
1126 \code
1127  typedef Sacado::Fad::DFad<double> DFadType;
1128 
1129 \endcode
1130 Declare the variables \a eps_fad and \a phi_fad as the special Wrapper data types
1131 \code
1134 
1135 \endcode
1136 Declare the summary data type relating all the dofs and initialising them too
1137 \code
1138  Sacado_Wrapper::DoFs_summary<dim> DoFs_summary;
1139  DoFs_summary.init_set_dofs(eps_fad, eps, phi_fad, phi);
1140 
1141 \endcode
1142 The variables are outputted to give some insight into the storage of the values (derivatives still trivial).
1143 \code
1144  std::cout << "eps_fad=" << eps_fad << std::endl;
1145  std::cout << "phi_fad=" << phi_fad << std::endl;
1146 
1147 \endcode
1148 Compute eps² = eps_ij * eps_jk in index notation
1149 \code
1150  SymmetricTensor<2,dim, Sacado::Fad::DFad<DFadType> > eps_fad_squared;
1151  for ( unsigned int i=0; i<dim; ++i)
1152  for ( unsigned int k=0; k<dim; ++k )
1153  for ( unsigned int j=0; j<dim; ++j )
1154  if ( i>=k )
1155  eps_fad_squared[i][k] += eps_fad[i][j] * eps_fad[j][k];
1156 
1157 \endcode
1158 Compute the strain energy density
1159 \code
1160  Sacado::Fad::DFad<DFadType> energy;
1161  energy = lambda/2. * trace(eps_fad)*trace(eps_fad) + mu * trace(eps_fad_squared) + 25 * phi_fad * trace(eps_fad);
1162 
1163 \endcode
1164 The energy is outputted (formatted by hand) to give some insight into the storage of the values and derivatives. \n
1165 energy=399 [ 17.5 32 40 21.5 48 25.5 150 ] \n
1166 [ 17.5 [ 5 0 0 1 0 1 25 ] 32 [ 0 8 0 0 0 0 0 ] 40 [ 0 0 8 0 0 0 0 ] \n
1167 21.5 [ 1 0 0 5 0 1 25 ] 48 [ 0 0 0 0 8 0 0 ] 25.5 [ 1 0 0 1 0 5 25 ] \n
1168 150 [ 25 0 0 25 0 25 0 ] ]
1169 
1170 \code
1171  std::cout << "energy=" << energy << std::endl;
1172 
1173 \endcode
1174 Compute sigma as \f[ \boldsymbol{\sigma} = \frac{\partial \Psi}{\partial \boldsymbol{\varepsilon}} \f]
1175 \code
1176  SymmetricTensor<2,dim> sigma_Sac;
1177  eps_fad.get_tangent(sigma_Sac, energy);
1178  std::cout << "sigma_Sacado=" << sigma_Sac << std::endl;
1179 
1180  double d_energy_d_phi;
1181  phi_fad.get_tangent(d_energy_d_phi, energy);
1182  std::cout << "d_energy_d_phi=" << d_energy_d_phi << std::endl;
1183 
1184 
1185 \endcode
1186 Analytical stress tensor:
1187 \code
1188  SymmetricTensor<2,dim> sigma;
1189  sigma = lambda*trace(eps)*unit_symmetric_tensor<dim>() + 2. * mu * eps + 25 * phi * unit_symmetric_tensor<dim>();
1190  std::cout << "analy. sigma=" << sigma << std::endl;
1191 
1192 
1193 \endcode
1194 Sacado stress tangent (or eps curvature) as \f[ \frac{\partial^2 \Psi}{\partial \boldsymbol{\varepsilon}^2} \f]
1195 \code
1196  SymmetricTensor<4,dim> C_Sac;
1197  eps_fad.get_curvature(C_Sac, energy);
1198 
1199 \endcode
1200 Sacado phi curvature as \f[ \frac{\partial^2 \Psi}{\partial \varphi^2} \f]
1201 \code
1202  double d2_energy_d_phi_2;
1203  phi_fad.get_curvature(d2_energy_d_phi_2, energy);
1204  std::cout << "d2_energy_d_phi_2=" << d2_energy_d_phi_2 << std::endl;
1205 
1206 \endcode
1207 Sacado derivatives \f[ \frac{\partial^2 \Psi}{\partial \boldsymbol{\varepsilon} \partial \varphi} \f]
1208 \code
1209  SymmetricTensor<2,dim> d2_energy_d_eps_d_phi;
1210  DoFs_summary.get_curvature(d2_energy_d_eps_d_phi, energy, eps_fad, phi_fad);
1211  std::cout << "d2_energy_d_eps_d_phi=" << d2_energy_d_eps_d_phi << std::endl;
1212 
1213 \endcode
1214 Sacado derivatives \f[ \frac{\partial^2 \Psi}{\partial \varphi \partial \boldsymbol{\varepsilon}} \f]
1215 \code
1216  SymmetricTensor<2,dim> d2_energy_d_phi_d_eps;
1217  DoFs_summary.get_curvature(d2_energy_d_phi_d_eps, energy, phi_fad, eps_fad);
1218  std::cout << "d2_energy_d_phi_d_eps=" << d2_energy_d_phi_d_eps << std::endl;
1219 
1220 \endcode
1221 When you consider the output: \n
1222 d2_energy_d_eps_d_phi=25 0 0 0 25 0 0 0 25 \n
1223 d2_energy_d_phi_d_eps=25 0 0 0 25 0 0 0 25 \n
1224 in detail you will notice that both second derivatives are identical. This compplies with the Schwarz integrability condition (Symmetry of second derivatives)
1225 (ignoring all limitation and requirements), it holds
1226 \f[ \frac{\partial^2 \Psi}{\partial \boldsymbol{\varepsilon} \partial \varphi} = \frac{\partial^2 \Psi}{\partial \varphi \partial \boldsymbol{\varepsilon}} \f]
1227 
1228 Analytical stress tangent
1229 \code
1230  SymmetricTensor<4,dim> C_analy;
1231  C_analy = lambda * outer_product(unit_symmetric_tensor<dim>(), unit_symmetric_tensor<dim>()) + 2. * mu * identity_tensor<dim>();
1232 
1233 \endcode
1234 Compute the error for the stress tangent
1235 \code
1236  double error_Sacado_vs_analy=0;
1237  for (unsigned int i=0; i<dim; ++i)
1238  for ( unsigned int j=0; j<dim; ++j)
1239  for ( unsigned int k=0; k<dim; ++k)
1240  for ( unsigned int l=0; l<dim; ++l)
1241  error_Sacado_vs_analy += std::fabs(C_Sac[i][j][k][l] - C_analy[i][j][k][l]);
1242  std::cout << "Numerical error=" << error_Sacado_vs_analy << std::endl;
1243 }
1244 
1245 
1246 \endcode
1247 @section Ex9 9. Example: Why we sometimes need the factor of 0.5 in the derivatives and sometimes we don't
1248 \code
1249 void sacado_test_9 ()
1250 
1251 {
1252  const unsigned int dim=3;
1253 
1254  std::cout << "Test 9:" << std::endl;
1255 
1256  double kappa_param = 5;
1257  fad_double kappa (kappa_param);
1258  double mu = 2;
1259 
1260  SymmetricTensor<2,dim, fad_double> sigma, eps;
1261 
1262  std::map<unsigned int,std::pair<unsigned int,unsigned int>> std_map_indicies;
1263 
1264  eps[0][0] = 1;
1265  eps[1][1] = 2;
1266  eps[2][2] = 3;
1267 
1268  eps[0][1] = 4;
1269  eps[0][2] = 5;
1270  eps[1][2] = 6;
1271 
1272 
1273  eps[0][0].diff(0,6);
1274  eps[0][1].diff(1,6);
1275  eps[0][2].diff(2,6);
1276  eps[1][1].diff(3,6);
1277  eps[1][2].diff(4,6);
1278  eps[2][2].diff(5,6);
1279 
1280  std::pair<unsigned int, unsigned int> tmp_pair;
1281  tmp_pair.first=0; tmp_pair.second=0;
1282  std_map_indicies[0] = tmp_pair;
1283 
1284  tmp_pair.first=0; tmp_pair.second=1;
1285  std_map_indicies[1] = tmp_pair;
1286 
1287  tmp_pair.first=0; tmp_pair.second=2;
1288  std_map_indicies[2] = tmp_pair;
1289 
1290  tmp_pair.first=1; tmp_pair.second=1;
1291  std_map_indicies[3] = tmp_pair;
1292 
1293  tmp_pair.first=1; tmp_pair.second=2;
1294  std_map_indicies[4] = tmp_pair;
1295 
1296  tmp_pair.first=2; tmp_pair.second=2;
1297  std_map_indicies[5] = tmp_pair;
1298 
1299  SymmetricTensor<2,dim, fad_double> stdTensor_I (( unit_symmetric_tensor<dim,fad_double>()) );
1300 
1301  sigma = kappa * (trace(eps) * stdTensor_I);
1302  SymmetricTensor<2,dim,fad_double> tmp = deviator<dim,fad_double>(symmetrize<dim,fad_double>(eps)); tmp*=(mu*2);
1303  sigma += tmp;
1304 
1305  std::cout << "sigma=" << sigma << std::endl;
1306 
1307 \endcode
1308 Retrieve the values stored in \a sigma:
1309 \code
1310  SymmetricTensor<2,dim> sigma_d;
1311  for ( unsigned int i=0; i<dim; ++i)
1312  for ( unsigned int j=0; j<dim; ++j )
1313  sigma_d[i][j] = sigma[i][j].val();
1314 
1315  SymmetricTensor<4,dim> C_Sacado;
1316 
1317  for ( unsigned int i=0; i<dim; ++i)
1318  for ( unsigned int j=0; j<dim; ++j )
1319  {
1320  double *derivs = &sigma[i][j].fastAccessDx(0); // Access the derivatives of the (i,j)-th component of \a sigma
1321 
1322  for(unsigned int x=0;x<((dim==2)?3:6);++x)
1323  {
1324  unsigned int k=std_map_indicies[x].first;
1325  unsigned int l=std_map_indicies[x].second;
1326 
1327  if(k!=l)/*Compare to Voigt notation since only SymmetricTensor instead of Tensor*/
1328  {
1329  C_Sacado[i][j][k][l] = 0.5*derivs[x];
1330  C_Sacado[i][j][l][k] = 0.5*derivs[x];
1331  }
1332  else
1333  C_Sacado[i][j][k][l] = derivs[x];
1334  }
1335  }
1336 
1337  SymmetricTensor<2,dim> eps_d;
1338  eps_d[0][0] = 1;
1339  eps_d[1][1] = 2;
1340  eps_d[2][2] = 3;
1341 
1342  eps_d[0][1] = 4;
1343  eps_d[0][2] = 5;
1344  eps_d[1][2] = 6;
1345 
1346  SymmetricTensor<2,dim> stress_from_tangent = C_Sacado*eps_d;
1347 
1348  std::cout << "C_Sacado*eps_d=" << stress_from_tangent << std::endl;
1349  std::cout << "sigma= " << sigma_d << std::endl;
1350 
1351 \endcode
1352 @todo Add the link to the ac.nz site
1353 
1354 As you can see we obtain the correct stress tensor only by using the factor 0.5 onto the off-diagonal elements of the symmetric tensor. \n
1355 @note So, why do we need the factor 0.5 then? (based on [homepages.engineering.auckland.ac.nz ...])\n
1356 When using symmetric tensor one has to be aware of the definition of the independent variables. In 3d those are the
1357 following 6 components: \n
1358 - \f$ A_{11} \f$
1359 - \f$ A_{22} \f$
1360 - \f$ A_{33} \f$
1361 - \f$ \overline{A_{12}} \f$
1362 - \f$ \overline{A_{13}} \f$
1363 - \f$ \overline{A_{23}} \f$ \n
1364 Here the overlined components denote the averaged values, such that
1365 - \f$ \overline{A_{12}} = 0.5 \cdot [ A_{12} + A_{21} ] = A_{12} = A_{21} \f$
1366 - ... \n
1367 This looks rather boring, because the averaged value of the two identical components (symmetry) is identical to the components.
1368 But (and that is a crucial BUT), this only holds for the values NOT the derivatives.
1369 The derivative of a function \f$ \Phi \f$ with respect to a symmetric tensor \a A is
1370 \f[ \frac{\partial \Phi}{\partial \overline{A_{12}}} = \frac{\partial \Phi}{\partial A_{12}} \cdot \frac{\partial A_{12}}{\partial \overline{A_{12}}}
1371  + \frac{\partial \Phi}{\partial A_{21}} \cdot \frac{\partial A_{21}}{\partial \overline{A_{21}}}
1372  = \frac{\partial \Phi}{\partial A_{12}} + \frac{\partial \Phi}{\partial A_{21}}
1373  = 2 \cdot \frac{\partial \Phi}{\partial A_{12}} \f]
1374 Note how the derivatives with respect to the actual tensor components (here 12 and 21) are identical. \n
1375 However, our derivatives (computed via Sacado) were calculated with respect to the actually independent
1376 variables (overlined). Still, in the tensor with the derivatives we want the derivatives with respect to the normal components (not overlined). \n
1377 Summary: \n
1378 When you compute derivatives with respect to a symmetric tensor, you have to scale the off-diagonal elements by 0.5 to account for the fact that you actually work with
1379 the averaged values. \n
1380 As a consequence, we use the factor 0.5 for all off-diagonal components of derivatives with respect to symmetric tensors. For second derivatives with respect to
1381 symmetric tensors the factors 0.5 and 0.25 come into play.
1382 
1383 \code
1384 }
1385 
1386 \endcode
1387 @note Here we template <dim,Number> for a function "on the same level".
1388 If you template a member function with <Number> inside a class template <dim>,
1389 use \n
1390 @code
1391 template<int dim>
1392 template<typename Number>
1393 SymmetricTensor<2,dim,Number> ...
1394 @endcode
1395 For the latter the declaration of the member function in the class body looks like this:
1396 @code
1397  template<typename Number>
1398 SymmetricTensor<2,dim,Number> ...
1399 @endcode
1400 
1401 \code
1402 template<int dim, typename Number>
1403 SymmetricTensor<2,dim,Number> stress_strain_relation ( const SymmetricTensor<2,dim,Number> &eps, const double &kappa, const double &mu )
1404 {
1405  SymmetricTensor<2,dim,Number> sigma;
1406 
1407  SymmetricTensor<2,dim,Number> stdTensor_I (( unit_symmetric_tensor<dim,Number>()) );
1408 
1409 \endcode
1410 Our stress equation is now computed in index notation to simplify the use of the constants and
1411 especially the use of the \a deviator.
1412 \code
1413  for ( unsigned int i=0; i<dim; ++i)
1414  for ( unsigned int j=0; j<dim; ++j )
1415  sigma[i][j] = kappa * trace(eps) * stdTensor_I[i][j] + 2. * mu * deviator(eps)[i][j];
1416 
1417  return sigma;
1418 }
1419 
1420 
1421 \endcode
1422 @section Ex10 10. Example: Using the wrapper with templated functions
1423 \code
1424 void sacado_test_10 ()
1425 {
1426  std::cout << "Test 10:" << std::endl;
1427  const unsigned int dim=3;
1428 \endcode
1429 To measure the computation time of certain sections (see below)
1430 \code
1431  TimerOutput timer (std::cout, TimerOutput::summary, TimerOutput::cpu_times);
1432 
1433 \endcode
1434 The following declarations are usually input arguments. So you receive the strain tensor and the constants out of doubles.
1435 \code
1436  SymmetricTensor<2,dim> eps_d;
1437  eps_d[0][0] = 1;
1438  eps_d[1][1] = 2;
1439  eps_d[2][2] = 3;
1440 
1441  eps_d[0][1] = 4;
1442  eps_d[0][2] = 5;
1443  eps_d[1][2] = 6;
1444 
1445  double kappa = 5;
1446  double mu = 2;
1447 
1448 \endcode
1449 Now we start working with Sacado: \n
1450 When we use the index notation to compute e.g. our stress we do not need to declare our constants (here kappa, mu) as
1451 fad_double.
1452 
1453 We declare our strain tensor as the special data type Sacado_Wrapper::SymTensor from the file "Sacado_Wrapper.h"
1454 where this data type was derived from the SymmetricTensor<2,dim,fad_double>.
1455 \code
1456  Sacado_Wrapper::SymTensor<dim> eps;
1457 
1458 \endcode
1459 Next we initialize our Sacado strain tensor with the values of the inputed double strain tensor:
1460 \code
1461  eps.init(eps_d);
1462 
1463 \endcode
1464 We define all the entries in the symmetric tensor \a eps as the dofs. So we can later derive any variable
1465 with respect to the strain tensor \a eps.
1466 @note Keep the order first declare \a eps, then initialise it and afterwards set it as dofs
1467 
1468 \code
1469  eps.set_dofs();
1470 
1471 \endcode
1472 Now we call the templated function via the Sacado data types to get the resulting stress tensor \a sigma. \n
1473 To give you an impression of how much longer computations with the Sacado data type take,
1474 we also measure the time needed for the computation
1475 \code
1476  timer.enter_subsection("Stress via fad_double");
1477  SymmetricTensor<2,dim,fad_double> sigma = stress_strain_relation ( eps, kappa, mu );
1478  timer.leave_subsection();
1479  std::cout << "Stress via fad_double:" << extract_value_from_Sacado<dim> (sigma) << std::endl;
1480 
1481 \endcode
1482 Finally we declare our desired tangent as the fourth order tensor \a C_Sacado and compute the tangent via
1483 the command \a get_tangent.
1484 \code
1485  SymmetricTensor<4,dim> C_Sacado;
1486  eps.get_tangent(C_Sacado, sigma);
1487 
1488 \endcode
1489 For comparsion we can also compute the pure values of the stress tensor by calling
1490 the templated function with the normal data type \a double.
1491 \code
1492  timer.enter_subsection("Stress via double");
1493  SymmetricTensor<2,dim,double> sigma_d = stress_strain_relation ( eps_d, kappa, mu );
1494  timer.leave_subsection();
1495  std::cout << "Stress via double: " << sigma_d << std::endl;
1496 
1497 \endcode
1498 As you can see from the times shown in the command window after you ran the code.
1499 The computation via the Sacado data types takes here around 4...6 times longer than
1500 using the standard type \a double. However, you have to keep in mind that Sacado
1501 computes the tangent whilst computing the stress. Depending on how complicated and computationally demanding the
1502 analytical tangent is, Sacado can possibly compete.
1503 
1504 @note
1505 This also means that you should only use the Sacado data types when you are actually
1506 interested in the derivatives. Else the computation is just too expensive.
1507 
1508 \code
1509 }
1510 
1511 
1512 /*
1513  * The main function just calls all the examples and puts some space between the outputs.
1514  */
1515 int main ()
1516 {
1517  sacado_test_scalar ();
1518 
1519  std::cout << std::endl;
1520 
1521  sacado_test_2 ();
1522 
1523  std::cout << std::endl;
1524 
1525  sacado_test_3 ();
1526 
1527  std::cout << std::endl;
1528 
1529  sacado_test_3B ();
1530 
1531  std::cout << std::endl;
1532 
1533  sacado_test_4();
1534 
1535  std::cout << std::endl;
1536 
1537  sacado_test_5();
1538 
1539  std::cout << std::endl;
1540 
1541  sacado_test_6();
1542 
1543  std::cout << std::endl;
1544 
1545  sacado_test_7();
1546 
1547  std::cout << std::endl;
1548 
1549  sacado_test_8();
1550 
1551  std::cout << std::endl;
1552 
1553  sacado_test_9();
1554 
1555  std::cout << std::endl;
1556 
1557  sacado_test_10();
1558 }
1559 \endcode
1560 
1561 @section END The End
1562 
1563 Hosted via GitHub according to https://goseeky.wordpress.com/2017/07/22/documentation-101-doxygen-with-github-pages/ \n
1564 Design of the documentation inspired by the deal.ii tutorial programs and created as shown here: https://github.com/jfriedlein/Custom_Doxygen_Documentation
1565 
1566 @}
1567 */
void sacado_test_5()
Definition: Sacado_example.cc:572
void sacado_test_scalar()
Definition: Sacado_example.cc:43
void sacado_test_6()
Definition: Sacado_example.cc:605
Sacado::Fad::DFad< double > DFadType
Definition: Sacado_Wrapper.h:19
void sacado_test_8()
Definition: Sacado_example.cc:807
void sacado_test_4()
Definition: Sacado_example.cc:426
void sacado_test_2()
Definition: Sacado_example.cc:72
Definition: Sacado_Wrapper.h:395
void get_curvature(double &Curvature, Sacado::Fad::DFad< DFadType > &argument)
Definition: Sacado_Wrapper.h:456
Definition: Sacado_Wrapper.h:505
Sacado::Fad::DFad< double > fad_double
Definition: Sacado-auxiliary_functions.h:7
void get_tangent(double &Tangent, Sacado::Fad::DFad< DFadType > &argument)
Definition: Sacado_Wrapper.h:436
Definition: Sacado_Wrapper.h:63
Definition: Sacado_Wrapper.h:193
void sacado_test_7()
Definition: Sacado_example.cc:659
void init_set_dofs(SymTensor2< dim > &eps, SymmetricTensor< 2, dim > &eps_init, SW_double2< dim > &double_arg, double &double_init)
Definition: Sacado_Wrapper.h:533