ASPECT
assembly.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2016 - 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 doc/COPYING. If not see
18  <http://www.gnu.org/licenses/>.
19  */
20 
21 #ifndef _aspect_volume_of_fluid_assembly_h
22 #define _aspect_volume_of_fluid_assembly_h
23 
24 #include <aspect/global.h>
26 
27 #include <deal.II/fe/fe_values.h>
28 
29 namespace aspect
30 {
31  using namespace dealii;
32 
33  template <int dim>
34  class Simulator;
35 
36  template <int dim>
38 
39  namespace internal
40  {
41  namespace Assembly
42  {
43  namespace Scratch
44  {
48  template <int dim>
50  {
51  VolumeOfFluidSystem (const FiniteElement<dim> &finite_element,
52  const FiniteElement<dim> &volume_of_fluid_element,
53  const Mapping<dim> &mapping,
54  const Quadrature<dim> &quadrature,
55  const Quadrature<dim-1> &face_quadrature);
56 
58 
63  VolumeOfFluidSystem &operator=(const VolumeOfFluidSystem &) = default;
64 
65  FEValues<dim> finite_element_values;
67  FEFaceValues<dim> face_finite_element_values;
69  FESubfaceValues<dim> subface_finite_element_values;
70 
71  std::vector<types::global_dof_index> local_dof_indices;
72 
73  // Field for exact cell volume (needed in some of the face flux calculations)
74  double volume;
75 
84  std::vector<double> phi_field;
85  std::vector<double> face_phi_field;
86 
87  std::vector<double> old_field_values;
88  /* Vector for interface normal in the unit cell */
89  std::vector<Tensor<1,dim>> cell_i_n_values;
90  /* "Distance" from cell center to interface as d value for the interface in the form @f$\vec{n}\cdot\vec{x}=d@f$ */
91  std::vector<double> cell_i_d_values;
92 
93  std::vector<Tensor<1,dim>> face_current_velocity_values;
94  std::vector<Tensor<1,dim>> face_old_velocity_values;
95  std::vector<Tensor<1,dim>> face_old_old_velocity_values;
96 
97  std::vector<double> neighbor_old_values;
98  /* Vector for interface normal in the unit cell */
99  std::vector<Tensor<1,dim>> neighbor_i_n_values;
100  /* "Distance" from cell center to interface as d value for the interface in the form @f$\vec{n}\cdot\vec{x}=d@f$ */
101  std::vector<double> neighbor_i_d_values;
102  };
103  }
104 
105  namespace CopyData
106  {
110  template <int dim>
112  {
119  VolumeOfFluidSystem(const FiniteElement<dim> &finite_element);
121 
126  VolumeOfFluidSystem &operator=(const VolumeOfFluidSystem &) = default;
127 
132  FullMatrix<double> local_matrix;
133  Vector<double> local_rhs;
134 
142  std::vector<Vector<double>> local_face_rhs;
143  std::vector<FullMatrix<double>> local_face_matrices_ext_ext;
144 
154  std::vector<bool> face_contributions_mask;
155 
165  std::vector<types::global_dof_index> local_dof_indices;
166 
175  std::vector<std::vector<types::global_dof_index>> neighbor_dof_indices;
176  };
177  }
178  }
179  }
180 
181  namespace Assemblers
182  {
183 
187  template <int dim>
189  {
190  public:
195  void local_assemble_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
196  const unsigned int calc_dir,
197  const bool update_from_old,
198  const typename DoFHandler<dim>::active_cell_iterator &cell,
201 
205  void local_assemble_boundary_face_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
206  const bool update_from_old,
207  const typename DoFHandler<dim>::active_cell_iterator &cell,
208  const unsigned int face_no,
211 
215  void local_assemble_internal_face_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
216  const bool update_from_old,
217  const typename DoFHandler<dim>::active_cell_iterator &cell,
218  const unsigned int face_no,
221 
225  void set_volume_fraction_threshold(const double value);
226  };
227  }
228 }
229 
230 
231 #endif
std::vector< std::vector< types::global_dof_index > > neighbor_dof_indices
Definition: assembly.h:175
std::vector< Vector< double > > local_face_rhs
Definition: assembly.h:142
std::vector< Tensor< 1, dim > > cell_i_n_values
Definition: assembly.h:89
std::vector< Tensor< 1, dim > > neighbor_i_n_values
Definition: assembly.h:99
std::vector< Tensor< 1, dim > > face_current_velocity_values
Definition: assembly.h:93
Definition: compat.h:59
std::vector< FullMatrix< double > > local_face_matrices_ext_ext
Definition: assembly.h:143
std::vector< Tensor< 1, dim > > face_old_old_velocity_values
Definition: assembly.h:95
std::vector< Tensor< 1, dim > > face_old_velocity_values
Definition: assembly.h:94
Definition: compat.h:42
std::vector< types::global_dof_index > local_dof_indices
Definition: assembly.h:165
std::vector< types::global_dof_index > local_dof_indices
Definition: assembly.h:71