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  template <int dim>
32  class Simulator;
33 
34  template <int dim>
36 
37  namespace internal
38  {
39  namespace Assembly
40  {
41  namespace Scratch
42  {
46  template <int dim>
48  {
49  VolumeOfFluidSystem (const FiniteElement<dim> &finite_element,
50  const FiniteElement<dim> &volume_of_fluid_element,
51  const Mapping<dim> &mapping,
52  const Quadrature<dim> &quadrature,
53  const Quadrature<dim-1> &face_quadrature);
54 
56 
62 
63  FEValues<dim> finite_element_values;
65  FEFaceValues<dim> face_finite_element_values;
67  FESubfaceValues<dim> subface_finite_element_values;
68 
69  std::vector<types::global_dof_index> local_dof_indices;
70 
71  // Field for exact cell volume (needed in some of the face flux calculations)
72  double volume;
73 
82  std::vector<double> phi_field;
83  std::vector<double> face_phi_field;
84 
85  std::vector<double> old_field_values;
86  /* Vector for interface normal in the unit cell */
87  std::vector<Tensor<1,dim>> cell_i_n_values;
88  /* "Distance" from cell center to interface as d value for the interface in the form @f$\vec{n}\cdot\vec{x}=d@f$ */
89  std::vector<double> cell_i_d_values;
90 
91  std::vector<Tensor<1,dim>> face_current_velocity_values;
92  std::vector<Tensor<1,dim>> face_old_velocity_values;
93  std::vector<Tensor<1,dim>> face_old_old_velocity_values;
94 
95  std::vector<double> neighbor_old_values;
96  /* Vector for interface normal in the unit cell */
97  std::vector<Tensor<1,dim>> neighbor_i_n_values;
98  /* "Distance" from cell center to interface as d value for the interface in the form @f$\vec{n}\cdot\vec{x}=d@f$ */
99  std::vector<double> neighbor_i_d_values;
100  };
101  }
102 
103  namespace CopyData
104  {
108  template <int dim>
110  {
117  VolumeOfFluidSystem(const FiniteElement<dim> &finite_element);
119 
125 
130  FullMatrix<double> local_matrix;
131  Vector<double> local_rhs;
132 
140  std::vector<Vector<double>> local_face_rhs;
141  std::vector<FullMatrix<double>> local_face_matrices_ext_ext;
142 
152  std::vector<bool> face_contributions_mask;
153 
163  std::vector<types::global_dof_index> local_dof_indices;
164 
173  std::vector<std::vector<types::global_dof_index>> neighbor_dof_indices;
174  };
175  }
176  }
177  }
178 
179  namespace Assemblers
180  {
181 
185  template <int dim>
187  {
188  public:
193  void local_assemble_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
194  const unsigned int calc_dir,
195  const bool update_from_old,
196  const typename DoFHandler<dim>::active_cell_iterator &cell,
199 
203  void local_assemble_boundary_face_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
204  const bool update_from_old,
205  const typename DoFHandler<dim>::active_cell_iterator &cell,
206  const unsigned int face_no,
209 
213  void local_assemble_internal_face_volume_of_fluid_system (const VolumeOfFluidField<dim> &field,
214  const bool update_from_old,
215  const typename DoFHandler<dim>::active_cell_iterator &cell,
216  const unsigned int face_no,
219 
223  void set_volume_fraction_threshold(const double value);
224  };
225  }
226 }
227 
228 
229 #endif
std::vector< std::vector< types::global_dof_index > > neighbor_dof_indices
Definition: assembly.h:173
VolumeOfFluidSystem & operator=(const VolumeOfFluidSystem &)=default
std::vector< Vector< double > > local_face_rhs
Definition: assembly.h:140
std::vector< Tensor< 1, dim > > cell_i_n_values
Definition: assembly.h:87
std::vector< Tensor< 1, dim > > neighbor_i_n_values
Definition: assembly.h:97
std::vector< Tensor< 1, dim > > face_current_velocity_values
Definition: assembly.h:91
Definition: compat.h:59
std::vector< FullMatrix< double > > local_face_matrices_ext_ext
Definition: assembly.h:141
std::vector< Tensor< 1, dim > > face_old_old_velocity_values
Definition: assembly.h:93
std::vector< Tensor< 1, dim > > face_old_velocity_values
Definition: assembly.h:92
VolumeOfFluidSystem(const FiniteElement< dim > &finite_element, const FiniteElement< dim > &volume_of_fluid_element, const Mapping< dim > &mapping, const Quadrature< dim > &quadrature, const Quadrature< dim-1 > &face_quadrature)
std::vector< types::global_dof_index > local_dof_indices
Definition: assembly.h:163
std::vector< types::global_dof_index > local_dof_indices
Definition: assembly.h:69