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  template <int dim>
33  class Simulator;
34 
35  struct AdvectionField;
36 
44  namespace internal
45  {
46  namespace Assembly
47  {
48  namespace Scratch
49  {
56  template <int dim>
57  struct ScratchBase
58  {
60  :
61  cell(),
62  face_number(numbers::invalid_unsigned_int)
63  {}
64 
65  ScratchBase(const ScratchBase &scratch)
66  :
67  cell(scratch.cell),
68  face_number(scratch.face_number)
69  {}
70 
71  virtual ~ScratchBase () = default;
72 
76  typename DoFHandler<dim>::active_cell_iterator cell;
77 
84  unsigned face_number;
85  };
86 
91  template <int dim>
92  struct StokesPreconditioner: public ScratchBase<dim>
93  {
94  StokesPreconditioner (const FiniteElement<dim> &finite_element,
95  const Quadrature<dim> &quadrature,
96  const Mapping<dim> &mapping,
97  const UpdateFlags update_flags,
98  const unsigned int n_compositional_fields,
99  const unsigned int stokes_dofs_per_cell,
100  const bool add_compaction_pressure,
101  const bool rebuild_matrix,
102  const bool use_bfbt);
103  StokesPreconditioner (const StokesPreconditioner &scratch);
104 
105  ~StokesPreconditioner () override;
106 
107  FEValues<dim> finite_element_values;
108 
109  void reinit (const typename DoFHandler<dim>::active_cell_iterator &cell_ref);
110 
111  std::vector<types::global_dof_index> local_dof_indices;
112  std::vector<unsigned int> dof_component_indices;
113  std::vector<SymmetricTensor<2,dim>> grads_phi_u;
114  std::vector<double> div_phi_u;
115  std::vector<double> phi_p;
116  std::vector<Tensor<1,dim>> phi_u;
117  std::vector<double> phi_p_c;
118  std::vector<Tensor<1,dim>> grad_phi_p;
119 
126 
133  };
134 
135 
136 
147  template <int dim>
148  struct StokesSystem : public StokesPreconditioner<dim>
149  {
150  StokesSystem (const FiniteElement<dim> &finite_element,
151  const Mapping<dim> &mapping,
152  const Quadrature<dim> &quadrature,
153  const Quadrature<dim-1> &face_quadrature,
154  const UpdateFlags update_flags,
155  const UpdateFlags face_update_flags,
156  const unsigned int n_compositional_fields,
157  const unsigned int stokes_dofs_per_cell,
158  const bool add_compaction_pressure,
159  const bool use_reference_density_profile,
160  const bool rebuild_stokes_matrix,
161  const bool rebuild_newton_stokes_matrix,
162  const bool use_bfbt);
163 
164  StokesSystem (const StokesSystem<dim> &scratch);
165 
166  FEFaceValues<dim> face_finite_element_values;
167 
169 
170  void reinit (const typename DoFHandler<dim>::active_cell_iterator &cell_ref,
171  const unsigned face_number_ref);
172 
173  std::vector<Tensor<1,dim>> phi_u;
174  std::vector<Tensor<1,dim>> velocity_values;
175  std::vector<double> velocity_divergence;
176  std::vector<Tensor<1,dim>> temperature_gradients;
177 
188 
198  std::vector<double> reference_densities;
200 
207  };
208 
209 
210 
215  template <int dim>
216  struct AdvectionSystem: public ScratchBase<dim>
217  {
218  AdvectionSystem (const FiniteElement<dim> &finite_element,
219  const FiniteElement<dim> &advection_element,
220  const Mapping<dim> &mapping,
221  const Quadrature<dim> &quadrature,
222  const Quadrature<dim-1> &face_quadrature,
223  const UpdateFlags update_flags,
224  const UpdateFlags face_update_flags,
225  const unsigned int n_compositional_fields,
226  const typename Simulator<dim>::AdvectionField &field);
227  AdvectionSystem (const AdvectionSystem &scratch);
228 
229  FEValues<dim> finite_element_values;
230 
231  void reinit (const typename DoFHandler<dim>::active_cell_iterator &cell_ref);
232 
233  std::unique_ptr<FEFaceValues<dim>> face_finite_element_values;
234  std::unique_ptr<FEFaceValues<dim>> neighbor_face_finite_element_values;
235  std::unique_ptr<FESubfaceValues<dim>> subface_finite_element_values;
236 
237  std::vector<types::global_dof_index> local_dof_indices;
238 
247  std::vector<double> phi_field;
248  std::vector<Tensor<1,dim>> grad_phi_field;
249  std::vector<double> laplacian_phi_field;
250  std::vector<double> face_phi_field;
251  std::vector<Tensor<1,dim>> face_grad_phi_field;
252  std::vector<double> neighbor_face_phi_field;
253  std::vector<Tensor<1,dim>> neighbor_face_grad_phi_field;
254 
255  std::vector<Tensor<1,dim>> old_velocity_values;
256  std::vector<Tensor<1,dim>> old_old_velocity_values;
257 
258  std::vector<double> old_pressure;
259  std::vector<double> old_old_pressure;
260  std::vector<Tensor<1,dim>> old_pressure_gradients;
261  std::vector<Tensor<1,dim>> old_old_pressure_gradients;
262 
263  std::vector<SymmetricTensor<2,dim>> old_strain_rates;
264  std::vector<SymmetricTensor<2,dim>> old_old_strain_rates;
265 
266  std::vector<double> old_temperature_values;
267  std::vector<double> old_old_temperature_values;
268 
269  std::vector<double> old_field_values;
270  std::vector<double> old_old_field_values;
271  std::vector<Tensor<1,dim>> old_field_grads;
272  std::vector<Tensor<1,dim>> old_old_field_grads;
273  std::vector<double> old_field_laplacians;
274  std::vector<double> old_old_field_laplacians;
275 
276  std::vector<std::vector<double>> old_composition_values;
277  std::vector<std::vector<double>> old_old_composition_values;
278 
279  std::vector<double> current_temperature_values;
280  std::vector<Tensor<1,dim>> current_velocity_values;
281  std::vector<Tensor<1,dim>> face_current_velocity_values;
282  std::vector<Tensor<1,dim>> mesh_velocity_values;
283  std::vector<Tensor<1,dim>> face_mesh_velocity_values;
284 
285  std::vector<SymmetricTensor<2,dim>> current_strain_rates;
286  std::vector<std::vector<double>> current_composition_values;
287  std::vector<double> current_velocity_divergences;
288 
295 
298 
301 
312 
321 
327  };
328  }
329 
330 
331 
344  namespace CopyData
345  {
351  template <int dim>
353  {
354  virtual ~CopyDataBase () = default;
355  };
356 
362  template <int dim>
363  struct StokesPreconditioner: public CopyDataBase<dim>
364  {
365  StokesPreconditioner (const unsigned int stokes_dofs_per_cell);
366 
367  StokesPreconditioner (const StokesPreconditioner &data);
368 
369  ~StokesPreconditioner () override = default;
370  StokesPreconditioner<dim> &operator= (const StokesPreconditioner<dim> &data) = default;
371 
372  FullMatrix<double> local_matrix;
374  std::vector<types::global_dof_index> local_dof_indices;
375 
381  void extract_stokes_dof_indices(const std::vector<types::global_dof_index> &all_dof_indices,
382  const Introspection<dim> &introspection,
383  const FiniteElement<dim> &finite_element);
384  };
385 
391  template <int dim>
392  struct StokesSystem : public StokesPreconditioner<dim>
393  {
394  StokesSystem (const unsigned int stokes_dofs_per_cell,
395  const bool do_pressure_rhs_compatibility_modification);
396  StokesSystem (const StokesSystem<dim> &data);
397 
398  ~StokesSystem () override = default;
399  StokesSystem<dim> &operator= (const StokesSystem<dim> &data) = default;
400 
401  Vector<double> local_rhs;
403  };
404 
410  template <int dim>
411  struct AdvectionSystem: public CopyDataBase<dim>
412  {
421  AdvectionSystem (const FiniteElement<dim> &finite_element,
422  const bool field_is_discontinuous);
423 
428  FullMatrix<double> local_matrix;
429 
441  std::vector<FullMatrix<double>> local_matrices_int_ext;
442  std::vector<FullMatrix<double>> local_matrices_ext_int;
443  std::vector<FullMatrix<double>> local_matrices_ext_ext;
444 
449  Vector<double> local_rhs;
450 
457  std::vector<bool> assembled_matrices;
458 
468  std::vector<types::global_dof_index> local_dof_indices;
469 
478  std::vector<std::vector<types::global_dof_index>> neighbor_dof_indices;
479  };
480  }
481  }
482  }
483 
488  namespace Assemblers
489  {
490 
500  unsigned int
501  n_interface_matrices (const ReferenceCell &reference_cell);
502 
508  unsigned int
509  nth_interface_matrix (const ReferenceCell &reference_cell,
510  const unsigned int face);
511 
517  unsigned int
518  nth_interface_matrix (const ReferenceCell &reference_cell,
519  const unsigned int face,
520  const unsigned int sub_face);
521 
536  template <int dim>
537  class Interface : public Plugins::InterfaceBase
538  {
539  public:
553  virtual
554  void
557 
591  virtual
592  void
593  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &) const;
594 
608  virtual
609  std::vector<double>
610  compute_residual(internal::Assembly::Scratch::ScratchBase<dim> &) const;
611  };
612 
613 
614 
622  template <int dim>
624  {
625  public:
627 
637  virtual
638  std::vector<double>
639  advection_prefactors(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const;
640 
650  virtual
651  std::vector<double>
652  diffusion_prefactors(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const;
653  };
654 
655 
656 
683  template <int dim>
684  class Manager
685  {
686  public:
687 
691  void reset ();
692 
697  std::vector<std::unique_ptr<Assemblers::Interface<dim>>> stokes_preconditioner;
698 
704  std::vector<std::unique_ptr<Assemblers::Interface<dim>>> stokes_system;
705 
713  std::vector<std::unique_ptr<Assemblers::Interface<dim>>> stokes_system_on_boundary_face;
714 
720  std::vector<std::vector<std::unique_ptr<Assemblers::Interface<dim>>>> advection_system;
721 
729  std::vector<std::vector<std::unique_ptr<Assemblers::Interface<dim>>>> advection_system_on_boundary_face;
730 
738  std::vector<std::vector<std::unique_ptr<Assemblers::Interface<dim>>>> advection_system_on_interior_face;
739 
759  struct Properties
760  {
765  Properties ();
766 
775 
784 
790 
797  UpdateFlags needed_update_flags;
798  };
799 
809  std::vector<Properties> advection_system_assembler_properties;
811  };
812  }
813 }
814 
815 
816 #endif
std::vector< SymmetricTensor< 2, dim > > old_strain_rates
Definition: interface.h:263
std::vector< Tensor< 1, dim > > old_field_grads
Definition: interface.h:271
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:468
std::vector< Tensor< 1, dim > > old_velocity_values
Definition: interface.h:255
std::vector< std::unique_ptr< Assemblers::Interface< dim > > > stokes_system
Definition: interface.h:704
std::vector< std::vector< std::unique_ptr< Assemblers::Interface< dim > > > > advection_system_on_interior_face
Definition: interface.h:738
Properties stokes_preconditioner_assembler_properties
Definition: interface.h:806
MaterialModel::MaterialModelInputs< dim > face_material_model_inputs
Definition: interface.h:186
std::vector< Tensor< 1, dim > > face_mesh_velocity_values
Definition: interface.h:283
std::vector< Tensor< 1, dim > > grad_phi_field
Definition: interface.h:248
std::vector< std::vector< std::unique_ptr< Assemblers::Interface< dim > > > > advection_system
Definition: interface.h:720
std::vector< Tensor< 1, dim > > face_grad_phi_field
Definition: interface.h:251
MaterialModel::MaterialModelOutputs< dim > material_model_outputs
Definition: interface.h:294
MaterialModel::MaterialModelInputs< dim > material_model_inputs
Definition: interface.h:293
HeatingModel::HeatingModelOutputs face_heating_model_outputs
Definition: interface.h:310
MaterialModel::MaterialModelOutputs< dim > neighbor_face_material_model_outputs
Definition: interface.h:300
std::vector< Tensor< 1, dim > > mesh_velocity_values
Definition: interface.h:282
const Simulator< dim >::AdvectionField * advection_field
Definition: interface.h:320
std::vector< Tensor< 1, dim > > old_pressure_gradients
Definition: interface.h:260
std::vector< SymmetricTensor< 2, dim > > old_old_strain_rates
Definition: interface.h:264
std::vector< Tensor< 1, dim > > old_old_field_grads
Definition: interface.h:272
unsigned int n_interface_matrices(const ReferenceCell &reference_cell)
std::vector< Tensor< 1, dim > > temperature_gradients
Definition: interface.h:176
std::vector< SymmetricTensor< 2, dim > > current_strain_rates
Definition: interface.h:285
ScratchBase(const ScratchBase &scratch)
Definition: interface.h:65
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:280
std::vector< Properties > advection_system_assembler_properties
Definition: interface.h:809
MaterialModel::MaterialModelOutputs< dim > face_material_model_outputs
Definition: interface.h:297
HeatingModel::HeatingModelOutputs neighbor_face_heating_model_outputs
Definition: interface.h:311
std::vector< Tensor< 1, dim > > neighbor_face_grad_phi_field
Definition: interface.h:253
std::unique_ptr< FEFaceValues< dim > > neighbor_face_finite_element_values
Definition: interface.h:234
std::vector< std::vector< types::global_dof_index > > neighbor_dof_indices
Definition: interface.h:478
std::vector< FullMatrix< double > > local_matrices_int_ext
Definition: interface.h:441
std::vector< Tensor< 1, dim > > old_old_velocity_values
Definition: interface.h:256
std::vector< FullMatrix< double > > local_matrices_ext_int
Definition: interface.h:442
std::vector< std::vector< double > > old_composition_values
Definition: interface.h:276
Definition: compat.h:59
std::vector< SymmetricTensor< 2, dim > > grads_phi_u
Definition: interface.h:113
Properties stokes_system_assembler_on_boundary_face_properties
Definition: interface.h:808
DoFHandler< dim >::active_cell_iterator cell
Definition: interface.h:76
std::vector< std::vector< double > > old_old_composition_values
Definition: interface.h:277
MaterialModel::MaterialModelOutputs< dim > face_material_model_outputs
Definition: interface.h:187
std::vector< Tensor< 1, dim > > old_old_pressure_gradients
Definition: interface.h:261
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:374
std::vector< Tensor< 1, dim > > velocity_values
Definition: interface.h:174
HeatingModel::HeatingModelOutputs heating_model_outputs
Definition: interface.h:309
std::vector< Properties > advection_system_assembler_on_face_properties
Definition: interface.h:810
std::vector< double > reference_densities_depth_derivative
Definition: interface.h:199
MaterialModel::MaterialModelInputs< dim > neighbor_face_material_model_inputs
Definition: interface.h:299
std::unique_ptr< FESubfaceValues< dim > > subface_finite_element_values
Definition: interface.h:235
std::vector< std::unique_ptr< Assemblers::Interface< dim > > > stokes_system_on_boundary_face
Definition: interface.h:713
std::vector< Tensor< 1, dim > > face_current_velocity_values
Definition: interface.h:281
std::vector< FullMatrix< double > > local_matrices_ext_ext
Definition: interface.h:443
std::vector< Tensor< 1, dim > > phi_u
Definition: interface.h:173
std::vector< std::vector< std::unique_ptr< Assemblers::Interface< dim > > > > advection_system_on_boundary_face
Definition: interface.h:729
MaterialModel::MaterialModelInputs< dim > face_material_model_inputs
Definition: interface.h:296
MaterialModel::MaterialModelInputs< dim > material_model_inputs
Definition: interface.h:124
std::vector< std::vector< double > > current_composition_values
Definition: interface.h:286
std::unique_ptr< FEFaceValues< dim > > face_finite_element_values
Definition: interface.h:233
std::vector< std::unique_ptr< Assemblers::Interface< dim > > > stokes_preconditioner
Definition: interface.h:697
MaterialModel::MaterialModelOutputs< dim > material_model_outputs
Definition: interface.h:125
Properties stokes_system_assembler_properties
Definition: interface.h:807
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:237
std::vector< types::global_dof_index > local_dof_indices
Definition: interface.h:111