]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Renamed field_data as cell_data and changed type from std::vector<double> to Vector...
authorVaishnavi Kale <vaishnav@ualberta.ca>
Mon, 9 Dec 2024 22:25:46 +0000 (15:25 -0700)
committerVaishnavi Kale <vaishnav@ualberta.ca>
Mon, 9 Dec 2024 22:25:46 +0000 (15:25 -0700)
include/deal.II/grid/grid_in.h
source/grid/grid_in.cc
tests/grid/grid_in_vtk_2d_field_data.cc

index 20ecab5d366211a717991a57c4f062e0f6a31f0f..90b33327acfe4530a21a5928472c426fbe64e077 100644 (file)
@@ -21,6 +21,7 @@
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/observer_pointer.h>
 #include <deal.II/base/point.h>
+#include <deal.II/lac/vector.h>
 
 #include <iostream>
 #include <map>
@@ -769,22 +770,23 @@ public:
   get_format_names();
 
   /**
-   * Return a map containing field data associated with the elements of an
+   * Return a map containing cell data associated with the elements of an
    * external vtk format mesh imported using read_vtk().
    * The format of the returned map is as
    * follows:
    * - std::string stores the name of the field data (identifier) as specified
    * in the external mesh
-   * - std::vector<double> stores value for the given identifier in each cell.
-   * To access the value, use field_data[name_field][cell_id]. For example, if
-   * the vtk mesh contains field data 'Density' defined at the cells,
-   * field_data['Density'][0] would be the density defined at cell ID '0',
+   * - Vector<double> stores value for the given identifier in each cell.
+   * To access the value, use cell_data[name_field][cell->active_cell_index()].
+   *  
+   * For example, if the vtk mesh contains field data 'Density' defined at the 
+   * cells, cell_data['Density'][0] would be the density defined at cell ID '0',
    * which corresponds to index 0 of the vector. The length of the vector in
-   * field_data['Density'] would equal the number of elements in the coarse
+   * cell_data['Density'] would equal the number of elements in the coarse
    * mesh.
    */
-  const std::map<std::string, std::vector<double>> &
-  get_field_data() const;
+  const std::map<std::string, Vector<double>> &
+  get_cell_data() const;
 
   /**
    * Exception
@@ -973,11 +975,11 @@ private:
    * The format is as follows:
    * - std::string stores the name of the field data (identifier) as specified
    * in the external mesh
-   * - std::vector<double> stores value for the given identifier in each
-   * cell_id.
-   * To access the value use field_data[name_field][cell_id].
+   * - Vector<double> stores value for the given identifier in each cell id.
+   *
+   * To access the value use cell_data[name_field][cell->active_cell_index()].
    */
-  std::map<std::string, std::vector<double>> field_data;
+  std::map<std::string, Vector<double>> cell_data;
 };
 
 /* -------------- declaration of explicit specializations ------------- */
index 84054b8b1925477d54bba510597e9d69f6046b33..01523b61401a3c2ec644c9989f21b1840a0fb1ee 100644 (file)
@@ -625,15 +625,15 @@ GridIn<dim, spacedim>::read_vtk(std::istream &in)
                         subcelldata.boundary_lines.size()) +
                       ") in 2d."));
                   in >> data_type;
-                  std::vector<double> temp_data;
-                  temp_data.resize(n_ids);
+                  Vector<double> temp_data;
+                  temp_data.reinit(n_ids);
                   for (unsigned int j = 0; j < n_ids; ++j)
                     {
                       in >> data;
                       if (j < cells.size())
                         temp_data[j] = data;
                     }
-                  this->field_data[section_name] = std::move(temp_data);
+                  this->cell_data[section_name] = std::move(temp_data);
                 }
             }
           else
@@ -654,10 +654,10 @@ GridIn<dim, spacedim>::read_vtk(std::istream &in)
 }
 
 template <int dim, int spacedim>
-const std::map<std::string, std::vector<double>> &
-GridIn<dim, spacedim>::get_field_data() const
+const std::map<std::string, Vector<double>> &
+GridIn<dim, spacedim>::get_cell_data() const
 {
-  return this->field_data;
+  return this->cell_data;
 }
 
 template <int dim, int spacedim>
index 152f30a8f054352ec3d858fbf3c85728788c075c..537fb26ddf59f14c8e87df00d401bd74d9325977 100644 (file)
@@ -45,16 +45,16 @@ check_file(const std::string name, typename GridIn<dim>::Format format)
   gridin.read(name, format);
   deallog << '\t' << triangulation.n_vertices() << '\t'
           << triangulation.n_cells() << std::endl;
-  std::map<std::string, std::vector<double>> field_data;
+  std::map<std::string, Vector<double>> cell_data;
 
-  field_data = gridin.get_field_data();
-  if (field_data.size() > 0)
+  cell_data = gridin.get_cell_data();
+  if (cell_data.size() > 0)
     {
-      std::map<std::string, std::vector<double>>::const_iterator iter =
-        field_data.find("Density");
-      if (iter != field_data.end())
+      std::map<std::string, Vector<double>>::const_iterator iter =
+        cell_data.find("Density");
+      if (iter != cell_data.end())
         {
-          std::vector<double> cell_density = iter->second;
+          Vector<double> cell_density = iter->second;
           for (typename Triangulation<dim>::active_cell_iterator cell =
                  triangulation.begin_active();
                cell != triangulation.end();

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.