]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added template Number to DerivativeForms. 936/head
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 15 May 2015 22:55:55 +0000 (00:55 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Fri, 15 May 2015 23:41:52 +0000 (01:41 +0200)
doc/news/changes.h
include/deal.II/base/derivative_form.h
tests/base/derivative_forms_01.cc [new file with mode: 0644]
tests/base/derivative_forms_01.with_trilinos=on.output [new file with mode: 0644]

index 32fa7b998c2263d2c11825b9261c5a0f75fbf1b8..1d9097ce2931ed68d34d9ff38b60bb65698a8841 100644 (file)
@@ -377,6 +377,12 @@ inconvenience this causes.
 
 
 <ol>
+  <li> New: DerivativeForm() now takes an additional optional
+  template argument specifying the type, similarly to Tensor() classes.
+  <br>
+  (Luca Heltai, 2015/05/16)
+  </li>
+
   <li> New: Utilities::MPI::min() functions.
   <br>
   (Timo Heister, 2015/05/12)
index e6a4b57af11b9eb09eb3f111eb765cd1f9df6059..69d3820fabb33abe99a1f05988474968704b9fe0 100644 (file)
@@ -29,11 +29,11 @@ DEAL_II_NAMESPACE_OPEN
  * R}^{\text{spacedim}}$, the second derivative a bilinear map from  ${\mathbb
  * R}^{\text{dim}} \times  {\mathbb R}^{\text{dim}}$ to ${\mathbb
  * R}^{\text{spacedim}}$ and so on.  In deal.II we represent these derivaties
- * using objects of type DerivativeForm<1,dim,spacedim>,
- * DerivativeForm<2,dim,spacedim> and so on.
- * @author Sebastian Pauletti, 2011
+ * using objects of type DerivativeForm<1,dim,spacedim,Number>,
+ * DerivativeForm<2,dim,spacedim,Number> and so on.
+ * @author Sebastian Pauletti, 2011, Luca Heltai, 2015
  */
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number=double>
 class DerivativeForm
 {
 public:
@@ -45,50 +45,50 @@ public:
   /**
    * Constructor from a second order tensor.
    */
-  DerivativeForm (const Tensor<2,dim> &);
+  DerivativeForm (const Tensor<2,dim,Number> &);
 
   /**
    * Read-Write access operator.
    */
-  Tensor<order,dim> &operator [] (const unsigned int i);
+  Tensor<order,dim,Number> &operator [] (const unsigned int i);
 
   /**
    * Read-only access operator.
    */
-  const Tensor<order,dim> &operator [] (const unsigned int i) const;
+  const Tensor<order,dim,Number> &operator [] (const unsigned int i) const;
 
   /**
    * Assignment operator.
    */
-  DerivativeForm   &operator = (const DerivativeForm <order, dim, spacedim> &);
+  DerivativeForm   &operator = (const DerivativeForm <order, dim, spacedim, Number> &);
 
 
   /**
    * Assignment operator.
    */
-  DerivativeForm   &operator = (const Tensor<2,dim> &);
+  DerivativeForm   &operator = (const Tensor<2,dim, Number> &);
 
   /**
    * Assignment operator.
    */
-  DerivativeForm   &operator = (const Tensor<1,dim> &);
+  DerivativeForm   &operator = (const Tensor<1,dim, Number> &);
 
   /**
-   * Converts a DerivativeForm <1,dim, dim> to Tensor<2,dim>. If the
+   * Converts a DerivativeForm <1,dim, dim> to Tensor<2,dim,Number>. If the
    * derivative is the Jacobian of F, then Tensor[i] = grad(F^i).
    */
-  operator Tensor<2,dim>() const;
+  operator Tensor<2,dim,Number>() const;
 
   /**
-   * Converts a DerivativeForm <1, dim, 1> to Tensor<1,dim>.
+   * Converts a DerivativeForm <1, dim, 1> to Tensor<1,dim,Number>.
    */
-  operator Tensor<1,dim>() const;
+  operator Tensor<1,dim,Number>() const;
 
   /**
    * Return the transpose of a rectangular DerivativeForm, that is to say
    * viewed as a two dimensional matrix.
    */
-  DerivativeForm<1, spacedim, dim> transpose () const;
+  DerivativeForm<1, spacedim, dim, Number> transpose () const;
 
 
   /**
@@ -104,7 +104,7 @@ public:
    * covariant matrix, namely $DF*G^{-1}$, where $G = DF^{t}*DF$. If $DF$ is
    * square, covariant from gives $DF^{-t}$.
    */
-  DerivativeForm<1, dim, spacedim> covariant_form() const;
+  DerivativeForm<1, dim, spacedim, Number> covariant_form() const;
 
 
 
@@ -127,14 +127,14 @@ private:
   /**
    * Auxiliary function that computes (*this) * T^{t}
    */
-  DerivativeForm<1, dim, spacedim> times_T_t (Tensor<2,dim> T) const;
+  DerivativeForm<1, dim, spacedim, Number> times_T_t (Tensor<2,dim,Number> T) const;
 
 
 private:
   /**
    * Array of tensors holding the subelements.
    */
-  Tensor<order,dim> tensor[spacedim];
+  Tensor<order,dim,Number> tensor[spacedim];
 
 
 };
@@ -147,9 +147,9 @@ private:
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim>::DerivativeForm  ()
+DerivativeForm<order, dim, spacedim, Number>::DerivativeForm  ()
 {
 // default constructor. not specifying an initializer list calls
 // the default constructor of the subobjects, which initialize them
@@ -158,9 +158,9 @@ DerivativeForm<order, dim, spacedim>::DerivativeForm  ()
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim>::DerivativeForm(const Tensor<2,dim> &T)
+DerivativeForm<order, dim, spacedim, Number>::DerivativeForm(const Tensor<2,dim,Number> &T)
 {
   Assert( (dim == spacedim) && (order==1),
           ExcMessage("Only allowed for square tensors."));
@@ -171,11 +171,11 @@ DerivativeForm<order, dim, spacedim>::DerivativeForm(const Tensor<2,dim> &T)
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim> &
-DerivativeForm<order, dim, spacedim>::
-operator = (const DerivativeForm<order, dim, spacedim> &ta)
+DerivativeForm<order, dim, spacedim, Number> &
+DerivativeForm<order, dim, spacedim, Number>::
+operator = (const DerivativeForm<order, dim, spacedim, Number> &ta)
 {
   for (unsigned int j=0; j<spacedim; ++j)
     (*this)[j] = ta[j];
@@ -184,10 +184,10 @@ operator = (const DerivativeForm<order, dim, spacedim> &ta)
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim> &DerivativeForm<order, dim, spacedim>::
-operator = (const Tensor<2,dim> &ta)
+DerivativeForm<order, dim, spacedim, Number> &DerivativeForm<order, dim, spacedim, Number>::
+operator = (const Tensor<2,dim,Number> &ta)
 {
   Assert( (dim == spacedim) && (order==1),
           ExcMessage("Only allowed for square tensors."));
@@ -201,10 +201,10 @@ operator = (const Tensor<2,dim> &ta)
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim> &DerivativeForm<order, dim, spacedim>::
-operator = (const Tensor<1,dim> &T)
+DerivativeForm<order, dim, spacedim, Number> &DerivativeForm<order, dim, spacedim, Number>::
+operator = (const Tensor<1,dim,Number> &T)
 {
   Assert( (1 == spacedim) && (order==1),
           ExcMessage("Only allowed for spacedim==1 and order==1."));
@@ -217,9 +217,9 @@ operator = (const Tensor<1,dim> &T)
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-Tensor<order,dim> &DerivativeForm<order, dim, spacedim>::
+Tensor<order,dim,Number> &DerivativeForm<order, dim, spacedim, Number>::
 operator[] (const unsigned int i)
 {
   Assert (i<spacedim, ExcIndexRange(i, 0, spacedim));
@@ -229,9 +229,9 @@ operator[] (const unsigned int i)
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-const Tensor<order,dim> &DerivativeForm<order, dim, spacedim>::
+const Tensor<order,dim,Number> &DerivativeForm<order, dim, spacedim, Number>::
 operator[] (const unsigned int i) const
 {
   Assert (i<spacedim, ExcIndexRange(i, 0, spacedim));
@@ -241,9 +241,9 @@ operator[] (const unsigned int i) const
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim>::operator Tensor<1,dim>() const
+DerivativeForm<order, dim, spacedim, Number>::operator Tensor<1,dim,Number>() const
 {
   Assert( (1 == spacedim) && (order==1),
           ExcMessage("Only allowed for spacedim==1."));
@@ -254,16 +254,14 @@ DerivativeForm<order, dim, spacedim>::operator Tensor<1,dim>() const
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<order, dim, spacedim>::operator Tensor<2,dim>() const
+DerivativeForm<order, dim, spacedim, Number>::operator Tensor<2,dim,Number>() const
 {
   Assert( (dim == spacedim) && (order==1),
           ExcMessage("Only allowed for square tensors."));
 
-
-
-  Tensor<2,dim> t;
+  Tensor<2,dim,Number> t;
 
   if ((dim == spacedim) && (order==1))
     for (unsigned int j=0; j<dim; ++j)
@@ -275,14 +273,14 @@ DerivativeForm<order, dim, spacedim>::operator Tensor<2,dim>() const
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<1,spacedim,dim>
-DerivativeForm<order,dim,spacedim>::
+DerivativeForm<1,spacedim,dim,Number>
+DerivativeForm<order,dim,spacedim,Number>::
 transpose () const
 {
   Assert(order==1, ExcMessage("Only for rectangular DerivativeForm."));
-  DerivativeForm<1,spacedim,dim> tt;
+  DerivativeForm<1,spacedim,dim,Number> tt;
 
   for (unsigned int i=0; i<spacedim; ++i)
     for (unsigned int j=0; j<dim; ++j)
@@ -293,13 +291,13 @@ transpose () const
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<1, dim, spacedim>
-DerivativeForm<order,dim,spacedim>::times_T_t (Tensor<2,dim> T) const
+DerivativeForm<1, dim, spacedim,Number>
+DerivativeForm<order,dim,spacedim,Number>::times_T_t (Tensor<2,dim,Number> T) const
 {
   Assert( order==1, ExcMessage("Only for order == 1."));
-  DerivativeForm<1,dim, spacedim> dest;
+  DerivativeForm<1,dim, spacedim,Number> dest;
   for (unsigned int i=0; i<spacedim; ++i)
     for (unsigned int j=0; j<dim; ++j)
       dest[i][j] = (*this)[i] * T[j];
@@ -308,22 +306,22 @@ DerivativeForm<order,dim,spacedim>::times_T_t (Tensor<2,dim> T) const
 }
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
 double
-DerivativeForm<order,dim,spacedim>::determinant () const
+DerivativeForm<order,dim,spacedim,Number>::determinant () const
 {
   Assert( order==1, ExcMessage("Only for order == 1."));
   if (dim == spacedim)
     {
-      Tensor<2,dim> T = (Tensor<2,dim>)( (*this) );
+      Tensor<2,dim,Number> T = (Tensor<2,dim,Number>)( (*this) );
       return dealii::determinant(T);
     }
   else
     {
       Assert( spacedim>dim, ExcMessage("Only for spacedim>dim."));
       DerivativeForm<1,spacedim,dim> DF_t = this->transpose();
-      Tensor<2, dim> G; //First fundamental form
+      Tensor<2,dim,Number> G; //First fundamental form
       for (unsigned int i=0; i<dim; ++i)
         for (unsigned int j=0; j<dim; ++j)
           G[i][j] = DF_t[i] * DF_t[j];
@@ -336,16 +334,16 @@ DerivativeForm<order,dim,spacedim>::determinant () const
 
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
-DerivativeForm<1,dim,spacedim>
-DerivativeForm<order,dim,spacedim>::covariant_form() const
+DerivativeForm<1,dim,spacedim,Number>
+DerivativeForm<order,dim,spacedim,Number>::covariant_form() const
 {
 
   if (dim == spacedim)
     {
 
-      Tensor<2,dim> DF_t (dealii::transpose(invert(  (Tensor<2,dim>)(*this)   )));
+      Tensor<2,dim,Number> DF_t (dealii::transpose(invert(  (Tensor<2,dim,Number>)(*this)   )));
       DerivativeForm<1,dim, spacedim> result = DF_t;
       return (result);
     }
@@ -353,7 +351,7 @@ DerivativeForm<order,dim,spacedim>::covariant_form() const
     {
 
       DerivativeForm<1,spacedim,dim> DF_t = this->transpose();
-      Tensor<2, dim> G; //First fundamental form
+      Tensor<2,dim,Number> G; //First fundamental form
       for (unsigned int i=0; i<dim; ++i)
         for (unsigned int j=0; j<dim; ++j)
           G[i][j] = DF_t[i] * DF_t[j];
@@ -365,12 +363,12 @@ DerivativeForm<order,dim,spacedim>::covariant_form() const
 }
 
 
-template <int order, int dim, int spacedim>
+template <int order, int dim, int spacedim, typename Number>
 inline
 std::size_t
-DerivativeForm<order, dim, spacedim>::memory_consumption ()
+DerivativeForm<order, dim, spacedim, Number>::memory_consumption ()
 {
-  return sizeof(DerivativeForm<order, dim, spacedim>);
+  return sizeof(DerivativeForm<order, dim, spacedim, Number>);
 }
 
 #endif // DOXYGEN
@@ -387,13 +385,13 @@ DerivativeForm<order, dim, spacedim>::memory_consumption ()
  * @relates DerivativeForm
  * @author Sebastian Pauletti, 2011
  */
-template <int spacedim, int dim>
+template <int spacedim, int dim, typename Number>
 inline
-Tensor<1, spacedim>
-apply_transformation (const DerivativeForm<1,dim,spacedim> &DF,
-                      const Tensor<1,dim>               &T)
+Tensor<1,spacedim,Number>
+apply_transformation (const DerivativeForm<1,dim,spacedim,Number> &DF,
+                      const Tensor<1,dim,Number>               &T)
 {
-  Tensor<1, spacedim> dest;
+  Tensor<1,spacedim,Number> dest;
   for (unsigned int i=0; i<spacedim; ++i)
     dest[i] = DF[i] * T;
   return dest;
@@ -408,11 +406,11 @@ apply_transformation (const DerivativeForm<1,dim,spacedim> &DF,
  * @author Sebastian Pauletti, 2011
  */
 //rank=2
-template <int spacedim, int dim>
+template <int spacedim, int dim, typename Number>
 inline
 DerivativeForm<1, spacedim, dim>
-apply_transformation (const DerivativeForm<1,dim,spacedim> &DF,
-                      const Tensor<2,dim>               &T)
+apply_transformation (const DerivativeForm<1,dim,spacedim,Number> &DF,
+                      const Tensor<2,dim,Number>               &T)
 {
 
   DerivativeForm<1, spacedim, dim> dest;
@@ -428,13 +426,13 @@ apply_transformation (const DerivativeForm<1,dim,spacedim> &DF,
  * @relates DerivativeForm
  * @author Sebastian Pauletti, 2011
  */
-template <int spacedim, int dim>
+template <int spacedim, int dim, typename Number>
 inline
-Tensor<2, spacedim>
-apply_transformation (const DerivativeForm<1,dim,spacedim> &DF1,
-                      const DerivativeForm<1,dim,spacedim> &DF2)
+Tensor<2,spacedim,Number>
+apply_transformation (const DerivativeForm<1,dim,spacedim,Number> &DF1,
+                      const DerivativeForm<1,dim,spacedim,Number> &DF2)
 {
-  Tensor<2, spacedim> dest;
+  Tensor<2,spacedim,Number> dest;
 
   for (unsigned int i=0; i<spacedim; ++i)
     dest[i] = apply_transformation(DF1, DF2[i]);
@@ -450,12 +448,12 @@ apply_transformation (const DerivativeForm<1,dim,spacedim> &DF1,
  * @relates DerivativeForm
  * @author Sebastian Pauletti, 2011
  */
-template <int dim, int spacedim>
+template <int dim, int spacedim, typename Number>
 inline
-DerivativeForm<1,spacedim,dim>
-transpose (const DerivativeForm<1,dim,spacedim> &DF)
+DerivativeForm<1,spacedim,dim,Number>
+transpose (const DerivativeForm<1,dim,spacedim,Number> &DF)
 {
-  DerivativeForm<1,spacedim,dim> tt;
+  DerivativeForm<1,spacedim,dim,Number> tt;
   tt = DF.transpose();
   return tt;
 }
diff --git a/tests/base/derivative_forms_01.cc b/tests/base/derivative_forms_01.cc
new file mode 100644 (file)
index 0000000..1f9994d
--- /dev/null
@@ -0,0 +1,77 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2000 - 2014 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <deal.II/base/derivative_form.h>
+#include "Sacado.hpp"
+
+// Compute the derivative of F: R^dim -> R^spacedim
+// We construct the function F_i, i=[0,spacedim)
+// and we compute the first and second derivatives using 
+// Sacado
+
+typedef typename Sacado::Fad::DFad<double> Sdouble;
+typedef typename Sacado::Fad::DFad<Sdouble> SSdouble;
+
+template<int dim, int spacedim>
+void test() {
+  Tensor<1,dim,SSdouble> x;
+  for(unsigned int j=0; j<dim; ++j) {
+    // This is for the first derivative
+    Sdouble jv = numbers::PI*(j+1)/2.0;
+    jv.diff(j, dim);
+    // Now the second
+    x[j] = jv;
+    x[j].diff(j, dim);
+  }
+  
+  Tensor<1,spacedim,SSdouble> F;
+  for(unsigned int i=0; i<spacedim; ++i)
+    for(unsigned int j=0; j<dim; ++j)
+      F[i] += std::sin(x[j]+numbers::PI*i/2);
+
+  DerivativeForm<1,dim,spacedim,Sdouble> dF;
+  for(unsigned int i=0; i<spacedim; ++i)
+    for(unsigned int j=0; j<dim; ++j)
+      dF[i][j] = F[i].dx(j);
+
+  DerivativeForm<2,dim,spacedim, double> ddF;
+  for(unsigned int i=0; i<spacedim; ++i)
+    for(unsigned int j=0; j<dim; ++j)
+      for(unsigned int k=0; k<dim; ++k)
+       ddF[i][j][k] = dF[i][j].dx(k);
+  
+  deallog << "dim = " << dim
+         << ", spacedim = " << spacedim << std::endl;
+  deallog << "x  : " << x << std::endl;
+  for(unsigned int i=0; i<spacedim; ++i)
+    deallog << "F[" << i << "] : " << F[i] << std::endl;
+  for(unsigned int i=0; i<spacedim; ++i)
+    deallog << "dF[" << i << "] : " << dF[i] << std::endl;
+  for(unsigned int i=0; i<spacedim; ++i)
+    deallog << "ddF[" << i << "] : " << ddF[i] << std::endl;
+}
+
+int main()
+{
+  initlog();
+  test<1,1>();
+  test<1,2>();
+  test<1,3>();
+  test<2,2>();
+  test<2,3>();
+  test<3,3>();
+}
diff --git a/tests/base/derivative_forms_01.with_trilinos=on.output b/tests/base/derivative_forms_01.with_trilinos=on.output
new file mode 100644 (file)
index 0000000..0fa1613
--- /dev/null
@@ -0,0 +1,55 @@
+
+DEAL::dim = 1, spacedim = 1
+DEAL::x  : 1.57080 [ 1.00000 ] [ 1.00000 [ ] ]
+DEAL::F[0] : 1.00000 [ 6.12323e-17 ] [ 6.12323e-17 [ -1.00000 ] ]
+DEAL::dF[0] : 6.12323e-17 [ -1.00000 ]
+DEAL::ddF[0] : -1.00000
+DEAL::dim = 1, spacedim = 2
+DEAL::x  : 1.57080 [ 1.00000 ] [ 1.00000 [ ] ]
+DEAL::F[0] : 1.00000 [ 6.12323e-17 ] [ 6.12323e-17 [ -1.00000 ] ]
+DEAL::F[1] : 1.22465e-16 [ -1.00000 ] [ -1.00000 [ -1.22465e-16 ] ]
+DEAL::dF[0] : 6.12323e-17 [ -1.00000 ]
+DEAL::dF[1] : -1.00000 [ -1.22465e-16 ]
+DEAL::ddF[0] : -1.00000
+DEAL::ddF[1] : -1.22465e-16
+DEAL::dim = 1, spacedim = 3
+DEAL::x  : 1.57080 [ 1.00000 ] [ 1.00000 [ ] ]
+DEAL::F[0] : 1.00000 [ 6.12323e-17 ] [ 6.12323e-17 [ -1.00000 ] ]
+DEAL::F[1] : 1.22465e-16 [ -1.00000 ] [ -1.00000 [ -1.22465e-16 ] ]
+DEAL::F[2] : -1.00000 [ -1.83697e-16 ] [ -1.83697e-16 [ 1.00000 ] ]
+DEAL::dF[0] : 6.12323e-17 [ -1.00000 ]
+DEAL::dF[1] : -1.00000 [ -1.22465e-16 ]
+DEAL::dF[2] : -1.83697e-16 [ 1.00000 ]
+DEAL::ddF[0] : -1.00000
+DEAL::ddF[1] : -1.22465e-16
+DEAL::ddF[2] : 1.00000
+DEAL::dim = 2, spacedim = 2
+DEAL::x  : 1.57080 [ 1.00000 0.00000 ] [ 1.00000 [ ] 0.00000 [ ] ] 3.14159 [ 0.00000 1.00000 ] [ 0.00000 [ ] 1.00000 [ ] ]
+DEAL::F[0] : 1.00000 [ 6.12323e-17 -1.00000 ] [ 6.12323e-17 [ -1.00000 0.00000 ] -1.00000 [ 0.00000 -1.22465e-16 ] ]
+DEAL::F[1] : -1.00000 [ -1.00000 -1.83697e-16 ] [ -1.00000 [ -1.22465e-16 0.00000 ] -1.83697e-16 [ 0.00000 1.00000 ] ]
+DEAL::dF[0] : 6.12323e-17 [ -1.00000 0.00000 ] -1.00000 [ 0.00000 -1.22465e-16 ]
+DEAL::dF[1] : -1.00000 [ -1.22465e-16 0.00000 ] -1.83697e-16 [ 0.00000 1.00000 ]
+DEAL::ddF[0] : -1.00000 0.00000 0.00000 -1.22465e-16
+DEAL::ddF[1] : -1.22465e-16 0.00000 0.00000 1.00000
+DEAL::dim = 2, spacedim = 3
+DEAL::x  : 1.57080 [ 1.00000 0.00000 ] [ 1.00000 [ ] 0.00000 [ ] ] 3.14159 [ 0.00000 1.00000 ] [ 0.00000 [ ] 1.00000 [ ] ]
+DEAL::F[0] : 1.00000 [ 6.12323e-17 -1.00000 ] [ 6.12323e-17 [ -1.00000 0.00000 ] -1.00000 [ 0.00000 -1.22465e-16 ] ]
+DEAL::F[1] : -1.00000 [ -1.00000 -1.83697e-16 ] [ -1.00000 [ -1.22465e-16 0.00000 ] -1.83697e-16 [ 0.00000 1.00000 ] ]
+DEAL::F[2] : -1.00000 [ -1.83697e-16 1.00000 ] [ -1.83697e-16 [ 1.00000 0.00000 ] 1.00000 [ 0.00000 2.44929e-16 ] ]
+DEAL::dF[0] : 6.12323e-17 [ -1.00000 0.00000 ] -1.00000 [ 0.00000 -1.22465e-16 ]
+DEAL::dF[1] : -1.00000 [ -1.22465e-16 0.00000 ] -1.83697e-16 [ 0.00000 1.00000 ]
+DEAL::dF[2] : -1.83697e-16 [ 1.00000 0.00000 ] 1.00000 [ 0.00000 2.44929e-16 ]
+DEAL::ddF[0] : -1.00000 0.00000 0.00000 -1.22465e-16
+DEAL::ddF[1] : -1.22465e-16 0.00000 0.00000 1.00000
+DEAL::ddF[2] : 1.00000 0.00000 0.00000 2.44929e-16
+DEAL::dim = 3, spacedim = 3
+DEAL::x  : 1.57080 [ 1.00000 0.00000 0.00000 ] [ 1.00000 [ ] 0.00000 [ ] 0.00000 [ ] ] 3.14159 [ 0.00000 1.00000 0.00000 ] [ 0.00000 [ ] 1.00000 [ ] 0.00000 [ ] ] 4.71239 [ 0.00000 0.00000 1.00000 ] [ 0.00000 [ ] 0.00000 [ ] 1.00000 [ ] ]
+DEAL::F[0] : 2.22045e-16 [ 6.12323e-17 -1.00000 -1.83697e-16 ] [ 6.12323e-17 [ -1.00000 0.00000 0.00000 ] -1.00000 [ 0.00000 -1.22465e-16 0.00000 ] -1.83697e-16 [ 0.00000 0.00000 1.00000 ] ]
+DEAL::F[1] : -1.00000 [ -1.00000 -1.83697e-16 1.00000 ] [ -1.00000 [ -1.22465e-16 0.00000 0.00000 ] -1.83697e-16 [ 0.00000 1.00000 0.00000 ] 1.00000 [ 0.00000 0.00000 2.44929e-16 ] ]
+DEAL::F[2] : -2.22045e-16 [ -1.83697e-16 1.00000 3.06162e-16 ] [ -1.83697e-16 [ 1.00000 0.00000 0.00000 ] 1.00000 [ 0.00000 2.44929e-16 0.00000 ] 3.06162e-16 [ 0.00000 0.00000 -1.00000 ] ]
+DEAL::dF[0] : 6.12323e-17 [ -1.00000 0.00000 0.00000 ] -1.00000 [ 0.00000 -1.22465e-16 0.00000 ] -1.83697e-16 [ 0.00000 0.00000 1.00000 ]
+DEAL::dF[1] : -1.00000 [ -1.22465e-16 0.00000 0.00000 ] -1.83697e-16 [ 0.00000 1.00000 0.00000 ] 1.00000 [ 0.00000 0.00000 2.44929e-16 ]
+DEAL::dF[2] : -1.83697e-16 [ 1.00000 0.00000 0.00000 ] 1.00000 [ 0.00000 2.44929e-16 0.00000 ] 3.06162e-16 [ 0.00000 0.00000 -1.00000 ]
+DEAL::ddF[0] : -1.00000 0.00000 0.00000 0.00000 -1.22465e-16 0.00000 0.00000 0.00000 1.00000
+DEAL::ddF[1] : -1.22465e-16 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 2.44929e-16
+DEAL::ddF[2] : 1.00000 0.00000 0.00000 0.00000 2.44929e-16 0.00000 0.00000 0.00000 -1.00000

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.