ASPECT
visualization.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2022 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>
124  class Interface
125  {
126  public:
146  explicit Interface (const std::string &physical_units = "");
147 
152  virtual ~Interface () = default;
153 
157  virtual void initialize ();
158 
162  virtual void update();
163 
175  virtual
176  std::string
177  get_physical_units () const;
178 
192  static
193  void
194  declare_parameters (ParameterHandler &prm);
195 
202  virtual
203  void
204  parse_parameters (ParameterHandler &prm);
205 
220  virtual
221  std::list<std::string>
223 
245  virtual
246  void save (std::map<std::string, std::string> &status_strings) const;
247 
259  virtual
260  void load (const std::map<std::string, std::string> &status_strings);
261 
262  private:
266  const std::string physical_units;
267  };
268 
269 
270 
283  template <int dim>
284  class CellDataVectorCreator : public Interface<dim>
285  {
286  public:
292  explicit CellDataVectorCreator (const std::string &physical_units = "");
293 
297  ~CellDataVectorCreator () override = default;
298 
317  virtual
318  std::pair<std::string, std::unique_ptr<Vector<float>>>
319  execute () const = 0;
320  };
321 
322 
329  template <int dim>
331  {
332  public:
338  virtual
339  ~SurfaceOnlyVisualization () = default;
340  };
341  }
342 
343 
359  template <int dim>
360  class Visualization : public Interface<dim>, public ::aspect::SimulatorAccess<dim>
361  {
362  public:
366  Visualization ();
367 
371  std::pair<std::string,std::string>
372  execute (TableHandler &statistics) override;
373 
377  void
378  update () override;
379 
398  static
399  void
400  register_visualization_postprocessor (const std::string &name,
401  const std::string &description,
402  void (*declare_parameters_function) (ParameterHandler &),
403  std::unique_ptr<VisualizationPostprocessors::Interface<dim>>(*factory_function) ());
404 
412  std::list<std::string>
413  required_other_postprocessors () const override;
414 
418  static
419  void
420  declare_parameters (ParameterHandler &prm);
421 
425  void
426  parse_parameters (ParameterHandler &prm) override;
427 
431  void save (std::map<std::string, std::string> &status_strings) const override;
432 
436  void load (const std::map<std::string, std::string> &status_strings) override;
437 
442  template <class Archive>
443  void serialize (Archive &ar, const unsigned int version);
444 
454  static
455  void
456  write_plugin_graph (std::ostream &output_stream);
457 
462  bool output_pointwise_stress_and_strain() const;
463 
467  DeclException1 (ExcPostprocessorNameNotFound,
468  std::string,
469  << "Could not find entry <"
470  << arg1
471  << "> among the names of registered postprocessors.");
472 
473  private:
480 
486 
494 
499  unsigned int last_output_timestep;
500 
505  unsigned int output_file_number;
506 
510  std::string output_format;
511 
518  unsigned int group_files;
519 
528 
540 
549 
556 
569 
577 
586 
594 
602 
611  void set_last_output_time (const double current_time);
612 
617  void mesh_changed_signal ();
618 
625  static
626  void writer (const std::string &filename,
627  const std::string &temporary_filename,
628  const std::string &file_contents);
629 
634  std::list<std::unique_ptr<VisualizationPostprocessors::Interface<dim>>> postprocessors;
635 
643  {
647  OutputHistory ();
648 
654  ~OutputHistory ();
655 
660  template <class Archive>
661  void serialize (Archive &ar, const unsigned int version);
662 
668 
673  std::string last_mesh_file_name;
674 
683  std::vector<std::pair<double,std::string>> times_and_pvtu_names;
684 
690  std::vector<std::vector<std::string>> output_file_names_by_timestep;
691 
698  std::vector<XDMFEntry> xdmf_entries;
699 
705  std::thread background_thread;
706  };
707 
713 
719 
737  template <typename DataOutType>
738  void write_master_files (const DataOutType &data_out,
739  const std::string &solution_file_prefix,
740  const std::vector<std::string> &filenames,
741  OutputHistory &output_history) const;
742 
743 
753  template <typename DataOutType>
754  std::string write_data_out_data(DataOutType &data_out,
755  OutputHistory &output_history,
756  const std::map<std::string,std::string> &visualization_field_names_and_units) const;
757  };
758  }
759 
760 
767 #define ASPECT_REGISTER_VISUALIZATION_POSTPROCESSOR(classname,name,description) \
768  template class classname<2>; \
769  template class classname<3>; \
770  namespace ASPECT_REGISTER_VISUALIZATION_POSTPROCESSOR_ ## classname \
771  { \
772  aspect::internal::Plugins::RegisterHelper<aspect::Postprocess::VisualizationPostprocessors::Interface<2>,classname<2>> \
773  dummy_ ## classname ## _2d (&aspect::Postprocess::Visualization<2>::register_visualization_postprocessor, \
774  name, description); \
775  aspect::internal::Plugins::RegisterHelper<aspect::Postprocess::VisualizationPostprocessors::Interface<3>,classname<3>> \
776  dummy_ ## classname ## _3d (&aspect::Postprocess::Visualization<3>::register_visualization_postprocessor, \
777  name, description); \
778  }
779 }
780 
781 
782 #endif
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
void average(const AveragingOperation operation, const typename DoFHandler< dim >::active_cell_iterator &cell, const Quadrature< dim > &quadrature_formula, const Mapping< dim > &mapping, MaterialModelOutputs< dim > &values_out)
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
static void declare_parameters(ParameterHandler &prm)
virtual void parse_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.")