ASPECT
structured_data.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2014 - 2023 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 
22 #ifndef _aspect_structured_data_h
23 #define _aspect_structured_data_h
24 
25 #include <aspect/global.h>
27 
28 #include <array>
29 
30 namespace aspect
31 {
32  namespace Utilities
33  {
34  using namespace dealii;
35 
55  template <int dim>
57  {
58  public:
68  StructuredDataLookup(const unsigned int n_components,
69  const double scale_factor);
70 
77  explicit StructuredDataLookup(const double scale_factor);
78 
114  void reinit(const std::vector<std::string> &column_names,
115  std::vector<std::vector<double>> &&coordinate_values,
116  std::vector<Table<dim,double>> &&data_table,
117  const MPI_Comm mpi_communicator = MPI_COMM_SELF,
118  const unsigned int root_process = numbers::invalid_unsigned_int);
119 
137  void
138  load_ascii(const std::string &filename,
139  const MPI_Comm communicator);
140 
149  void
150  load_netcdf(const std::string &filename, const std::vector<std::string> &data_column_names = {});
151 
152 
167  void
168  load_file(const std::string &filename,
169  const MPI_Comm communicator);
170 
181  double
182  get_data(const Point<dim> &position,
183  const unsigned int component) const;
184 
194  Tensor<1,dim>
195  get_gradients(const Point<dim> &position,
196  const unsigned int component);
197 
204  std::vector<std::string>
205  get_column_names() const;
206 
212  bool
213  has_equidistant_coordinates() const;
214 
227  const std::vector<double> &
228  get_interpolation_point_coordinates(const unsigned int dimension) const;
229 
235  unsigned int
236  get_column_index_from_name(const std::string &column_name) const;
237 
243  std::string
244  get_column_name_from_index(const unsigned int column_index) const;
245 
249  double get_maximum_component_value(const unsigned int component) const;
250 
251  private:
255  unsigned int n_components;
256 
262  std::vector<std::string> data_component_names;
263 
269  std::vector<std::unique_ptr<Function<dim>>> data;
270 
274  std::array<std::vector<double>,dim> coordinate_values;
275 
279  std::vector<double> maximum_component_value;
280 
284  TableIndices<dim> table_points;
285 
290  const double scale_factor;
291 
297 
302  TableIndices<dim>
303  compute_table_indices(const TableIndices<dim> &sizes, const std::size_t idx) const;
304 
305  };
306 
311  template <int dim>
313  {
314  public:
318  AsciiDataBase();
319 
323  static
324  void
325  declare_parameters (ParameterHandler &prm,
326  const std::string &default_directory,
327  const std::string &default_filename,
328  const std::string &subsection_name = "Ascii data model");
329 
333  void
334  parse_parameters (ParameterHandler &prm,
335  const std::string &subsection_name = "Ascii data model");
336 
340  std::string data_directory;
341 
347  std::string data_file_name;
348 
355  double scale_factor;
356  };
357 
362  template <int dim>
364  {
365  public:
370 
375  virtual
376  void
377  initialize (const std::set<types::boundary_id> &boundary_ids,
378  const unsigned int components);
379 
386  void
387  update();
388 
392  double
393  get_data_component (const types::boundary_id boundary_indicator,
394  const Point<dim> &position,
395  const unsigned int component) const;
396 
400  double
401  get_maximum_component_value (const types::boundary_id boundary_indicator,
402  const unsigned int component) const;
403 
407  Tensor<1,dim-1>
408  vector_gradient(const types::boundary_id boundary_indicator,
409  const Point<dim> &p,
410  const unsigned int component) const;
411 
428  static
429  void
430  declare_parameters (ParameterHandler &prm,
431  const std::string &default_directory,
432  const std::string &default_filename,
433  const std::string &subsection_name = "Ascii data model",
434  const bool declare_time_dependent_parameters = true);
435 
449  void
450  parse_parameters (ParameterHandler &prm,
451  const std::string &subsection_name = "Ascii data model",
452  const bool parse_time_dependent_parameters = true);
453 
454  protected:
460 
466 
474 
480 
485  double time_weight;
486 
493 
498  std::map<types::boundary_id,
500 
504  std::map<types::boundary_id,
506 
510  void
511  update_data (const types::boundary_id boundary_id,
512  const bool reload_both_files);
513 
518  void
519  end_time_dependence ();
520 
524  std::string
525  create_filename (const int filenumber,
526  const types::boundary_id boundary_id) const;
527  };
528 
529 
530 
535  template <int dim>
537  {
538  public:
543 
548  virtual
549  void
550  initialize (const unsigned int components);
551 
552 
556  double
557  get_data_component (const Point<dim> &position,
558  const unsigned int component) const;
559 
563  static
564  void
565  declare_parameters (ParameterHandler &prm,
566  const std::string &default_directory,
567  const std::string &default_filename,
568  const std::string &subsection_name = "Ascii data model");
569 
573  void
574  parse_parameters (ParameterHandler &prm,
575  const std::string &subsection_name = "Ascii data model");
576 
577  protected:
582  std::unique_ptr<aspect::Utilities::StructuredDataLookup<dim>> lookup;
583 
589  std::unique_ptr<aspect::Utilities::StructuredDataLookup<3>> slice_lookup;
590 
596 
602  Tensor<2,3> rotation_matrix;
603  };
604 
605 
610  template <int dim>
612  {
613  public:
618 
623  virtual
624  void
625  initialize (const unsigned int components);
626 
627 
631  double
632  get_data_component (const Point<dim> &position,
633  const unsigned int component) const;
634 
635 
639  static
640  void
641  declare_parameters (ParameterHandler &prm,
642  const std::string &default_directory,
643  const std::string &default_filename,
644  const std::string &subsection_name = "Ascii data model");
645 
649  void
650  parse_parameters (ParameterHandler &prm,
651  const std::string &subsection_name = "Ascii data model");
652 
653  protected:
658  std::vector<std::unique_ptr<aspect::Utilities::StructuredDataLookup<dim-1>>> lookups;
659 
660  private:
661 
665  std::string data_directory;
666 
670  std::vector<std::string> data_file_names;
671 
676 
680  std::string interpolation_scheme;
681 
682 
683  };
684 
685 
689  template <int dim>
691  {
692  public:
697 
702  virtual
703  void
704  initialize (const MPI_Comm communicator);
705 
706 
710  double
711  get_data_component (const Point<1> &position,
712  const unsigned int component) const;
713 
720  std::vector<std::string>
721  get_column_names() const;
722 
737  const std::vector<double> &
738  get_interpolation_point_coordinates() const;
739 
745  unsigned int
746  get_column_index_from_name(const std::string &column_name) const;
747 
753  unsigned int
754  maybe_get_column_index_from_name(const std::string &column_name) const;
755 
761  std::string
762  get_column_name_from_index(const unsigned int column_index) const;
763  protected:
768  std::unique_ptr<aspect::Utilities::StructuredDataLookup<1>> lookup;
769  };
770 
771 
772 
773  template<int dim>
775  }
776 }
777 
778 #endif
std::unique_ptr< aspect::Utilities::StructuredDataLookup< 3 > > slice_lookup
std::map< types::boundary_id, std::unique_ptr< aspect::Utilities::StructuredDataLookup< dim-1 > > > lookups
std::vector< std::unique_ptr< Function< dim > > > data
std::vector< std::string > data_component_names
void declare_parameters(ParameterHandler &prm)
std::vector< std::string > data_file_names
std::vector< double > maximum_component_value
std::vector< std::unique_ptr< aspect::Utilities::StructuredDataLookup< dim-1 > > > lookups
std::unique_ptr< aspect::Utilities::StructuredDataLookup< 1 > > lookup
std::map< types::boundary_id, std::unique_ptr< aspect::Utilities::StructuredDataLookup< dim-1 > > > old_lookups
std::array< std::vector< double >, dim > coordinate_values
Definition: compat.h:42
std::unique_ptr< aspect::Utilities::StructuredDataLookup< dim > > lookup