]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Tensor<rank,dim,Number> - Implement operator* with TensorAccessors::contract
authorMatthias Maier <tamiko@43-1.org>
Sat, 12 Sep 2015 02:36:59 +0000 (21:36 -0500)
committerMatthias Maier <tamiko@43-1.org>
Sat, 12 Sep 2015 20:55:14 +0000 (15:55 -0500)
include/deal.II/base/tensor.h

index 7d231bc49d609447b8649ddeed29b1b12ace1ce1..b86ea17fb9220adb902029cc5a9689f2809a6ade 100644 (file)
@@ -244,13 +244,13 @@ public:
   template <class Archive>
   void serialize(Archive &ar, const unsigned int version);
 
-private:
   /**
    * Internal type declaration that is used to specialize the return type
    * of operator[]() for Tensor<1,dim,Number>
    */
   typedef Number tensor_type;
 
+private:
   /**
    * The value of this scalar object.
    */
@@ -557,13 +557,13 @@ public:
                   << arg1
                   << ", but this is not possible for tensors of the current type.");
 
-private:
   /**
    * Internal type declaration that is used to specialize the return type
    * of operator[]() for Tensor<1,dim,Number>
    */
   typedef Tensor<rank_, dim, Number> tensor_type;
 
+private:
   /**
    * Array of tensors holding the subelements.
    */
@@ -1164,21 +1164,6 @@ Tensor<rank_,dim,Number>::serialize(Archive &ar, const unsigned int)
 
 /* ----------------- Non-member functions operating on tensors. ------------ */
 
-
-#ifndef DEAL_II_WITH_CXX11
-template <typename T, typename U, int rank, int dim>
-struct ProductType<T,Tensor<rank,dim,U> >
-{
-  typedef Tensor<rank,dim,typename ProductType<T,U>::type> type;
-};
-
-template <typename T, typename U, int rank, int dim>
-struct ProductType<Tensor<rank,dim,T>,U>
-{
-  typedef Tensor<rank,dim,typename ProductType<T,U>::type> type;
-};
-#endif
-
 /**
  * @name Output functions for Tensor objects
  */
@@ -1227,35 +1212,80 @@ std::ostream &operator << (std::ostream &out, const Tensor<0,dim,Number> &p)
  */
 //@{
 
+
+#ifndef DEAL_II_WITH_CXX11
+template <typename T, typename U, int rank, int dim>
+struct ProductType<T,Tensor<rank,dim,U> >
+{
+  typedef Tensor<rank,dim,typename ProductType<T,U>::type> type;
+};
+
+template <typename T, typename U, int rank, int dim>
+struct ProductType<Tensor<rank,dim,T>,U>
+{
+  typedef Tensor<rank,dim,typename ProductType<T,U>::type> type;
+};
+#endif
+
+
+
 /**
- * Scalar multiplication of a tensor of rank 0 with a scalar from the left.
+ * Scalar multiplication of a tensor of rank 0 with an object from the
+ * left.
+ *
+ * This function unwraps the underlying @p Number stored in the Tensor and
+ * multiplies @p object with it.
  *
  * @relates Tensor<0,dim,Number>
  * @relates EnableIfScalar
  */
-template <int dim, typename Number, typename OtherNumber>
+template <int dim, typename Number, typename Other>
 inline
-Tensor<0,dim,typename ProductType<typename EnableIfScalar<OtherNumber>::type, Number>::type>
-operator * (const OtherNumber           factor,
+typename ProductType<Other, Number>::type
+operator * (const Other                 object,
             const Tensor<0,dim,Number> &t)
 {
-  return factor * static_cast<const Number &>(t);
+  return object * static_cast<const Number &>(t);
 }
 
 
 /**
- * Scalar multiplication of a tensor of rank 0 with a scalar from the right.
+ * Scalar multiplication of a tensor of rank 0 with an object from the
+ * right.
+ *
+ * This function unwraps the underlying @p Number stored in the Tensor and
+ * multiplies @p object with it.
  *
  * @relates Tensor<0,dim,Number>
  * @relates EnableIfScalar
  */
-template <int dim, typename Number, typename OtherNumber>
+template <int dim, typename Number, typename Other>
 inline
-Tensor<0,dim,typename ProductType<Number, typename EnableIfScalar<OtherNumber>::type>::type>
+typename ProductType<Number, Other>::type
 operator * (const Tensor<0,dim,Number> &t,
-            const OtherNumber           factor)
+            const Other                 object)
 {
-  return static_cast<const Number &>(t) * factor;
+  return static_cast<const Number &>(t) * object;
+}
+
+
+/**
+ * Scalar multiplication of two tensors of rank 0.
+ *
+ * This function unwraps the underlying objects of type @p Number and @p
+ * OtherNumber that are stored within the Tensor and multiplies them.
+ * It returns an unwrapped number of product type.
+ *
+ * @relates Tensor<0,dim,Number>
+ */
+template <int dim, typename Number, typename OtherNumber>
+inline
+typename ProductType<Number, OtherNumber>::type // FIXME: TEST!
+operator * (const Tensor<0, dim, Number>      &src1,
+            const Tensor<0, dim, OtherNumber> &src2)
+{
+  return static_cast<const Number &>(src1) *
+         static_cast<const OtherNumber &>(src2);
 }
 
 
@@ -1307,24 +1337,9 @@ operator- (const Tensor<0,dim,Number> &p, const Tensor<0,dim,OtherNumber> &q)
  * Multiplication of a tensor of general rank with a scalar number from the
  * right.
  *
- * The purpose of this operator is to enable only multiplication of a tensor
- * by a scalar number (i.e., a floating point number, a complex floating point
- * number, etc.). The function is written in a way that only allows the
- * compiler to consider the function if the second argument is indeed a scalar
- * number -- in other words, @p OtherNumber will not match, for example
- * <code>std::vector@<double@></code> as the product of a tensor and a vector
- * clearly would make no sense. The mechanism by which the compiler is
- * prohibited of considering this operator for multiplication with non-scalar
- * types are explained in the documentation of the EnableIfScalar class.
- *
- * The return type of the function is chosen so that it matches the types of
- * both the tensor and the scalar argument. For example, if you multiply a
- * <code>Tensor@<1,dim,double@></code> by <code>std::complex@<double@></code>,
- * then the result will be a
- * <code>Tensor@<1,dim,std::complex@<double@>@></code>. In other words, the
- * type with which the returned tensor stores its components equals the type
- * you would get if you multiplied an individual component of the input tensor
- * by the scalar factor.
+ * Only multiplication with a scalar number type (i.e., a floating point
+ * number, a complex floating point number, etc.), see the documentation of
+ * EnableIfScalar for details.
  *
  * @relates Tensor
  * @relates EnableIfScalar
@@ -1450,20 +1465,54 @@ operator- (const Tensor<rank,dim,Number> &p, const Tensor<rank,dim,OtherNumber>
  */
 //@{
 
+
 /**
- * Returns the contraction of two Tensors of rank 0.
+ * The dot product (single contraction) for tensors: Return a tensor of
+ * rank $(\text{rank\_1} + \text{rank\_2} - 2)$ that is the contraction of
+ * the last index of a tensor @p src1 of rank @p rank_1 with the first
+ * index of a tensor @p src2 of rank @p rank_2:
+ * @f[
+ *   \text{result}_{i_1,..,i_{r1},j_1,..,j_{r2}}
+ *   = \sum_{k}
+ *     \text{left}_{i_1,..,i_{r1}, k}
+ *     \text{right}_{j_1,..,j_{r2}, k}
+ * @f]
+ *
+ * @note For the <tt>Tensor</tt> class, the multiplication operator only
+ * performs a contraction over a single pair of indices. This is in
+ * contrast to the multiplication operator for symmetric tensors, which
+ * does the double contraction.
+ *
+ * @note In case the contraction yields tensor of rank 0 the scalar
+ * number is returned as an unwrapped number type
  *
- * @relates Tensor<0,dim,Number>
+ * @relates Tensor
  */
-template <int dim, typename Number, typename OtherNumber>
+template <int rank_1, int rank_2, int dim,
+          typename Number, typename OtherNumber,
+          typename = typename std::enable_if<rank_1 != 0>::type,
+          typename = typename std::enable_if<rank_2 != 0>::type>
+
 inline
-typename ProductType<Number, OtherNumber>::type
-operator* (const Tensor<0,dim,Number> &p, const Tensor<0,dim,OtherNumber> &q)
+typename Tensor<rank_1 + rank_2 - 2, dim, typename ProductType<Number, OtherNumber>::type>::tensor_type
+operator * (const Tensor<rank_1, dim, Number> &src1,
+            const Tensor<rank_2, dim, OtherNumber> &src2)
 {
-  return static_cast<const Number &>(p) * static_cast<const OtherNumber &>(q);
+  typename Tensor<rank_1 + rank_2 - 2, dim, typename ProductType<Number, OtherNumber>::type>::tensor_type result;
+
+  TensorAccessors::internal::ReorderedIndexView<0, rank_2, const Tensor<rank_2, dim, OtherNumber> >
+  reordered = TensorAccessors::reordered_index_view<0, rank_2>(src2);
+  TensorAccessors::contract<1, rank_1, rank_2, dim>(result, src1, reordered);
+
+  return result;
 }
 
+
 //@}
+/**
+ * @name To be refactored
+ */
+//@{
 
 
 /**
@@ -1487,30 +1536,6 @@ contract (const Tensor<1,dim,Number> &src1,
 }
 
 
-/**
- * Multiplication operator performing a contraction of the last index of the
- * first argument and the first index of the second argument. This function
- * therefore does the same as the corresponding <tt>contract</tt> function,
- * but returns the result as a return value, rather than writing it into the
- * reference given as the first argument to the <tt>contract</tt> function.
- *
- * Note that for the <tt>Tensor</tt> class, the multiplication operator only
- * performs a contraction over a single pair of indices. This is in contrast
- * to the multiplication operator for symmetric tensors, which does the double
- * contraction.
- *
- * @relates Tensor
- */
-template <int dim, typename Number, typename OtherNumber>
-inline
-typename ProductType<Number,OtherNumber>::type
-operator * (const Tensor<1,dim,Number> &src1,
-            const Tensor<1,dim,OtherNumber> &src2)
-{
-  return contract(src1, src2);
-}
-
-
 /**
  * Double contract two tensors of rank 2, thus computing the Frobenius inner
  * product <tt> sum<sub>i,j</sub> src1[i][j]*src2[i][j]</tt>.
@@ -1553,37 +1578,6 @@ void contract (Tensor<1,dim,Number>       &dest,
 }
 
 
-/**
- * Multiplication operator performing a contraction of the last index of the
- * first argument and the first index of the second argument. This function
- * therefore does the same as the corresponding <tt>contract</tt> function,
- * but returns the result as a return value, rather than writing it into the
- * reference given as the first argument to the <tt>contract</tt> function.
- *
- * Note that for the <tt>Tensor</tt> class, the multiplication operator only
- * performs a contraction over a single pair of indices. This is in contrast
- * to the multiplication operator for symmetric tensors, which does the double
- * contraction.
- *
- * @relates Tensor
- * @author Wolfgang Bangerth, 2005
- */
-template <int dim, typename Number>
-Tensor<1,dim,Number>
-operator * (const Tensor<2,dim,Number> &src1,
-            const Tensor<1,dim,Number> &src2)
-{
-  Tensor<1,dim,Number> dest;
-  for (unsigned int i=0; i<dim; ++i)
-    {
-      dest[i] = src1[i][0] * src2[0];
-      for (unsigned int j=1; j<dim; ++j)
-        dest[i] += src1[i][j] * src2[j];
-    }
-  return dest;
-}
-
-
 /**
  * Contract a tensor of rank 1 with a tensor of rank 2. The result is
  * <tt>dest[i] = sum_j src1[j] src2[j][i]</tt>.
@@ -1606,38 +1600,6 @@ void contract (Tensor<1,dim,Number>       &dest,
 }
 
 
-/**
- * Multiplication operator performing a contraction of the last index of the
- * first argument and the first index of the second argument. This function
- * therefore does the same as the corresponding <tt>contract</tt> function,
- * but returns the result as a return value, rather than writing it into the
- * reference given as the first argument to the <tt>contract</tt> function.
- *
- * Note that for the <tt>Tensor</tt> class, the multiplication operator only
- * performs a contraction over a single pair of indices. This is in contrast
- * to the multiplication operator for symmetric tensors, which does the double
- * contraction.
- *
- * @relates Tensor
- * @author Wolfgang Bangerth, 2005
- */
-template <int dim, typename Number>
-inline
-Tensor<1,dim,Number>
-operator * (const Tensor<1,dim,Number> &src1,
-            const Tensor<2,dim,Number> &src2)
-{
-  Tensor<1,dim,Number> dest;
-  for (unsigned int i=0; i<dim; ++i)
-    {
-      dest[i] = src1[0] * src2[0][i];
-      for (unsigned int j=1; j<dim; ++j)
-        dest[i] += src1[j] * src2[j][i];
-    }
-  return dest;
-}
-
-
 /**
  * Contract a tensor of rank 2 with a tensor of rank 2. The result is
  * <tt>dest[i][k] = sum_j src1[i][j] src2[j][k]</tt>.
@@ -1661,37 +1623,6 @@ void contract (Tensor<2,dim,Number>       &dest,
 }
 
 
-
-/**
- * Multiplication operator performing a contraction of the last index of the
- * first argument and the first index of the second argument. This function
- * therefore does the same as the corresponding <tt>contract</tt> function,
- * but returns the result as a return value, rather than writing it into the
- * reference given as the first argument to the <tt>contract</tt> function.
- *
- * Note that for the <tt>Tensor</tt> class, the multiplication operator only
- * performs a contraction over a single pair of indices. This is in contrast
- * to the multiplication operator for symmetric tensors, which does the double
- * contraction.
- *
- * @relates Tensor
- * @author Wolfgang Bangerth, 2005
- */
-template <int dim, typename Number>
-inline
-Tensor<2,dim,Number>
-operator * (const Tensor<2,dim,Number> &src1,
-            const Tensor<2,dim,Number> &src2)
-{
-  Tensor<2,dim,Number> dest;
-  for (unsigned int i=0; i<dim; ++i)
-    for (unsigned int j=0; j<dim; ++j)
-      for (unsigned int k=0; k<dim; ++k)
-        dest[i][j] += src1[i][k] * src2[k][j];
-  return dest;
-}
-
-
 /**
  * Contract a tensor of rank 2 with a tensor of rank 2. The contraction is
  * performed over index <tt>index1</tt> of the first tensor, and
@@ -1934,37 +1865,6 @@ void contract (Tensor<3,dim,Number>       &dest,
 }
 
 
-/**
- * Multiplication operator performing a contraction of the last index of the
- * first argument and the first index of the second argument. This function
- * therefore does the same as the corresponding <tt>contract</tt> function,
- * but returns the result as a return value, rather than writing it into the
- * reference given as the first argument to the <tt>contract</tt> function.
- *
- * Note that for the <tt>Tensor</tt> class, the multiplication operator only
- * performs a contraction over a single pair of indices. This is in contrast
- * to the multiplication operator for symmetric tensors, which does the double
- * contraction.
- *
- * @relates Tensor
- * @author Wolfgang Bangerth, 2005
- */
-template <int dim, typename Number>
-inline
-Tensor<3,dim,Number>
-operator * (const Tensor<3,dim,Number> &src1,
-            const Tensor<2,dim,Number> &src2)
-{
-  Tensor<3,dim,Number> dest;
-  for (unsigned int i=0; i<dim; ++i)
-    for (unsigned int j=0; j<dim; ++j)
-      for (unsigned int k=0; k<dim; ++k)
-        for (unsigned int l=0; l<dim; ++l)
-          dest[i][j][k] += src1[i][j][l] * src2[l][k];
-  return dest;
-}
-
-
 /**
  * Contract a tensor of rank 2 with a tensor of rank 3. The result is
  * <tt>dest[i][j][l] = sum_k src1[i][k] src2[k][j][l]</tt>.
@@ -1987,61 +1887,6 @@ void contract (Tensor<3,dim,Number>       &dest,
 }
 
 
-/**
- * Multiplication operator performing a contraction of the last index of the
- * first argument and the first index of the second argument. This function
- * therefore does the same as the corresponding <tt>contract</tt> function,
- * but returns the result as a return value, rather than writing it into the
- * reference given as the first argument to the <tt>contract</tt> function.
- *
- * Note that for the <tt>Tensor</tt> class, the multiplication operator only
- * performs a contraction over a single pair of indices. This is in contrast
- * to the multiplication operator for symmetric tensors, which does the double
- * contraction.
- *
- * @relates Tensor
- * @author Wolfgang Bangerth, 2005
- */
-template <int dim, typename Number>
-inline
-Tensor<3,dim,Number>
-operator * (const Tensor<2,dim,Number> &src1,
-            const Tensor<3,dim,Number> &src2)
-{
-  Tensor<3,dim,Number> dest;
-  for (unsigned int i=0; i<dim; ++i)
-    for (unsigned int j=0; j<dim; ++j)
-      for (unsigned int k=0; k<dim; ++k)
-        for (unsigned int l=0; l<dim; ++l)
-          dest[i][j][k] += src1[i][l] * src2[l][j][k];
-  return dest;
-}
-
-
-/**
- * Contract a tensor of rank 3 with a tensor of rank 3. The result is
- * <tt>dest[i][j][k][l] = sum_m src1[i][j][m] src2[m][k][l]</tt>.
- *
- * @relates Tensor
- * @author Wolfgang Bangerth, 1998
- */
-template <int dim, typename Number>
-inline
-Tensor<4,dim,Number>
-operator * (const Tensor<3,dim,Number> &src1,
-            const Tensor<3,dim,Number> &src2)
-{
-  Tensor<4,dim,Number> dest;
-  for (unsigned int i=0; i<dim; ++i)
-    for (unsigned int j=0; j<dim; ++j)
-      for (unsigned int k=0; k<dim; ++k)
-        for (unsigned int l=0; l<dim; ++l)
-          for (unsigned int m=0; m<dim; ++m)
-            dest[i][j][k][l] += src1[i][j][m] * src2[m][k][l];
-  return dest;
-}
-
-
 /**
  * Contract the last two indices of <tt>src1</tt> with the two indices
  * <tt>src2</tt>, creating a rank-2 tensor. This is the matrix-vector product
@@ -2233,7 +2078,6 @@ void outer_product (Tensor<1,dim,Number>       &dst,
 }
 
 
-
 /**
  * Form the outer product of two tensors of rank 1 and 0, i.e. <tt>dst[i] =
  * src1[i] * src2</tt>. Of course, this is only a scaling of <tt>src1</tt>,
@@ -2254,55 +2098,6 @@ void outer_product (Tensor<1,dim,Number>       &dst,
 }
 
 
-/**
- * Cross-product in 2d. This is just a rotation by 90 degrees clockwise to
- * compute the outer normal from a tangential vector. This function is defined
- * for all space dimensions to allow for dimension independent programming
- * (e.g. within switches over the space dimension), but may only be called if
- * the actual dimension of the arguments is two (e.g. from the <tt>dim==2</tt>
- * case in the switch).
- *
- * @relates Tensor
- * @author Guido Kanschat, 2001
- */
-template <int dim, typename Number>
-inline
-void
-cross_product (Tensor<1,dim,Number>       &dst,
-               const Tensor<1,dim,Number> &src)
-{
-  Assert (dim==2, ExcInternalError());
-
-  dst[0] = src[1];
-  dst[1] = -src[0];
-}
-
-
-/**
- * Cross-product of 2 vectors in 3d. This function is defined for all space
- * dimensions to allow for dimension independent programming (e.g. within
- * switches over the space dimension), but may only be called if the actual
- * dimension of the arguments is three (e.g. from the <tt>dim==3</tt> case in
- * the switch).
- *
- * @relates Tensor
- * @author Guido Kanschat, 2001
- */
-template <int dim, typename Number>
-inline
-void
-cross_product (Tensor<1,dim,Number>       &dst,
-               const Tensor<1,dim,Number> &src1,
-               const Tensor<1,dim,Number> &src2)
-{
-  Assert (dim==3, ExcInternalError());
-
-  dst[0] = src1[1]*src2[2] - src1[2]*src2[1];
-  dst[1] = src1[2]*src2[0] - src1[0]*src2[2];
-  dst[2] = src1[0]*src2[1] - src1[1]*src2[0];
-}
-
-
 /**
  * Compute the scalar product $a:b=\sum_{i,j} a_{ij}b_{ij}$ between two
  * tensors $a,b$ of rank 2. We don't use <code>operator*</code> for this
@@ -2342,7 +2137,6 @@ Number determinant (const Tensor<rank,1,Number> &t)
 }
 
 
-
 /**
  * Compute the determinant of a tensor of rank one and dimension one. Since
  * this is a number, the return value is, of course, the number itself.
@@ -2358,6 +2152,62 @@ Number determinant (const Tensor<1,1,Number> &t)
 }
 
 
+//@}
+/**
+ * @name Special operations on tensors of rank 1
+ */
+//@{
+
+
+/**
+ * Cross-product in 2d. This is just a rotation by 90 degrees clockwise to
+ * compute the outer normal from a tangential vector. This function is defined
+ * for all space dimensions to allow for dimension independent programming
+ * (e.g. within switches over the space dimension), but may only be called if
+ * the actual dimension of the arguments is two (e.g. from the <tt>dim==2</tt>
+ * case in the switch).
+ *
+ * @relates Tensor
+ * @author Guido Kanschat, 2001
+ */
+template <int dim, typename Number>
+inline
+void
+cross_product (Tensor<1,dim,Number>       &dst,
+               const Tensor<1,dim,Number> &src)
+{
+  Assert (dim==2, ExcInternalError());
+
+  dst[0] = src[1];
+  dst[1] = -src[0];
+}
+
+
+/**
+ * Cross-product of 2 vectors in 3d. This function is defined for all space
+ * dimensions to allow for dimension independent programming (e.g. within
+ * switches over the space dimension), but may only be called if the actual
+ * dimension of the arguments is three (e.g. from the <tt>dim==3</tt> case in
+ * the switch).
+ *
+ * @relates Tensor
+ * @author Guido Kanschat, 2001
+ */
+template <int dim, typename Number>
+inline
+void
+cross_product (Tensor<1,dim,Number>       &dst,
+               const Tensor<1,dim,Number> &src1,
+               const Tensor<1,dim,Number> &src2)
+{
+  Assert (dim==3, ExcInternalError());
+
+  dst[0] = src1[1]*src2[2] - src1[2]*src2[1];
+  dst[1] = src1[2]*src2[0] - src1[0]*src2[2];
+  dst[2] = src1[0]*src2[1] - src1[1]*src2[0];
+}
+
+
 /**
  * Compute the determinant of a tensor of rank two and dimension one. Since
  * this is a number, the return value is, of course, the number itself.
@@ -2373,7 +2223,6 @@ Number determinant (const Tensor<2,1,Number> &t)
 }
 
 
-
 /**
  * Compute the determinant of a tensor or rank 2, here for <tt>dim==2</tt>.
  *
@@ -2453,7 +2302,6 @@ Number determinant (const Tensor<2,dim,Number> &t)
 }
 
 
-
 /**
  * Compute and return the trace of a tensor of rank 2, i.e. the sum of its
  * diagonal entries.
@@ -2471,7 +2319,6 @@ Number trace (const Tensor<2,dim,Number> &d)
 }
 
 
-
 /**
  * Compute and return the inverse of the given tensor. Since the compiler can
  * perform the return value optimization, and since the size of the return
@@ -2540,7 +2387,6 @@ invert (const Tensor<2,dim,Number> &t)
 }
 
 
-
 /**
  * Return the transpose of the given tensor. Since the compiler can perform
  * the return value optimization, and since the size of the return object is
@@ -2587,8 +2433,6 @@ transpose (const Tensor<2,1,Number> &t)
 }
 
 
-
-
 /**
  * Return the transpose of the given tensor. This is the specialization of the
  * general template for <tt>dim==2</tt>.
@@ -2606,8 +2450,6 @@ transpose (const Tensor<2,2,Number> &t)
 }
 
 
-
-
 /**
  * Return the transpose of the given tensor. This is the specialization of the
  * general template for <tt>dim==3</tt>.
@@ -2657,7 +2499,6 @@ l1_norm (const Tensor<2,dim,Number> &t)
 }
 
 
-
 /**
  * Return the $l_\infty$ norm of the given rank-2 tensor, where $||t||_\infty
  * = \max_i \sum_j |t_{ij}|$ (maximum of the sums over rows).
@@ -2684,7 +2525,7 @@ linfty_norm (const Tensor<2,dim,Number> &t)
   return max;
 }
 
-
+//@}
 
 DEAL_II_NAMESPACE_CLOSE
 

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.