ASPECT
melt.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2016 - 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 LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19  */
20 
21 
22 #ifndef _aspect_melt_h
23 #define _aspect_melt_h
24 
26 #include <aspect/global.h>
30 
31 namespace aspect
32 {
33  using namespace dealii;
34 
35  namespace MaterialModel
36  {
42  template <int dim>
44  {
45  public:
54  MeltInputs (const unsigned int n_points);
55 
59  std::vector<double> compaction_pressures;
60 
65  std::vector<Tensor<1,dim>> fluid_velocities;
66 
70  void fill (const LinearAlgebra::BlockVector &solution,
71  const FEValuesBase<dim> &fe_values,
72  const Introspection<dim> &introspection) override;
73  };
74 
75  template <int dim>
77  {
78  public:
79  MeltOutputs (const unsigned int n_points,
80  const unsigned int /*n_comp*/)
81  {
82  compaction_viscosities.resize(n_points);
83  fluid_viscosities.resize(n_points);
84  permeabilities.resize(n_points);
85  fluid_densities.resize(n_points);
86  fluid_density_gradients.resize(n_points, Tensor<1,dim>());
87  }
88 
94  std::vector<double> compaction_viscosities;
95 
99  std::vector<double> fluid_viscosities;
100 
104  std::vector<double> permeabilities;
105 
109  std::vector<double> fluid_densities;
110 
117  std::vector<Tensor<1,dim>> fluid_density_gradients;
118 
124  void average (const MaterialAveraging::AveragingOperation operation,
125  const FullMatrix<double> &projection_matrix,
126  const FullMatrix<double> &expansion_matrix) override;
127  };
128 
140  template <int dim>
142  {
143  public:
148  virtual ~MeltFractionModel () = default;
149 
158  virtual void melt_fractions (const MaterialModel::MaterialModelInputs<dim> &in,
159  std::vector<double> &melt_fractions) const = 0;
160 
168  template <typename ModelType>
169  static
170  bool is_melt_fraction_model (const ModelType &model_object);
171 
177  template <typename ModelType>
178  static
179  const MeltFractionModel<dim> &
180  as_melt_fraction_model (const ModelType &model_object);
181  };
182 
183 
184 
185  template <int dim>
186  template <typename ModelType>
187  inline
188  bool
189  MeltFractionModel<dim>::is_melt_fraction_model (const ModelType &model_object)
190  {
191  return (dynamic_cast<const MeltFractionModel<dim>*>(&model_object)
192  != nullptr);
193  }
194 
195 
196  template <int dim>
197  template <typename ModelType>
198  inline
199  const MeltFractionModel<dim> &
200  MeltFractionModel<dim>::as_melt_fraction_model (const ModelType &model_object)
201  {
202  Assert (is_melt_fraction_model(model_object) == true,
203  ExcMessage ("This function can only be called for model objects "
204  "whose types are derived from MeltFractionModel."));
205 
206  return dynamic_cast<const MeltFractionModel<dim>&>(model_object);
207  }
208 
209 
210 
214  template <int dim>
216  {
217  public:
222  virtual double reference_darcy_coefficient () const = 0;
223 
241  double p_c_scale (const MaterialModel::MaterialModelInputs<dim> &inputs,
243  const MeltHandler<dim> &melt_handler,
244  const bool consider_is_melt_cell) const;
245  };
246 
247 
248  }
249 
250 
251 
252  namespace Assemblers
253  {
259  template <int dim>
261  public SimulatorAccess<dim>
262  {
263  public:
269  void
270  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
271  };
272 
277  template <int dim>
279  {
280  public:
281  void
282  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
283  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
284  };
285 
290  template <int dim>
291  class MeltStokesSystem : public MeltInterface<dim>
292  {
293  public:
294  void
295  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
296  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
297  };
298 
299 
305  template <int dim>
307  {
308  public:
309  void
310  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
311  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
312  };
313 
318  template <int dim>
320  {
321  public:
322  void
323  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
324  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
325 
330  std::vector<double>
331  compute_residual(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const override;
332  };
333 
339  template <int dim>
341  {
342  public:
343  void
344  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
345  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
346  };
347 
351  template <int dim>
353  {
354  public:
355  void
356  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
357  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
358  };
359  }
360 
361 
362  namespace Melt
363  {
364  template <int dim>
365  struct Parameters
366  {
371  static void declare_parameters (ParameterHandler &prm);
372 
382  void parse_parameters (ParameterHandler &prm);
383 
395 
405 
410 
420  };
421  }
422 
423 
430  template <int dim>
431  class MeltHandler: public SimulatorAccess<dim>
432  {
433  public:
434  MeltHandler(ParameterHandler &prm);
435 
442  static void create_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &output);
443 
449  void edit_finite_element_variables(const Parameters<dim> &parameters,
450  std::vector<VariableDeclaration<dim>> &variables);
451 
457  void set_assemblers (Assemblers::Manager<dim> &assemblers) const;
458 
459 
466  void initialize() const;
467 
471  void initialize_simulator (const Simulator<dim> &simulator_object) override;
472 
487  void compute_melt_variables(LinearAlgebra::BlockSparseMatrix &system_matrix,
488  LinearAlgebra::BlockVector &solution,
489  LinearAlgebra::BlockVector &system_rhs);
490 
494  bool is_porosity (const typename Simulator<dim>::AdvectionField &advection_field) const;
495 
500  void apply_free_surface_stabilization_with_melt (const double free_surface_theta,
501  const typename DoFHandler<dim>::active_cell_iterator &cell,
504 
511  void add_current_constraints(AffineConstraints<double> &constraints);
512 
518  bool is_melt_cell(const typename DoFHandler<dim>::active_cell_iterator &cell) const;
519 
527  double limited_darcy_coefficient(const double K_D,
528  const bool is_melt_cell) const;
529 
534  get_boundary_fluid_pressure () const;
535 
541 
542  private:
548  const std::unique_ptr<aspect::BoundaryFluidPressure::Interface<dim>> boundary_fluid_pressure;
549 
557  std::vector<bool> is_melt_cell_vector;
558 
565  AffineConstraints<double> current_constraints;
566 
567  };
568 
569 }
570 
571 #endif
bool average_melt_velocity
Definition: melt.h:419
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:239
std::vector< bool > is_melt_cell_vector
Definition: melt.h:557
Melt::Parameters< dim > melt_parameters
Definition: melt.h:540
void average(const AveragingOperation operation, const typename DoFHandler< dim >::active_cell_iterator &cell, const Quadrature< dim > &quadrature_formula, const Mapping< dim > &mapping, MaterialModelOutputs< dim > &values_out)
std::vector< Tensor< 1, dim > > fluid_density_gradients
Definition: melt.h:117
::TrilinosWrappers::BlockSparseMatrix BlockSparseMatrix
Definition: global.h:250
void declare_parameters(ParameterHandler &prm)
bool heat_advection_by_melt
Definition: melt.h:404
std::vector< Tensor< 1, dim > > fluid_velocities
Definition: melt.h:65
std::vector< double > fluid_viscosities
Definition: melt.h:99
std::vector< double > compaction_viscosities
Definition: melt.h:94
std::vector< double > permeabilities
Definition: melt.h:104
double melt_scaling_factor_threshold
Definition: melt.h:394
const std::unique_ptr< aspect::BoundaryFluidPressure::Interface< dim > > boundary_fluid_pressure
Definition: melt.h:548
std::vector< double > compaction_pressures
Definition: melt.h:59
std::vector< double > fluid_densities
Definition: melt.h:109
bool use_discontinuous_p_c
Definition: melt.h:409
MeltOutputs(const unsigned int n_points, const unsigned int)
Definition: melt.h:79
Definition: compat.h:42
AffineConstraints< double > current_constraints
Definition: melt.h:565