ASPECT
handler.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2021 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 doc/COPYING. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef _aspect_volume_of_fluid_handler_h
22 #define _aspect_volume_of_fluid_handler_h
23 
24 #include <aspect/simulator.h>
28 
29 #include <boost/serialization/map.hpp>
30 
31 namespace aspect
32 {
39  template <int dim>
40  class VolumeOfFluidHandler : public SimulatorAccess<dim>
41  {
42  public:
46  VolumeOfFluidHandler(Simulator<dim> &simulator, ParameterHandler &prm);
47 
52  void edit_finite_element_variables (std::vector<VariableDeclaration<dim>> &vars);
53 
57  static
58  void declare_parameters (ParameterHandler &prm);
59 
63  void parse_parameters (ParameterHandler &prm);
64 
68  unsigned int get_n_fields() const;
69 
73  const std::string name_for_field_index(unsigned int i) const;
74 
79  const VolumeOfFluidField<dim> &field_struct_for_field_index(unsigned int i) const;
80 
84  double get_volume_fraction_threshold() const;
85 
89  unsigned int field_index_for_name(const std::string &fieldname) const;
90 
95  void initialize (ParameterHandler &prm);
96 
101 
105  void initialize_from_composition_field (const VolumeOfFluidField<dim> &field);
106 
110  void initialize_from_level_set (const VolumeOfFluidField<dim> &field);
111 
115  void update_volume_of_fluid_normals (const VolumeOfFluidField<dim> &field,
116  LinearAlgebra::BlockVector &solution);
117 
123  void update_volume_of_fluid_composition (const typename Simulator<dim>::AdvectionField &composition_field,
124  const VolumeOfFluidField<dim> &volume_of_fluid_field,
125  LinearAlgebra::BlockVector &solution);
126 
130  void do_volume_of_fluid_update (const typename Simulator<dim>::AdvectionField &advection_field);
131 
140  void assemble_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
141  const unsigned int calculation_dim,
142  const bool update_from_old_solution);
143 
148  void solve_volume_of_fluid_system (const VolumeOfFluidField<dim> &field);
149 
150 
151  private:
156 
162 
167 
172 
177  std::vector<VolumeOfFluidField<dim>> data;
178 
184 
188  static constexpr double volume_of_fluid_reconstruct_epsilon = 1e-13;
189 
194 
199  unsigned int n_init_samples;
200 
205  std::vector<std::string> volume_of_fluid_field_names;
206 
210  std::map<unsigned int, unsigned int> volume_of_fluid_composition_map_index;
211 
217  std::vector<VolumeOfFluid::VolumeOfFluidInputType::Kind> initialization_data_type;
218 
219  friend class Simulator<dim>;
220  };
221 
222 }
223 
224 #endif
static constexpr double volume_of_fluid_reconstruct_epsilon
Definition: handler.h:188
unsigned int get_n_fields() const
std::vector< std::string > volume_of_fluid_field_names
Definition: handler.h:205
void initialize_from_level_set(const VolumeOfFluidField< dim > &field)
const Simulator< dim > * simulator
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:268
void assemble_volume_of_fluid_system(const VolumeOfFluidField< dim > &field, const unsigned int calculation_dim, const bool update_from_old_solution)
std::vector< VolumeOfFluidField< dim > > data
Definition: handler.h:177
friend class Simulator< dim >
Definition: handler.h:219
void update_volume_of_fluid_composition(const typename Simulator< dim >::AdvectionField &composition_field, const VolumeOfFluidField< dim > &volume_of_fluid_field, LinearAlgebra::BlockVector &solution)
static void declare_parameters(ParameterHandler &prm)
Assemblers::VolumeOfFluidAssembler< dim > assembler
Definition: handler.h:166
void copy_local_to_global_volume_of_fluid_system(const internal::Assembly::CopyData::VolumeOfFluidSystem< dim > &data)
void initialize_from_composition_field(const VolumeOfFluidField< dim > &field)
const VolumeOfFluidField< dim > & field_struct_for_field_index(unsigned int i) const
std::vector< VolumeOfFluid::VolumeOfFluidInputType::Kind > initialization_data_type
Definition: handler.h:217
unsigned int field_index_for_name(const std::string &fieldname) const
Simulator< dim > & sim
Definition: handler.h:155
void edit_finite_element_variables(std::vector< VariableDeclaration< dim >> &vars)
void do_volume_of_fluid_update(const typename Simulator< dim >::AdvectionField &advection_field)
void initialize(ParameterHandler &prm)
double volume_of_fluid_solver_tolerance
Definition: handler.h:193
unsigned int n_volume_of_fluid_fields
Definition: handler.h:171
Definition: compat.h:59
std::map< unsigned int, unsigned int > volume_of_fluid_composition_map_index
Definition: handler.h:210
VolumeOfFluidHandler(Simulator< dim > &simulator, ParameterHandler &prm)
double get_volume_fraction_threshold() const
void update_volume_of_fluid_normals(const VolumeOfFluidField< dim > &field, LinearAlgebra::BlockVector &solution)
unsigned int n_init_samples
Definition: handler.h:199
void parse_parameters(ParameterHandler &prm)
const std::string name_for_field_index(unsigned int i) const
void solve_volume_of_fluid_system(const VolumeOfFluidField< dim > &field)