]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Provide the cell via DataOutFaces. Add tests.
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 21 Dec 2016 16:28:59 +0000 (09:28 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Sun, 25 Dec 2016 12:27:35 +0000 (05:27 -0700)
source/numerics/data_out_faces.cc
tests/numerics/data_out_faces_postprocessor_scalar_02.cc [new file with mode: 0644]
tests/numerics/data_out_faces_postprocessor_scalar_02.output [new file with mode: 0644]
tests/numerics/data_out_faces_postprocessor_vector_02.cc [new file with mode: 0644]
tests/numerics/data_out_faces_postprocessor_vector_02.output [new file with mode: 0644]
tests/numerics/data_out_postprocessor_scalar_02.cc

index cc999ba3f0dc546367908175d22d6756521f9de8..a81fd310bc04d1ec00217c9bbbb67c82914f6db1 100644 (file)
@@ -169,6 +169,12 @@ build_one_patch (const FaceDescriptor *cell_and_face,
                   if (update_flags & update_normal_vectors)
                     data.patch_values_scalar.normals = this_fe_patch_values.get_all_normal_vectors();
 
+                  const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(),
+                                                                              cell_and_face->first->level(),
+                                                                              cell_and_face->first->index(),
+                                                                              this->dof_data[dataset]->dof_handler);
+                  data.patch_values_scalar.template set_cell<DoFHandlerType> (dh_cell);
+
                   postprocessor->
                   evaluate_scalar_field(data.patch_values_scalar,
                                         data.postprocessed_values[dataset]);
@@ -194,6 +200,12 @@ build_one_patch (const FaceDescriptor *cell_and_face,
                   if (update_flags & update_normal_vectors)
                     data.patch_values_system.normals = this_fe_patch_values.get_all_normal_vectors();
 
+                  const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(),
+                                                                              cell_and_face->first->level(),
+                                                                              cell_and_face->first->index(),
+                                                                              this->dof_data[dataset]->dof_handler);
+                  data.patch_values_system.template set_cell<DoFHandlerType> (dh_cell);
+
                   postprocessor->
                   evaluate_vector_field(data.patch_values_system,
                                         data.postprocessed_values[dataset]);
diff --git a/tests/numerics/data_out_faces_postprocessor_scalar_02.cc b/tests/numerics/data_out_faces_postprocessor_scalar_02.cc
new file mode 100644 (file)
index 0000000..2345b90
--- /dev/null
@@ -0,0 +1,125 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+
+// tests DataPostprocessor: access the cell we are currently working
+// on for a scalar finite element field and DataOutFaces
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/base/function.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/numerics/matrix_tools.h>
+#include <deal.II/lac/vector.h>
+
+#include <deal.II/numerics/data_out_faces.h>
+#include <deal.II/numerics/data_postprocessor.h>
+#include <fstream>
+
+#include <deal.II/base/logstream.h>
+
+
+std::ofstream logfile("output");
+
+
+template <int dim>
+class MyPostprocessor : public DataPostprocessorScalar<dim>
+{
+public:
+  MyPostprocessor ()
+    :
+    DataPostprocessorScalar<dim> ("data", update_values)
+  {}
+
+  virtual
+  void
+  evaluate_scalar_field (const DataPostprocessorInputs::Scalar<dim> &input_data,
+                         std::vector<Vector<double> >               &computed_quantities) const
+  {
+    for (unsigned int q=0; q<input_data.solution_values.size(); ++q)
+      {
+        Assert (computed_quantities[q].size() == 1,
+                ExcInternalError());
+
+        // get the cell this all belongs to
+        typename DoFHandler<dim>::cell_iterator
+        cell = input_data.template get_cell<DoFHandler<dim> >();
+
+        Assert (input_data.solution_values[q]
+                ==
+                double(cell->active_cell_index()),
+                ExcInternalError());
+
+        computed_quantities[q][0] = input_data.solution_values[q];
+      }
+  }
+};
+
+
+
+template <int dim>
+void test ()
+{
+  Triangulation<dim>   triangulation;
+  FE_DGQ<dim>          fe(1);
+  DoFHandler<dim>      dof_handler(triangulation);
+
+  GridGenerator::hyper_cube (triangulation, 0, 1);
+  triangulation.refine_global (1);
+  triangulation.begin_active()->set_refine_flag ();
+  triangulation.execute_coarsening_and_refinement ();
+
+  dof_handler.distribute_dofs (fe);
+
+  // create a vector with piecewise constants, and set the elements of
+  // the vector so that on each cell the field has values equal to the
+  // cell's active_cell_index
+  Vector<double> solution (dof_handler.n_dofs());
+  for (typename DoFHandler<dim>::active_cell_iterator cell = dof_handler.begin_active();
+       cell != dof_handler.end(); ++cell)
+    {
+      std::vector<types::global_dof_index> local_dof_indices (cell->get_fe().dofs_per_cell);
+      cell->get_dof_indices (local_dof_indices);
+      for (unsigned int i=0; i<local_dof_indices.size(); ++i)
+        solution(local_dof_indices[i]) = cell->active_cell_index();
+    }
+
+  MyPostprocessor<dim> p;
+  DataOutFaces<dim> data_out;
+  data_out.attach_dof_handler (dof_handler);
+  data_out.add_data_vector (solution, p);
+  data_out.build_patches ();
+  data_out.write_gnuplot (logfile);
+}
+
+
+int main ()
+{
+  logfile << std::setprecision(2);
+  deallog << std::setprecision(2);
+
+  test<2>();
+  test<3>();
+
+  return 0;
+}
diff --git a/tests/numerics/data_out_faces_postprocessor_scalar_02.output b/tests/numerics/data_out_faces_postprocessor_scalar_02.output
new file mode 100644 (file)
index 0000000..6836551
--- /dev/null
@@ -0,0 +1,285 @@
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <data> 
+1 0 0 
+1 0.5 0 
+
+
+0.5 0 0 
+1 0 0 
+
+
+0 0.5 1 
+0 1 1 
+
+
+0 1 1 
+0.5 1 1 
+
+
+1 0.5 2 
+1 1 2 
+
+
+0.5 1 2 
+1 1 2 
+
+
+0 0 3 
+0 0.25 3 
+
+
+0 0 3 
+0.25 0 3 
+
+
+0.25 0 4 
+0.5 0 4 
+
+
+0 0.25 5 
+0 0.5 5 
+
+
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <z> <data> 
+1 0 0 0 
+1 0.5 0 0 
+
+1 0 0.5 0 
+1 0.5 0.5 0 
+
+
+0.5 0 0 0 
+0.5 0 0.5 0 
+
+1 0 0 0 
+1 0 0.5 0 
+
+
+0.5 0 0 0 
+1 0 0 0 
+
+0.5 0.5 0 0 
+1 0.5 0 0 
+
+
+0 0.5 0 1 
+0 1 0 1 
+
+0 0.5 0.5 1 
+0 1 0.5 1 
+
+
+0 1 0 1 
+0 1 0.5 1 
+
+0.5 1 0 1 
+0.5 1 0.5 1 
+
+
+0 0.5 0 1 
+0.5 0.5 0 1 
+
+0 1 0 1 
+0.5 1 0 1 
+
+
+1 0.5 0 2 
+1 1 0 2 
+
+1 0.5 0.5 2 
+1 1 0.5 2 
+
+
+0.5 1 0 2 
+0.5 1 0.5 2 
+
+1 1 0 2 
+1 1 0.5 2 
+
+
+0.5 0.5 0 2 
+1 0.5 0 2 
+
+0.5 1 0 2 
+1 1 0 2 
+
+
+0 0 0.5 3 
+0 0.5 0.5 3 
+
+0 0 1 3 
+0 0.5 1 3 
+
+
+0 0 0.5 3 
+0 0 1 3 
+
+0.5 0 0.5 3 
+0.5 0 1 3 
+
+
+0 0 1 3 
+0.5 0 1 3 
+
+0 0.5 1 3 
+0.5 0.5 1 3 
+
+
+1 0 0.5 4 
+1 0.5 0.5 4 
+
+1 0 1 4 
+1 0.5 1 4 
+
+
+0.5 0 0.5 4 
+0.5 0 1 4 
+
+1 0 0.5 4 
+1 0 1 4 
+
+
+0.5 0 1 4 
+1 0 1 4 
+
+0.5 0.5 1 4 
+1 0.5 1 4 
+
+
+0 0.5 0.5 5 
+0 1 0.5 5 
+
+0 0.5 1 5 
+0 1 1 5 
+
+
+0 1 0.5 5 
+0 1 1 5 
+
+0.5 1 0.5 5 
+0.5 1 1 5 
+
+
+0 0.5 1 5 
+0.5 0.5 1 5 
+
+0 1 1 5 
+0.5 1 1 5 
+
+
+1 0.5 0.5 6 
+1 1 0.5 6 
+
+1 0.5 1 6 
+1 1 1 6 
+
+
+0.5 1 0.5 6 
+0.5 1 1 6 
+
+1 1 0.5 6 
+1 1 1 6 
+
+
+0.5 0.5 1 6 
+1 0.5 1 6 
+
+0.5 1 1 6 
+1 1 1 6 
+
+
+0 0 0 7 
+0 0.25 0 7 
+
+0 0 0.25 7 
+0 0.25 0.25 7 
+
+
+0 0 0 7 
+0 0 0.25 7 
+
+0.25 0 0 7 
+0.25 0 0.25 7 
+
+
+0 0 0 7 
+0.25 0 0 7 
+
+0 0.25 0 7 
+0.25 0.25 0 7 
+
+
+0.25 0 0 8 
+0.25 0 0.25 8 
+
+0.5 0 0 8 
+0.5 0 0.25 8 
+
+
+0.25 0 0 8 
+0.5 0 0 8 
+
+0.25 0.25 0 8 
+0.5 0.25 0 8 
+
+
+0 0.25 0 9 
+0 0.5 0 9 
+
+0 0.25 0.25 9 
+0 0.5 0.25 9 
+
+
+0 0.25 0 9 
+0.25 0.25 0 9 
+
+0 0.5 0 9 
+0.25 0.5 0 9 
+
+
+0.25 0.25 0 10 
+0.5 0.25 0 10 
+
+0.25 0.5 0 10 
+0.5 0.5 0 10 
+
+
+0 0 0.25 11 
+0 0.25 0.25 11 
+
+0 0 0.5 11 
+0 0.25 0.5 11 
+
+
+0 0 0.25 11 
+0 0 0.5 11 
+
+0.25 0 0.25 11 
+0.25 0 0.5 11 
+
+
+0.25 0 0.25 12 
+0.25 0 0.5 12 
+
+0.5 0 0.25 12 
+0.5 0 0.5 12 
+
+
+0 0.25 0.25 13 
+0 0.5 0.25 13 
+
+0 0.25 0.5 13 
+0 0.5 0.5 13 
+
+
diff --git a/tests/numerics/data_out_faces_postprocessor_vector_02.cc b/tests/numerics/data_out_faces_postprocessor_vector_02.cc
new file mode 100644 (file)
index 0000000..53edf76
--- /dev/null
@@ -0,0 +1,128 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+
+// tests DataPostprocessor: access the cell we are currently working
+// on for a vector finite element field and DataOutFaces
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/base/function.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/numerics/matrix_tools.h>
+#include <deal.II/lac/vector.h>
+
+#include <deal.II/numerics/data_out_faces.h>
+#include <deal.II/numerics/data_postprocessor.h>
+#include <fstream>
+
+#include <deal.II/base/logstream.h>
+
+
+std::ofstream logfile("output");
+
+
+template <int dim>
+class MyPostprocessor : public DataPostprocessorVector<dim>
+{
+public:
+  MyPostprocessor ()
+    :
+    DataPostprocessorVector<dim> ("data", update_values)
+  {}
+
+  virtual
+  void
+  evaluate_vector_field (const DataPostprocessorInputs::Vector<dim> &input_data,
+                         std::vector<Vector<double> >                    &computed_quantities) const
+  {
+    for (unsigned int q=0; q<input_data.solution_values.size(); ++q)
+      {
+        // we only have one scalar field to output
+        Assert (input_data.solution_values[q].size() == 2,
+                ExcInternalError());
+        Assert (computed_quantities[q].size() == dim,
+                ExcInternalError());
+
+        // get the cell this all belongs to
+        typename DoFHandler<dim>::cell_iterator
+        cell = input_data.template get_cell<DoFHandler<dim> >();
+
+        Assert (input_data.solution_values[q](0)
+                ==
+                double(cell->active_cell_index()),
+                ExcInternalError());
+
+        computed_quantities[q][0] = input_data.solution_values[q](0);
+      }
+  }
+};
+
+
+
+template <int dim>
+void test ()
+{
+  Triangulation<dim>   triangulation;
+  FESystem<dim>        fe(FE_DGQ<dim>(0),2);
+  DoFHandler<dim>      dof_handler(triangulation);
+
+  GridGenerator::hyper_cube (triangulation, 0, 1);
+  triangulation.refine_global (1);
+  triangulation.begin_active()->set_refine_flag ();
+  triangulation.execute_coarsening_and_refinement ();
+
+  dof_handler.distribute_dofs (fe);
+
+  // create a vector with piecewise constants, and set the elements of
+  // the vector so that on each cell the field has values equal to the
+  // cell's active_cell_index
+  Vector<double> solution (dof_handler.n_dofs());
+  for (typename DoFHandler<dim>::active_cell_iterator cell = dof_handler.begin_active();
+       cell != dof_handler.end(); ++cell)
+    {
+      std::vector<types::global_dof_index> local_dof_indices (cell->get_fe().dofs_per_cell);
+      cell->get_dof_indices (local_dof_indices);
+      for (unsigned int i=0; i<local_dof_indices.size(); ++i)
+        solution(local_dof_indices[i]) = cell->active_cell_index();
+    }
+
+  MyPostprocessor<dim> p;
+  DataOutFaces<dim> data_out;
+  data_out.attach_dof_handler (dof_handler);
+  data_out.add_data_vector (solution, p);
+  data_out.build_patches ();
+  data_out.write_gnuplot (logfile);
+}
+
+
+int main ()
+{
+  logfile << std::setprecision(2);
+  deallog << std::setprecision(2);
+
+  test<2>();
+  test<3>();
+
+  return 0;
+}
diff --git a/tests/numerics/data_out_faces_postprocessor_vector_02.output b/tests/numerics/data_out_faces_postprocessor_vector_02.output
new file mode 100644 (file)
index 0000000..85f2587
--- /dev/null
@@ -0,0 +1,285 @@
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <data> <data> 
+1 0 0 0 
+1 0.5 0 0 
+
+
+0.5 0 0 0 
+1 0 0 0 
+
+
+0 0.5 1 0 
+0 1 1 0 
+
+
+0 1 1 0 
+0.5 1 1 0 
+
+
+1 0.5 2 0 
+1 1 2 0 
+
+
+0.5 1 2 0 
+1 1 2 0 
+
+
+0 0 3 0 
+0 0.25 3 0 
+
+
+0 0 3 0 
+0.25 0 3 0 
+
+
+0.25 0 4 0 
+0.5 0 4 0 
+
+
+0 0.25 5 0 
+0 0.5 5 0 
+
+
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <z> <data> <data> <data> 
+1 0 0 0 0 0 
+1 0.5 0 0 0 0 
+
+1 0 0.5 0 0 0 
+1 0.5 0.5 0 0 0 
+
+
+0.5 0 0 0 0 0 
+0.5 0 0.5 0 0 0 
+
+1 0 0 0 0 0 
+1 0 0.5 0 0 0 
+
+
+0.5 0 0 0 0 0 
+1 0 0 0 0 0 
+
+0.5 0.5 0 0 0 0 
+1 0.5 0 0 0 0 
+
+
+0 0.5 0 1 0 0 
+0 1 0 1 0 0 
+
+0 0.5 0.5 1 0 0 
+0 1 0.5 1 0 0 
+
+
+0 1 0 1 0 0 
+0 1 0.5 1 0 0 
+
+0.5 1 0 1 0 0 
+0.5 1 0.5 1 0 0 
+
+
+0 0.5 0 1 0 0 
+0.5 0.5 0 1 0 0 
+
+0 1 0 1 0 0 
+0.5 1 0 1 0 0 
+
+
+1 0.5 0 2 0 0 
+1 1 0 2 0 0 
+
+1 0.5 0.5 2 0 0 
+1 1 0.5 2 0 0 
+
+
+0.5 1 0 2 0 0 
+0.5 1 0.5 2 0 0 
+
+1 1 0 2 0 0 
+1 1 0.5 2 0 0 
+
+
+0.5 0.5 0 2 0 0 
+1 0.5 0 2 0 0 
+
+0.5 1 0 2 0 0 
+1 1 0 2 0 0 
+
+
+0 0 0.5 3 0 0 
+0 0.5 0.5 3 0 0 
+
+0 0 1 3 0 0 
+0 0.5 1 3 0 0 
+
+
+0 0 0.5 3 0 0 
+0 0 1 3 0 0 
+
+0.5 0 0.5 3 0 0 
+0.5 0 1 3 0 0 
+
+
+0 0 1 3 0 0 
+0.5 0 1 3 0 0 
+
+0 0.5 1 3 0 0 
+0.5 0.5 1 3 0 0 
+
+
+1 0 0.5 4 0 0 
+1 0.5 0.5 4 0 0 
+
+1 0 1 4 0 0 
+1 0.5 1 4 0 0 
+
+
+0.5 0 0.5 4 0 0 
+0.5 0 1 4 0 0 
+
+1 0 0.5 4 0 0 
+1 0 1 4 0 0 
+
+
+0.5 0 1 4 0 0 
+1 0 1 4 0 0 
+
+0.5 0.5 1 4 0 0 
+1 0.5 1 4 0 0 
+
+
+0 0.5 0.5 5 0 0 
+0 1 0.5 5 0 0 
+
+0 0.5 1 5 0 0 
+0 1 1 5 0 0 
+
+
+0 1 0.5 5 0 0 
+0 1 1 5 0 0 
+
+0.5 1 0.5 5 0 0 
+0.5 1 1 5 0 0 
+
+
+0 0.5 1 5 0 0 
+0.5 0.5 1 5 0 0 
+
+0 1 1 5 0 0 
+0.5 1 1 5 0 0 
+
+
+1 0.5 0.5 6 0 0 
+1 1 0.5 6 0 0 
+
+1 0.5 1 6 0 0 
+1 1 1 6 0 0 
+
+
+0.5 1 0.5 6 0 0 
+0.5 1 1 6 0 0 
+
+1 1 0.5 6 0 0 
+1 1 1 6 0 0 
+
+
+0.5 0.5 1 6 0 0 
+1 0.5 1 6 0 0 
+
+0.5 1 1 6 0 0 
+1 1 1 6 0 0 
+
+
+0 0 0 7 0 0 
+0 0.25 0 7 0 0 
+
+0 0 0.25 7 0 0 
+0 0.25 0.25 7 0 0 
+
+
+0 0 0 7 0 0 
+0 0 0.25 7 0 0 
+
+0.25 0 0 7 0 0 
+0.25 0 0.25 7 0 0 
+
+
+0 0 0 7 0 0 
+0.25 0 0 7 0 0 
+
+0 0.25 0 7 0 0 
+0.25 0.25 0 7 0 0 
+
+
+0.25 0 0 8 0 0 
+0.25 0 0.25 8 0 0 
+
+0.5 0 0 8 0 0 
+0.5 0 0.25 8 0 0 
+
+
+0.25 0 0 8 0 0 
+0.5 0 0 8 0 0 
+
+0.25 0.25 0 8 0 0 
+0.5 0.25 0 8 0 0 
+
+
+0 0.25 0 9 0 0 
+0 0.5 0 9 0 0 
+
+0 0.25 0.25 9 0 0 
+0 0.5 0.25 9 0 0 
+
+
+0 0.25 0 9 0 0 
+0.25 0.25 0 9 0 0 
+
+0 0.5 0 9 0 0 
+0.25 0.5 0 9 0 0 
+
+
+0.25 0.25 0 10 0 0 
+0.5 0.25 0 10 0 0 
+
+0.25 0.5 0 10 0 0 
+0.5 0.5 0 10 0 0 
+
+
+0 0 0.25 11 0 0 
+0 0.25 0.25 11 0 0 
+
+0 0 0.5 11 0 0 
+0 0.25 0.5 11 0 0 
+
+
+0 0 0.25 11 0 0 
+0 0 0.5 11 0 0 
+
+0.25 0 0.25 11 0 0 
+0.25 0 0.5 11 0 0 
+
+
+0.25 0 0.25 12 0 0 
+0.25 0 0.5 12 0 0 
+
+0.5 0 0.25 12 0 0 
+0.5 0 0.5 12 0 0 
+
+
+0 0.25 0.25 13 0 0 
+0 0.5 0.25 13 0 0 
+
+0 0.25 0.5 13 0 0 
+0 0.5 0.5 13 0 0 
+
+
index 560edbdedc82a49a4527d4dbc7c8378c10b9cc2e..be28e10de27c9a9dc932fe2e87123b3b3eff341f 100644 (file)
@@ -63,7 +63,7 @@ public:
 
         // get the cell this all belongs to
         typename DoFHandler<dim>::cell_iterator
-          cell = input_data.template get_cell<DoFHandler<dim> >();
+        cell = input_data.template get_cell<DoFHandler<dim> >();
 
         Assert (input_data.solution_values[q]
                 ==
@@ -103,7 +103,7 @@ void test ()
       for (unsigned int i=0; i<local_dof_indices.size(); ++i)
         solution(local_dof_indices[i]) = cell->active_cell_index();
     }
-  
+
   MyPostprocessor<dim> p;
   DataOut<dim> data_out;
   data_out.attach_dof_handler (dof_handler);

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.