ASPECT
visualization.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2024 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_postprocess_visualization_h
23 #define _aspect_postprocess_visualization_h
24 
27 #include <aspect/plugins.h>
28 
29 #include <deal.II/base/thread_management.h>
30 #include <deal.II/numerics/data_postprocessor.h>
31 #include <deal.II/base/data_out_base.h>
32 #include <deal.II/numerics/data_out.h>
33 
34 #include <thread>
35 
36 namespace aspect
37 {
38  namespace Postprocess
39  {
40  namespace VisualizationPostprocessors
41  {
45  inline void average_quantities(std::vector<Vector<double>> &quantities)
46  {
47  const unsigned int N = quantities.size();
48  const unsigned int M = quantities[0].size();
49  for (unsigned int m=0; m<M; ++m)
50  {
51  double sum = 0;
52  for (unsigned int q=0; q<N; ++q)
53  sum += quantities[q](m);
54 
55  const double average = sum/N;
56  for (unsigned int q=0; q<N; ++q)
57  quantities[q](m) = average;
58  }
59  }
60 
123  template <int dim>
125  {
126  public:
146  explicit Interface (const std::string &physical_units = "");
147 
159  virtual
160  std::string
161  get_physical_units () const;
162 
177  virtual
178  std::list<std::string>
180 
202  virtual
203  void save (std::map<std::string, std::string> &status_strings) const;
204 
216  virtual
217  void load (const std::map<std::string, std::string> &status_strings);
218 
219  private:
223  const std::string physical_units;
224  };
225 
226 
227 
240  template <int dim>
241  class CellDataVectorCreator : public Interface<dim>
242  {
243  public:
249  explicit CellDataVectorCreator (const std::string &physical_units = "");
250 
254  ~CellDataVectorCreator () override = default;
255 
274  virtual
275  std::pair<std::string, std::unique_ptr<Vector<float>>>
276  execute () const = 0;
277  };
278 
279 
286  template <int dim>
288  {
289  public:
295  virtual
296  ~SurfaceOnlyVisualization () = default;
297  };
298  }
299 
300 
316  template <int dim>
317  class Visualization : public Interface<dim>, public ::aspect::SimulatorAccess<dim>
318  {
319  public:
323  Visualization ();
324 
328  std::pair<std::string,std::string>
329  execute (TableHandler &statistics) override;
330 
334  void
335  update () override;
336 
355  static
356  void
357  register_visualization_postprocessor (const std::string &name,
358  const std::string &description,
359  void (*declare_parameters_function) (ParameterHandler &),
360  std::unique_ptr<VisualizationPostprocessors::Interface<dim>>(*factory_function) ());
361 
369  std::list<std::string>
370  required_other_postprocessors () const override;
371 
375  static
376  void
377  declare_parameters (ParameterHandler &prm);
378 
382  void
383  parse_parameters (ParameterHandler &prm) override;
384 
388  void save (std::map<std::string, std::string> &status_strings) const override;
389 
393  void load (const std::map<std::string, std::string> &status_strings) override;
394 
399  template <class Archive>
400  void serialize (Archive &ar, const unsigned int version);
401 
411  static
412  void
413  write_plugin_graph (std::ostream &output_stream);
414 
419  bool output_pointwise_stress_and_strain() const;
420 
424  DeclException1 (ExcPostprocessorNameNotFound,
425  std::string,
426  << "Could not find entry <"
427  << arg1
428  << "> among the names of registered postprocessors.");
429 
430  private:
437 
443 
451 
456  unsigned int last_output_timestep;
457 
462  unsigned int output_file_number;
463 
467  std::string output_format;
468 
475  unsigned int group_files;
476 
485 
497 
506 
513 
526 
534 
543 
551 
559 
567 
576  void set_last_output_time (const double current_time);
577 
582  void mesh_changed_signal ();
583 
590  static
591  void writer (const std::string &filename,
592  const std::string &temporary_filename,
593  const std::string &file_contents);
594 
599  std::list<std::unique_ptr<VisualizationPostprocessors::Interface<dim>>> postprocessors;
600 
608  {
612  OutputHistory ();
613 
619  ~OutputHistory ();
620 
625  template <class Archive>
626  void serialize (Archive &ar, const unsigned int version);
627 
633 
638  std::string last_mesh_file_name;
639 
648  std::vector<std::pair<double,std::string>> times_and_pvtu_names;
649 
655  std::vector<std::vector<std::string>> output_file_names_by_timestep;
656 
663  std::vector<XDMFEntry> xdmf_entries;
664 
670  std::thread background_thread;
671  };
672 
678 
684 
702  template <typename DataOutType>
703  void write_description_files (const DataOutType &data_out,
704  const std::string &solution_file_prefix,
705  const std::vector<std::string> &filenames,
706  OutputHistory &output_history) const;
707 
708 
718  template <typename DataOutType>
719  std::string write_data_out_data(DataOutType &data_out,
720  OutputHistory &output_history,
721  const std::map<std::string,std::string> &visualization_field_names_and_units) const;
722  };
723  }
724 
725 
732 #define ASPECT_REGISTER_VISUALIZATION_POSTPROCESSOR(classname,name,description) \
733  template class classname<2>; \
734  template class classname<3>; \
735  namespace ASPECT_REGISTER_VISUALIZATION_POSTPROCESSOR_ ## classname \
736  { \
737  aspect::internal::Plugins::RegisterHelper<aspect::Postprocess::VisualizationPostprocessors::Interface<2>,classname<2>> \
738  dummy_ ## classname ## _2d (&aspect::Postprocess::Visualization<2>::register_visualization_postprocessor, \
739  name, description); \
740  aspect::internal::Plugins::RegisterHelper<aspect::Postprocess::VisualizationPostprocessors::Interface<3>,classname<3>> \
741  dummy_ ## classname ## _3d (&aspect::Postprocess::Visualization<3>::register_visualization_postprocessor, \
742  name, description); \
743  }
744 }
745 
746 
747 #endif
virtual void parse_parameters(ParameterHandler &prm)
void write_plugin_graph(std::ostream &output_stream)
virtual void save(std::map< std::string, std::string > &status_strings) const
std::vector< std::pair< double, std::string > > times_and_pvtu_names
virtual std::list< std::string > required_other_postprocessors() const
std::list< std::unique_ptr< VisualizationPostprocessors::Interface< dim > > > postprocessors
Interface(const std::string &physical_units="")
void average_quantities(std::vector< Vector< double >> &quantities)
Definition: visualization.h:45
std::vector< std::vector< std::string > > output_file_names_by_timestep
Definition: compat.h:59
void average(const AveragingOperation operation, const typename DoFHandler< dim >::active_cell_iterator &cell, const Quadrature< dim > &quadrature_formula, const Mapping< dim > &mapping, const MaterialProperties::Property &requested_properties, MaterialModelOutputs< dim > &values_out)
static void declare_parameters(ParameterHandler &prm)
virtual void load(const std::map< std::string, std::string > &status_strings)
DeclException1(ProbabilityFunctionNegative, Point< dim >,<< "Your probability density function in the particle generator " "returned a negative probability density for the following position: "<< arg1<< ". Please check your function expression.")