]> https://gitweb.dealii.org/ - dealii.git/commitdiff
instantiated get_function_third_derivatives template functions 1696/head
authorMaien Hamed <tomaien@hotmail.com>
Wed, 30 Sep 2015 20:45:56 +0000 (22:45 +0200)
committerMaien Hamed <tomaien@hotmail.com>
Wed, 30 Sep 2015 20:56:07 +0000 (22:56 +0200)
source/fe/fe_values.inst.in
tests/fe/fe_values_view_30.cc [new file with mode: 0644]
tests/fe/fe_values_view_30.output [new file with mode: 0644]

index 8678b601a1f6ba10a4b48acf0422aac7cd01de07..e4842ec21e1a7c59188ca7267ca2f15618ce6408 100644 (file)
@@ -75,6 +75,10 @@ for (VEC : SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimensi
       void FEValuesViews::Scalar<deal_II_dimension, deal_II_space_dimension>
       ::get_function_laplacians<dealii::VEC>
       (const dealii::VEC&, std::vector<ProductType<dealii::VEC::value_type,double>::type> &) const;
+    template
+      void FEValuesViews::Scalar<deal_II_dimension, deal_II_space_dimension>
+      ::get_function_third_derivatives<dealii::VEC>
+      (const dealii::VEC&, std::vector<ProductType<dealii::VEC::value_type,dealii::Tensor<3,deal_II_space_dimension> >::type>&) const;
 
     template
       void FEValuesViews::Vector<deal_II_dimension, deal_II_space_dimension>
@@ -104,6 +108,10 @@ for (VEC : SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimensi
       void FEValuesViews::Vector<deal_II_dimension, deal_II_space_dimension>
       ::get_function_laplacians<dealii::VEC>
       (const dealii::VEC&, std::vector<ProductType<dealii::VEC::value_type,dealii::Tensor<1,deal_II_space_dimension> >::type >&) const;
+    template
+      void FEValuesViews::Vector<deal_II_dimension, deal_II_space_dimension>
+      ::get_function_third_derivatives<dealii::VEC>
+      (const dealii::VEC&, std::vector<ProductType<dealii::VEC::value_type,dealii::Tensor<4,deal_II_space_dimension> >::type >&) const;
 
     template
       void FEValuesViews::SymmetricTensor<2, deal_II_dimension, deal_II_space_dimension>
@@ -204,6 +212,22 @@ for (VEC : SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimensi
       void FEValuesBase<deal_II_dimension,deal_II_space_dimension>::get_function_laplacians<VEC>
       (const VEC&, const VectorSlice<const std::vector<types::global_dof_index> >&,
        std::vector<std::vector<VEC::value_type> > &, bool) const;
+       
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_space_dimension>::get_function_third_derivatives<VEC>
+      (const VEC&, std::vector<dealii::Tensor<3,deal_II_space_dimension,VEC::value_type> > &) const;
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_space_dimension>::get_function_third_derivatives<VEC>
+      (const VEC&, const VectorSlice<const std::vector<types::global_dof_index> >&,
+       std::vector<dealii::Tensor<3,deal_II_space_dimension,VEC::value_type> > &) const;
+
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_space_dimension>::get_function_third_derivatives<VEC>
+      (const VEC&, std::vector<std::vector<dealii::Tensor<3,deal_II_space_dimension,VEC::value_type> > > &, bool) const;
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_space_dimension>::get_function_third_derivatives<VEC>
+      (const VEC&, const VectorSlice<const std::vector<types::global_dof_index> >&,
+       VectorSlice<std::vector<std::vector<dealii::Tensor<3,deal_II_space_dimension,VEC::value_type> > > >, bool) const;
 
 #endif
   }
@@ -225,6 +249,9 @@ for (deal_II_dimension : DIMENSIONS)
     template
         void FEValuesViews::Scalar<deal_II_dimension>::get_function_laplacians<dealii::IndexSet>
         (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,double>::type>&) const;
+    template
+        void FEValuesViews::Scalar<deal_II_dimension>::get_function_third_derivatives<dealii::IndexSet>
+        (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,dealii::Tensor<3,deal_II_dimension> >::type>&) const;
 
     template
         void FEValuesViews::Vector<deal_II_dimension>::get_function_values<dealii::IndexSet>
@@ -247,6 +274,9 @@ for (deal_II_dimension : DIMENSIONS)
     template
         void FEValuesViews::Vector<deal_II_dimension>::get_function_laplacians<dealii::IndexSet>
         (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,dealii::Tensor<1,deal_II_dimension> >::type>&) const;
+    template
+        void FEValuesViews::Vector<deal_II_dimension>::get_function_third_derivatives<dealii::IndexSet>
+        (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,dealii::Tensor<4,deal_II_dimension> >::type>&) const;
 
     template
         void FEValuesViews::SymmetricTensor<2,deal_II_dimension,deal_II_dimension>::get_function_values<dealii::IndexSet>
@@ -280,6 +310,10 @@ for (deal_II_dimension : DIMENSIONS)
     void FEValuesViews::Scalar<deal_II_dimension, deal_II_dimension+1>
     ::get_function_laplacians<dealii::IndexSet>
     (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,double>::type>&) const;
+    template
+    void FEValuesViews::Scalar<deal_II_dimension, deal_II_dimension+1>
+    ::get_function_third_derivatives<dealii::IndexSet>
+    (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,dealii::Tensor<3,deal_II_dimension+1> >::type>&) const;
 
     template
     void FEValuesViews::Vector<deal_II_dimension, deal_II_dimension+1>
@@ -305,6 +339,10 @@ for (deal_II_dimension : DIMENSIONS)
     void FEValuesViews::Vector<deal_II_dimension, deal_II_dimension+1>
     ::get_function_laplacians<dealii::IndexSet>
     (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,dealii::Tensor<1,deal_II_dimension+1> >::type>&) const;
+    template
+    void FEValuesViews::Vector<deal_II_dimension, deal_II_dimension+1>
+    ::get_function_third_derivatives<dealii::IndexSet>
+    (const dealii::IndexSet&, std::vector<ProductType<IndexSet::value_type,dealii::Tensor<4,deal_II_dimension+1> >::type>&) const;
 
     template
     void FEValuesViews::SymmetricTensor<2, deal_II_dimension, deal_II_dimension+1>
@@ -407,6 +445,22 @@ for (deal_II_dimension : DIMENSIONS)
       (const IndexSet&, const VectorSlice<const std::vector<types::global_dof_index> >&,
        std::vector<std::vector<IndexSet::value_type> > &, bool) const;
 
+    template
+      void FEValuesBase<deal_II_dimension>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, std::vector<dealii::Tensor<3,deal_II_dimension,IndexSet::value_type> > &) const;
+    template
+      void FEValuesBase<deal_II_dimension>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, const VectorSlice<const std::vector<types::global_dof_index> >&,
+       std::vector<dealii::Tensor<3,deal_II_dimension,IndexSet::value_type> > &) const;
+
+    template
+      void FEValuesBase<deal_II_dimension>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, std::vector<std::vector<dealii::Tensor<3,deal_II_dimension,IndexSet::value_type> > > &, bool) const;
+    template
+      void FEValuesBase<deal_II_dimension>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, const VectorSlice<const std::vector<types::global_dof_index> >&,
+       VectorSlice<std::vector<std::vector<dealii::Tensor<3,deal_II_dimension,IndexSet::value_type> > > >, bool) const;
+
 
 #if deal_II_dimension != 3
 
@@ -484,5 +538,21 @@ for (deal_II_dimension : DIMENSIONS)
       (const IndexSet&, const VectorSlice<const std::vector<types::global_dof_index> >&,
        std::vector<std::vector<IndexSet::value_type> > &, bool) const;
 
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_dimension+1>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, std::vector<dealii::Tensor<3,deal_II_dimension+1,IndexSet::value_type> > &) const;
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_dimension+1>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, const VectorSlice<const std::vector<types::global_dof_index> >&,
+       std::vector<dealii::Tensor<3,deal_II_dimension+1,IndexSet::value_type> > &) const;
+
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_dimension+1>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, std::vector<std::vector<dealii::Tensor<3,deal_II_dimension+1,IndexSet::value_type> > > &, bool) const;
+    template
+      void FEValuesBase<deal_II_dimension,deal_II_dimension+1>::get_function_third_derivatives<IndexSet>
+      (const IndexSet&, const VectorSlice<const std::vector<types::global_dof_index> >&,
+       VectorSlice<std::vector<std::vector<dealii::Tensor<3,deal_II_dimension+1,IndexSet::value_type> > > >, bool) const;
+
 #endif
   }
diff --git a/tests/fe/fe_values_view_30.cc b/tests/fe/fe_values_view_30.cc
new file mode 100644 (file)
index 0000000..8469b16
--- /dev/null
@@ -0,0 +1,188 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2007 - 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// test the FEValues views and extractor classes. this test is for
+// get_function_hessians for vector components and a non-primitive element
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/function.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/dofs/dof_tools.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_nedelec.h>
+#include <deal.II/fe/fe_raviart_thomas.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/fe/mapping_q1.h>
+
+#include <fstream>
+
+template<int dim>
+class VectorFunction : public Function<dim>
+{
+public:
+  VectorFunction() : Function<dim>(dim) {}
+  virtual double value (const Point<dim> &p, const unsigned int component) const;
+  virtual void vector_value(const Point<dim> &p, Vector<double> &values) const;
+};
+
+template<>
+double VectorFunction<2>::value(const Point<2> &p, const unsigned int component) const
+{
+  Assert (component < 2,  ExcIndexRange (component, 0, 1));
+
+  const double PI = numbers::PI;
+  double val = 0.0;
+  switch (component)
+    {
+    case 0:
+      val = pow(p(0),3);
+      break;
+    case 1:
+      val = pow(p(1),2)*p(0);
+      break;
+    }
+  return val;
+}
+
+template<>
+double VectorFunction<3>::value(const Point<3> &p, const unsigned int component) const
+{
+  Assert (component < 3, ExcIndexRange (component, 0, 2));
+
+  const double PI = numbers::PI;
+  double val = 0.0;
+  switch (component)
+    {
+    case 0:
+      val = pow(p(0),3);
+      break;
+    case 1:
+      val = pow(p(1),2)*p(0);
+      break;
+    case 2:
+      val = p(2)*p(1)*p(0);
+      break;
+    }
+  return val;
+}
+
+template<int dim>
+void VectorFunction<dim>::vector_value(const Point<dim> &p, Vector<double> &values) const
+{
+  for (int i = 0; i < dim; ++i)
+    values(i) = value(p, i);
+}
+
+template<int dim>
+void test (const Triangulation<dim> &tr,
+           const FiniteElement<dim> &fe)
+{
+  deallog << "FE=" << fe.get_name()
+          << std::endl;
+
+  DoFHandler<dim> dof_handler(tr);
+  dof_handler.distribute_dofs(fe);
+
+  VectorFunction<dim> fe_function;
+
+  MappingQ1<dim> mapping;
+
+  const QGauss<dim> quadrature(2);
+  FEValues<dim> fe_values (mapping, fe, quadrature,
+                           update_values | update_gradients | update_3rd_derivatives);
+
+  ConstraintMatrix constraints;
+  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
+  constraints.close();
+
+  Vector<double> function_vals(dof_handler.n_dofs());
+  VectorTools::project(mapping, dof_handler, constraints, quadrature, fe_function, function_vals);
+
+  fe_values.reinit (dof_handler.begin_active());
+
+  std::vector<Tensor<4,dim> > selected_vector_values (quadrature.size());
+  std::vector<std::vector<Tensor<3,dim> > >
+  vector_values (quadrature.size(),
+                 std::vector<Tensor<3,dim> >(fe.n_components()));
+
+  fe_values.get_function_third_derivatives (function_vals, vector_values);
+
+  for (unsigned int c=0; c<fe.n_components(); ++c)
+    // use a vector extractor if there
+    // are sufficiently many components
+    // left after the current component
+    // 'c'
+    if (c+dim <= fe.n_components())
+      {
+        FEValuesExtractors::Vector vector_components (c);
+        fe_values[vector_components].get_function_third_derivatives (function_vals,
+            selected_vector_values);
+        deallog << "component=" << c << std::endl;
+
+        for (unsigned int q=0; q<fe_values.n_quadrature_points; ++q)
+          for (unsigned int d=0; d<dim; ++d)
+            {
+              for (unsigned int e=0; e<dim; ++e)
+                for (unsigned int f=0; f<dim; ++f)
+                  for (unsigned int g=0; g<dim; ++g)
+                    deallog << selected_vector_values[q][d][e][f][g] << (e<dim-1 && f<dim-1 && g<dim-1? ", " : "; ");
+              deallog << std::endl;
+              Assert ((selected_vector_values[q][d] - vector_values[q][c+d]).norm()
+                      <= 1e-12 * selected_vector_values[q][d].norm(),
+                      ExcInternalError());
+            }
+      }
+}
+
+
+
+template<int dim>
+void test_hyper_sphere()
+{
+  Triangulation<dim> tr;
+  GridGenerator::hyper_ball(tr);
+
+  static const HyperBallBoundary<dim> boundary;
+  tr.set_boundary (0, boundary);
+
+
+  const unsigned int order = 3;
+
+  test(tr, FESystem<dim> (FE_Q<dim>(order), dim) );
+}
+
+
+int main()
+{
+  std::ofstream logfile ("output");
+  deallog << std::setprecision (3);
+
+  deallog.attach(logfile);
+  deallog.depth_console (0);
+  deallog.threshold_double(1.e-5);
+
+  test_hyper_sphere<2>();
+  test_hyper_sphere<3>();
+}
diff --git a/tests/fe/fe_values_view_30.output b/tests/fe/fe_values_view_30.output
new file mode 100644 (file)
index 0000000..cf36f91
--- /dev/null
@@ -0,0 +1,37 @@
+
+DEAL::FE=FESystem<2>[FE_Q<2>(3)^2]
+DEAL::component=0
+DEAL::38.5, 319.; 319.; 13.1; 319.; 13.1; 13.1; -605.; 
+DEAL::16.0, 104.; 104.; 62.6; 104.; 62.6; 62.6; 30.3; 
+DEAL::38.5, -352.; -352.; 159.; -352.; 159.; 159.; 412.; 
+DEAL::16.0, -140.; -140.; 213.; -140.; 213.; 213.; 189.; 
+DEAL::43.3, -204.; -204.; 140.; -204.; 140.; 140.; -303.; 
+DEAL::12.3, -83.1; -83.1; 124.; -83.1; 124.; 124.; 112.; 
+DEAL::43.3, 172.; 172.; -185.; 172.; -185.; -185.; -707.; 
+DEAL::12.3, 56.7; 56.7; -186.; 56.7; -186.; -186.; -700.; 
+DEAL::FE=FESystem<3>[FE_Q<3>(3)^3]
+DEAL::component=0
+DEAL::473., 1.51e+03, 1.10e+03; 1.51e+03, 683., 1.14e+03; 1.10e+03; 1.14e+03; 453.; 1.51e+03, 683., 1.14e+03; 683., 40.1, -545.; 1.14e+03; -545.; -606.; 1.10e+03; 1.14e+03; 453.; 1.14e+03; -545.; -606.; 453.; -606.; -34.7; 
+DEAL::162., 595., 537.; 595., 293., 797.; 537.; 797.; 296.; 595., 293., 797.; 293., 20.9, 140.; 797.; 140.; 211.; 537.; 797.; 296.; 797.; 140.; 211.; 296.; 211.; 9.62; 
+DEAL::95.1, 122., 350.; 122., 168., 150.; 350.; 150.; 385.; 122., 168., 150.; 168., 137., 156.; 150.; 156.; 85.1; 350.; 150.; 385.; 150.; 156.; 85.1; 385.; 85.1; 92.1; 
+DEAL::473., -1.57e+03, -1.06e+03; -1.57e+03, 709., 1.85e+03; -1.06e+03; 1.85e+03; 415.; -1.57e+03, 709., 1.85e+03; 709., -26.7, -145.; 1.85e+03; -145.; -208.; -1.06e+03; 1.85e+03; 415.; 1.85e+03; -145.; -208.; 415.; -208.; 28.5; 
+DEAL::162., -743., -521.; -743., 452., 1.38e+03; -521.; 1.38e+03; 275.; -743., 452., 1.38e+03; 452., -18.0, -351.; 1.38e+03; -351.; -459.; -521.; 1.38e+03; 275.; 1.38e+03; -351.; -459.; 275.; -459.; -4.64; 
+DEAL::95.1, -183., -415.; -183., 296., 276.; -415.; 276.; 481.; -183., 296., 276.; 296., -120., -194.; 276.; -194.; -134.; -415.; 276.; 481.; 276.; -194.; -134.; 481.; -134.; -119.; 
+DEAL::428., -1.09e+03, 1.01e+03; -1.09e+03, 206., -999.; 1.01e+03; -999.; 434.; -1.09e+03, 206., -999.; 206., 40.1, 202.; -999.; 202.; 63.2; 1.01e+03; -999.; 434.; -999.; 202.; 63.2; 434.; 63.2; -0.126; 
+DEAL::118., -431., 404.; -431., 145., -684.; 404.; -684.; 264.; -431., 145., -684.; 145., 20.9, 416.; -684.; 416.; -389.; 404.; -684.; 264.; -684.; 416.; -389.; 264.; -389.; -2.49; 
+DEAL::-101., 282., -277.; 282., -230., 257.; -277.; 257.; -308.; 282., -230., 257.; -230., 137., -242.; 257.; -242.; 324.; -277.; 257.; -308.; 257.; -242.; 324.; -308.; 324.; -127.; 
+DEAL::428., 1.17e+03, -982.; 1.17e+03, 318., -1.27e+03; -982.; -1.27e+03; 337.; 1.17e+03, 318., -1.27e+03; 318., -26.7, -216.; -1.27e+03; -216.; 287.; -982.; -1.27e+03; 337.; -1.27e+03; -216.; 287.; 337.; 287.; 3.19; 
+DEAL::118., 571., -293.; 571., 274., -940.; -293.; -940.; 52.3; 571., 274., -940.; 274., -18.0, -387.; -940.; -387.; 478.; -293.; -940.; 52.3; -940.; -387.; 478.; 52.3; 478.; 7.87; 
+DEAL::-101., -292., 355.; -292., -306., 510.; 355.; 510.; -440.; -292., -306., 510.; -306., -120., 323.; 510.; 323.; -493.; 355.; 510.; -440.; 510.; 323.; -493.; -440.; -493.; 124.; 
+DEAL::424., 1.45e+03, -766.; 1.45e+03, 597., -791.; -766.; -791.; 233.; 1.45e+03, 597., -791.; 597., -14.1, 96.3; -791.; 96.3; 257.; -766.; -791.; 233.; -791.; 96.3; 257.; 233.; 257.; -34.7; 
+DEAL::149., 466., -357.; 466., 223., -609.; -357.; -609.; 111.; 466., 223., -609.; 223., 8.21, -272.; -609.; -272.; 300.; -357.; -609.; 111.; -609.; -272.; 300.; 111.; 300.; 9.62; 
+DEAL::-121., -127., 553.; -127., -101., 214.; 553.; 214.; -483.; -127., -101., 214.; -101., -99.3, 215.; 214.; 215.; -116.; 553.; 214.; -483.; 214.; 215.; -116.; -483.; -116.; 92.1; 
+DEAL::424., -1.53e+03, 736.; -1.53e+03, 750., -1.00e+03; 736.; -1.00e+03; 206.; -1.53e+03, 750., -1.00e+03; 750., -10.9, 83.5; -1.00e+03; 83.5; -53.5; 736.; -1.00e+03; 206.; -1.00e+03; 83.5; -53.5; 206.; -53.5; 28.5; 
+DEAL::149., -606., 343.; -606., 376., -803.; 343.; -803.; 100.; -606., 376., -803.; 376., -10.4, 270.; -803.; 270.; -173.; 343.; -803.; 100.; -803.; 270.; -173.; 100.; -173.; -4.64; 
+DEAL::-121., 172., -515.; 172., -304., 339.; -515.; 339.; -474.; 172., -304., 339.; -304., 114., -348.; 339.; -348.; 187.; -515.; 339.; -474.; 339.; -348.; 187.; -474.; 187.; -119.; 
+DEAL::408., -995., -849.; -995., 238., 665.; -849.; 665.; 228.; -995., 238., 665.; 238., -14.1, -39.1; 665.; -39.1; -11.3; -849.; 665.; 228.; 665.; -39.1; -11.3; 228.; -11.3; -0.126; 
+DEAL::109., -360., -314.; -360., 151., 490.; -314.; 490.; 92.8; -360., 151., 490.; 151., 8.21, -242.; 490.; -242.; -169.; -314.; 490.; 92.8; 490.; -242.; -169.; 92.8; -169.; -2.49; 
+DEAL::117., -292., -442.; -292., 135., 294.; -442.; 294.; 390.; -292., 135., 294.; 135., -99.3, -258.; 294.; -258.; -312.; -442.; 294.; 390.; 294.; -258.; -312.; 390.; -312.; -127.; 
+DEAL::408., 991., 883.; 991., 157., 763.; 883.; 763.; 330.; 991., 157., 763.; 157., -10.9, 107.; 763.; 107.; 3.25; 883.; 763.; 330.; 763.; 107.; 3.25; 330.; 3.25; 3.19; 
+DEAL::109., 472., 285.; 472., 250., 625.; 285.; 625.; 164.; 472., 250., 625.; 250., -10.4, 237.; 625.; 237.; 134.; 285.; 625.; 164.; 625.; 237.; 134.; 164.; 134.; 7.87; 
+DEAL::117., 424., 445.; 424., 425., 853.; 445.; 853.; 448.; 424., 425., 853.; 425., 114., 471.; 853.; 471.; 540.; 445.; 853.; 448.; 853.; 471.; 540.; 448.; 540.; 124.; 

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.