ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2017 - 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 #ifndef _aspect_simulator_assemblers_interface_h
22 #define _aspect_simulator_assemblers_interface_h
23 
24 #include <aspect/global.h>
27 
28 #include <deal.II/fe/fe_values.h>
29 
30 namespace aspect
31 {
32  using namespace dealii;
33 
34  template <int dim>
35  class Simulator;
36 
37  struct AdvectionField;
38 
46  namespace internal
47  {
48  namespace Assembly
49  {
50  namespace Scratch
51  {
58  template <int dim>
59  struct ScratchBase
60  {
62  :
63  cell(),
64  face_number(numbers::invalid_unsigned_int)
65  {}
66 
67  ScratchBase(const ScratchBase &scratch)
68  :
69  cell(scratch.cell),
70  face_number(scratch.face_number)
71  {}
72 
73  virtual ~ScratchBase () = default;
74 
78  typename DoFHandler<dim>::active_cell_iterator cell;
79 
86  unsigned face_number;
87  };
88 
93  template <int dim>
94  struct StokesPreconditioner: public ScratchBase<dim>
95  {
96  StokesPreconditioner (const FiniteElement<dim> &finite_element,
97  const Quadrature<dim> &quadrature,
98  const Mapping<dim> &mapping,
99  const UpdateFlags update_flags,
100  const unsigned int n_compositional_fields,
101  const unsigned int stokes_dofs_per_cell,
102  const bool add_compaction_pressure,
103  const bool rebuild_matrix,
104  const bool use_bfbt);
105  StokesPreconditioner (const StokesPreconditioner &scratch);
106 
107  ~StokesPreconditioner () override;
108 
109  FEValues<dim> finite_element_values;
110 
111  void reinit (const typename DoFHandler<dim>::active_cell_iterator &cell_ref);
112 
113  std::vector<types::global_dof_index> local_dof_indices;
114  std::vector<unsigned int> dof_component_indices;
115  std::vector<SymmetricTensor<2,dim>> grads_phi_u;
116  std::vector<double> div_phi_u;
117  std::vector<double> phi_p;
118  std::vector<Tensor<1,dim>> phi_u;
119  std::vector<double> phi_p_c;
120  std::vector<Tensor<1,dim>> grad_phi_p;
121 
128 
135  };
136 
137 
138 
149  template <int dim>
150  struct StokesSystem : public StokesPreconditioner<dim>
151  {
152  StokesSystem (const FiniteElement<dim> &finite_element,
153  const Mapping<dim> &mapping,
154  const Quadrature<dim> &quadrature,
155  const Quadrature<dim-1> &face_quadrature,
156  const UpdateFlags update_flags,
157  const UpdateFlags face_update_flags,
158  const unsigned int n_compositional_fields,
159  const unsigned int stokes_dofs_per_cell,
160  const bool add_compaction_pressure,
161  const bool use_reference_density_profile,
162  const bool rebuild_stokes_matrix,
163  const bool rebuild_newton_stokes_matrix,
164  const bool use_bfbt);
165 
166  StokesSystem (const StokesSystem<dim> &scratch);
167 
168  FEFaceValues<dim> face_finite_element_values;
169 
171 
172  void reinit (const typename DoFHandler<dim>::active_cell_iterator &cell_ref,
173  const unsigned face_number_ref);
174 
175  std::vector<Tensor<1,dim>> phi_u;
176  std::vector<Tensor<1,dim>> velocity_values;
177  std::vector<double> velocity_divergence;
178  std::vector<Tensor<1,dim>> temperature_gradients;
179 
190 
200  std::vector<double> reference_densities;
202 
209  };
210 
211 
212 
217  template <int dim>
218  struct AdvectionSystem: public ScratchBase<dim>
219  {
220  AdvectionSystem (const FiniteElement<dim> &finite_element,
221  const FiniteElement<dim> &advection_element,
222  const Mapping<dim> &mapping,
223  const Quadrature<dim> &quadrature,
224  const Quadrature<dim-1> &face_quadrature,
225  const UpdateFlags update_flags,
226  const UpdateFlags face_update_flags,
227  const unsigned int n_compositional_fields,
228  const typename Simulator<dim>::AdvectionField &field);
229  AdvectionSystem (const AdvectionSystem &scratch);
230 
231  FEValues<dim> finite_element_values;
232 
233  void reinit (const typename DoFHandler<dim>::active_cell_iterator &cell_ref);
234 
235  std::unique_ptr<FEFaceValues<dim>> face_finite_element_values;
236  std::unique_ptr<FEFaceValues<dim>> neighbor_face_finite_element_values;
237  std::unique_ptr<FESubfaceValues<dim>> subface_finite_element_values;
238 
239  std::vector<types::global_dof_index> local_dof_indices;
240 
249  std::vector<double> phi_field;
250  std::vector<Tensor<1,dim>> grad_phi_field;
251  std::vector<double> laplacian_phi_field;
252  std::vector<double> face_phi_field;
253  std::vector<Tensor<1,dim>> face_grad_phi_field;
254  std::vector<double> neighbor_face_phi_field;
255  std::vector<Tensor<1,dim>> neighbor_face_grad_phi_field;
256 
257  std::vector<Tensor<1,dim>> old_velocity_values;
258  std::vector<Tensor<1,dim>> old_old_velocity_values;
259 
260  std::vector<double> old_pressure;
261  std::vector<double> old_old_pressure;
262  std::vector<Tensor<1,dim>> old_pressure_gradients;
263  std::vector<Tensor<1,dim>> old_old_pressure_gradients;
264 
265  std::vector<SymmetricTensor<2,dim>> old_strain_rates;
266  std::vector<SymmetricTensor<2,dim>> old_old_strain_rates;
267 
268  std::vector<double> old_temperature_values;
269  std::vector<double> old_old_temperature_values;
270 
271  std::vector<double> old_field_values;
272  std::vector<double> old_old_field_values;
273  std::vector<Tensor<1,dim>> old_field_grads;
274  std::vector<Tensor<1,dim>> old_old_field_grads;
275  std::vector<double> old_field_laplacians;
276  std::vector<double> old_old_field_laplacians;
277 
278  std::vector<std::vector<double>> old_composition_values;
279  std::vector<std::vector<double>> old_old_composition_values;
280 
281  std::vector<double> current_temperature_values;
282  std::vector<Tensor<1,dim>> current_velocity_values;
283  std::vector<Tensor<1,dim>> face_current_velocity_values;
284  std::vector<Tensor<1,dim>> mesh_velocity_values;
285  std::vector<Tensor<1,dim>> face_mesh_velocity_values;
286 
287  std::vector<SymmetricTensor<2,dim>> current_strain_rates;
288  std::vector<std::vector<double>> current_composition_values;
289  std::vector<double> current_velocity_divergences;
290 
297 
300 
303 
314 
323 
329  };
330  }
331 
332 
333 
346  namespace CopyData
347  {
353  template <int dim>
355  {
356  virtual ~CopyDataBase () = default;
357  };
358 
364  template <int dim>
365  struct StokesPreconditioner: public CopyDataBase<dim>
366  {
367  StokesPreconditioner (const unsigned int stokes_dofs_per_cell);
368 
369  StokesPreconditioner (const StokesPreconditioner &data);
370 
371  ~StokesPreconditioner () override = default;
372  StokesPreconditioner<dim> &operator= (const StokesPreconditioner<dim> &data) = default;
373 
374  FullMatrix<double> local_matrix;
376  std::vector<types::global_dof_index> local_dof_indices;
377 
383  void extract_stokes_dof_indices(const std::vector<types::global_dof_index> &all_dof_indices,
384  const Introspection<dim> &introspection,
385  const FiniteElement<dim> &finite_element);
386  };
387 
393  template <int dim>
394  struct StokesSystem : public StokesPreconditioner<dim>
395  {
396  StokesSystem (const unsigned int stokes_dofs_per_cell,
397  const bool do_pressure_rhs_compatibility_modification);
398  StokesSystem (const StokesSystem<dim> &data);
399 
400  ~StokesSystem () override = default;
401  StokesSystem<dim> &operator= (const StokesSystem<dim> &data) = default;
402 
403  Vector<double> local_rhs;
405  };
406 
412  template <int dim>
413  struct AdvectionSystem: public CopyDataBase<dim>
414  {
423  AdvectionSystem (const FiniteElement<dim> &finite_element,
424  const bool field_is_discontinuous);
425 
430  FullMatrix<double> local_matrix;
431 
443  std::vector<FullMatrix<double>> local_matrices_int_ext;
444  std::vector<FullMatrix<double>> local_matrices_ext_int;
445  std::vector<FullMatrix<double>> local_matrices_ext_ext;
446 
451  Vector<double> local_rhs;
452 
459  std::vector<bool> assembled_matrices;
460 
470  std::vector<types::global_dof_index> local_dof_indices;
471 
480  std::vector<std::vector<types::global_dof_index>> neighbor_dof_indices;
481  };
482  }
483  }
484  }
485 
490  namespace Assemblers
491  {
492 
502  unsigned int
503  n_interface_matrices (const ReferenceCell &reference_cell);
504 
510  unsigned int
511  nth_interface_matrix (const ReferenceCell &reference_cell,
512  const unsigned int face);
513 
519  unsigned int
520  nth_interface_matrix (const ReferenceCell &reference_cell,
521  const unsigned int face,
522  const unsigned int sub_face);
523 
538  template <int dim>
539  class Interface : public Plugins::InterfaceBase
540  {
541  public:
555  virtual
556  void
559 
593  virtual
594  void
595  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &) const;
596 
610  virtual
611  std::vector<double>
612  compute_residual(internal::Assembly::Scratch::ScratchBase<dim> &) const;
613  };
614 
615 
616 
624  template <int dim>
626  {
627  public:
629 
639  virtual
640  std::vector<double>
641  advection_prefactors(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const;
642 
652  virtual
653  std::vector<double>
654  diffusion_prefactors(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const;
655  };
656 
657 
658 
685  template <int dim>
686  class Manager
687  {
688  public:
689 
693  void reset ();
694 
699  std::vector<std::unique_ptr<Assemblers::Interface<dim>>> stokes_preconditioner;
700 
706  std::vector<std::unique_ptr<Assemblers::Interface<dim>>> stokes_system;
707 
715  std::vector<std::unique_ptr<Assemblers::Interface<dim>>> stokes_system_on_boundary_face;
716 
722  std::vector<std::vector<std::unique_ptr<Assemblers::Interface<dim>>>> advection_system;
723 
731  std::vector<std::vector<std::unique_ptr<Assemblers::Interface<dim>>>> advection_system_on_boundary_face;
732 
740  std::vector<std::vector<std::unique_ptr<Assemblers::Interface<dim>>>> advection_system_on_interior_face;
741 
761  struct Properties
762  {
767  Properties ();
768 
777 
786 
792 
799  UpdateFlags needed_update_flags;
800  };
801 
811  std::vector<Properties> advection_system_assembler_properties;
813  };
814  }
815 }
816 
817 
818 #endif
std::vector< SymmetricTensor< 2, dim > > old_strain_rates
Definition: interface.h:265
std::vector< Tensor< 1, dim > > old_field_grads
Definition: interface.h:273
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:470
std::vector< Tensor< 1, dim > > old_velocity_values
Definition: interface.h:257
std::vector< std::unique_ptr< Assemblers::Interface< dim > > > stokes_system
Definition: interface.h:706
std::vector< std::vector< std::unique_ptr< Assemblers::Interface< dim > > > > advection_system_on_interior_face
Definition: interface.h:740
Properties stokes_preconditioner_assembler_properties
Definition: interface.h:808
MaterialModel::MaterialModelInputs< dim > face_material_model_inputs
Definition: interface.h:188
std::vector< Tensor< 1, dim > > face_mesh_velocity_values
Definition: interface.h:285
std::vector< Tensor< 1, dim > > grad_phi_field
Definition: interface.h:250
std::vector< std::vector< std::unique_ptr< Assemblers::Interface< dim > > > > advection_system
Definition: interface.h:722
std::vector< Tensor< 1, dim > > face_grad_phi_field
Definition: interface.h:253
MaterialModel::MaterialModelOutputs< dim > material_model_outputs
Definition: interface.h:296
MaterialModel::MaterialModelInputs< dim > material_model_inputs
Definition: interface.h:295
HeatingModel::HeatingModelOutputs face_heating_model_outputs
Definition: interface.h:312
MaterialModel::MaterialModelOutputs< dim > neighbor_face_material_model_outputs
Definition: interface.h:302
std::vector< Tensor< 1, dim > > mesh_velocity_values
Definition: interface.h:284
const Simulator< dim >::AdvectionField * advection_field
Definition: interface.h:322
std::vector< Tensor< 1, dim > > old_pressure_gradients
Definition: interface.h:262
std::vector< SymmetricTensor< 2, dim > > old_old_strain_rates
Definition: interface.h:266
std::vector< Tensor< 1, dim > > old_old_field_grads
Definition: interface.h:274
unsigned int n_interface_matrices(const ReferenceCell &reference_cell)
std::vector< Tensor< 1, dim > > temperature_gradients
Definition: interface.h:178
std::vector< SymmetricTensor< 2, dim > > current_strain_rates
Definition: interface.h:287
ScratchBase(const ScratchBase &scratch)
Definition: interface.h:67
unsigned int nth_interface_matrix(const ReferenceCell &reference_cell, const unsigned int face, const unsigned int sub_face)
std::vector< Tensor< 1, dim > > current_velocity_values
Definition: interface.h:282
std::vector< Properties > advection_system_assembler_properties
Definition: interface.h:811
MaterialModel::MaterialModelOutputs< dim > face_material_model_outputs
Definition: interface.h:299
HeatingModel::HeatingModelOutputs neighbor_face_heating_model_outputs
Definition: interface.h:313
std::vector< Tensor< 1, dim > > neighbor_face_grad_phi_field
Definition: interface.h:255
std::unique_ptr< FEFaceValues< dim > > neighbor_face_finite_element_values
Definition: interface.h:236
std::vector< std::vector< types::global_dof_index > > neighbor_dof_indices
Definition: interface.h:480
std::vector< FullMatrix< double > > local_matrices_int_ext
Definition: interface.h:443
std::vector< Tensor< 1, dim > > old_old_velocity_values
Definition: interface.h:258
std::vector< FullMatrix< double > > local_matrices_ext_int
Definition: interface.h:444
std::vector< std::vector< double > > old_composition_values
Definition: interface.h:278
Definition: compat.h:59
std::vector< SymmetricTensor< 2, dim > > grads_phi_u
Definition: interface.h:115
Properties stokes_system_assembler_on_boundary_face_properties
Definition: interface.h:810
DoFHandler< dim >::active_cell_iterator cell
Definition: interface.h:78
std::vector< std::vector< double > > old_old_composition_values
Definition: interface.h:279
MaterialModel::MaterialModelOutputs< dim > face_material_model_outputs
Definition: interface.h:189
std::vector< Tensor< 1, dim > > old_old_pressure_gradients
Definition: interface.h:263
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:376
std::vector< Tensor< 1, dim > > velocity_values
Definition: interface.h:176
HeatingModel::HeatingModelOutputs heating_model_outputs
Definition: interface.h:311
std::vector< Properties > advection_system_assembler_on_face_properties
Definition: interface.h:812
std::vector< double > reference_densities_depth_derivative
Definition: interface.h:201
MaterialModel::MaterialModelInputs< dim > neighbor_face_material_model_inputs
Definition: interface.h:301
std::unique_ptr< FESubfaceValues< dim > > subface_finite_element_values
Definition: interface.h:237
std::vector< std::unique_ptr< Assemblers::Interface< dim > > > stokes_system_on_boundary_face
Definition: interface.h:715
std::vector< Tensor< 1, dim > > face_current_velocity_values
Definition: interface.h:283
std::vector< FullMatrix< double > > local_matrices_ext_ext
Definition: interface.h:445
std::vector< Tensor< 1, dim > > phi_u
Definition: interface.h:175
std::vector< std::vector< std::unique_ptr< Assemblers::Interface< dim > > > > advection_system_on_boundary_face
Definition: interface.h:731
MaterialModel::MaterialModelInputs< dim > face_material_model_inputs
Definition: interface.h:298
MaterialModel::MaterialModelInputs< dim > material_model_inputs
Definition: interface.h:126
std::vector< std::vector< double > > current_composition_values
Definition: interface.h:288
std::unique_ptr< FEFaceValues< dim > > face_finite_element_values
Definition: interface.h:235
Definition: compat.h:42
std::vector< std::unique_ptr< Assemblers::Interface< dim > > > stokes_preconditioner
Definition: interface.h:699
MaterialModel::MaterialModelOutputs< dim > material_model_outputs
Definition: interface.h:127
Properties stokes_system_assembler_properties
Definition: interface.h:809
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:239
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:113