ASPECT
melt.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 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:
87  MeltOutputs (const unsigned int n_points,
88  const unsigned int n_comp);
89 
95  std::vector<double> compaction_viscosities;
96 
100  std::vector<double> fluid_viscosities;
101 
105  std::vector<double> permeabilities;
106 
110  std::vector<double> fluid_densities;
111 
118  std::vector<Tensor<1,dim>> fluid_density_gradients;
119 
125  void average (const MaterialAveraging::AveragingOperation operation,
126  const FullMatrix<double> &projection_matrix,
127  const FullMatrix<double> &expansion_matrix) override;
128  };
129 
141  template <int dim>
143  {
144  public:
149  virtual ~MeltFractionModel () = default;
150 
159  virtual void melt_fractions (const MaterialModel::MaterialModelInputs<dim> &in,
160  std::vector<double> &melt_fractions) const = 0;
161 
169  template <typename ModelType>
170  static
171  bool is_melt_fraction_model (const ModelType &model_object);
172 
178  template <typename ModelType>
179  static
180  const MeltFractionModel<dim> &
181  as_melt_fraction_model (const ModelType &model_object);
182  };
183 
184 
185 
186  template <int dim>
187  template <typename ModelType>
188  inline
189  bool
190  MeltFractionModel<dim>::is_melt_fraction_model (const ModelType &model_object)
191  {
192  return (dynamic_cast<const MeltFractionModel<dim>*>(&model_object)
193  != nullptr);
194  }
195 
196 
197  template <int dim>
198  template <typename ModelType>
199  inline
200  const MeltFractionModel<dim> &
201  MeltFractionModel<dim>::as_melt_fraction_model (const ModelType &model_object)
202  {
203  Assert (is_melt_fraction_model(model_object) == true,
204  ExcMessage ("This function can only be called for model objects "
205  "whose types are derived from MeltFractionModel."));
206 
207  return dynamic_cast<const MeltFractionModel<dim>&>(model_object);
208  }
209 
210 
211 
215  template <int dim>
217  {
218  public:
223  virtual double reference_darcy_coefficient () const = 0;
224 
242  double p_c_scale (const MaterialModel::MaterialModelInputs<dim> &inputs,
244  const MeltHandler<dim> &melt_handler,
245  const bool consider_is_melt_cell) const;
246  };
247 
248 
249  }
250 
251 
252 
253  namespace Assemblers
254  {
260  template <int dim>
262  public SimulatorAccess<dim>
263  {
264  public:
270  void
271  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
272  };
273 
278  template <int dim>
280  {
281  public:
282  void
283  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
284  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
285  };
286 
291  template <int dim>
292  class MeltStokesSystem : public MeltInterface<dim>
293  {
294  public:
295  void
296  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
297  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
298  };
299 
300 
306  template <int dim>
308  {
309  public:
310  void
311  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
312  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
313  };
314 
319  template <int dim>
321  {
322  public:
323  void
324  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
325  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
326 
331  std::vector<double>
332  compute_residual(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const override;
333  };
334 
340  template <int dim>
342  {
343  public:
344  void
345  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
346  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
347  };
348 
352  template <int dim>
354  {
355  public:
356  void
357  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
358  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
359  };
360  }
361 
362 
363  namespace Melt
364  {
365  template <int dim>
366  struct Parameters
367  {
372  static void declare_parameters (ParameterHandler &prm);
373 
383  void parse_parameters (ParameterHandler &prm);
384 
396 
406 
411 
421  };
422  }
423 
424 
431  template <int dim>
432  class MeltHandler: public SimulatorAccess<dim>
433  {
434  public:
435  MeltHandler(ParameterHandler &prm);
436 
443  static void create_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &output);
444 
450  void edit_finite_element_variables(const Parameters<dim> &parameters,
451  std::vector<VariableDeclaration<dim>> &variables);
452 
458  void set_assemblers (Assemblers::Manager<dim> &assemblers) const;
459 
460 
467  void initialize() const;
468 
472  void initialize_simulator (const Simulator<dim> &simulator_object) override;
473 
488  void compute_melt_variables(LinearAlgebra::BlockSparseMatrix &system_matrix,
489  LinearAlgebra::BlockVector &solution,
490  LinearAlgebra::BlockVector &system_rhs) const;
491 
495  bool is_porosity (const typename Simulator<dim>::AdvectionField &advection_field) const;
496 
501  void apply_free_surface_stabilization_with_melt (const double free_surface_theta,
502  const typename DoFHandler<dim>::active_cell_iterator &cell,
505 
512  void add_current_constraints(AffineConstraints<double> &constraints);
513 
519  bool is_melt_cell(const typename DoFHandler<dim>::active_cell_iterator &cell) const;
520 
528  double limited_darcy_coefficient(const double K_D,
529  const bool is_melt_cell) const;
530 
535  get_boundary_fluid_pressure () const;
536 
542 
543  private:
549  const std::unique_ptr<aspect::BoundaryFluidPressure::Interface<dim>> boundary_fluid_pressure;
550 
558  std::vector<bool> is_melt_cell_vector;
559 
566  AffineConstraints<double> current_constraints;
567 
568  };
569 
570 }
571 
572 #endif
bool average_melt_velocity
Definition: melt.h:420
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:261
std::vector< bool > is_melt_cell_vector
Definition: melt.h:558
Melt::Parameters< dim > melt_parameters
Definition: melt.h:541
std::vector< Tensor< 1, dim > > fluid_density_gradients
Definition: melt.h:118
::TrilinosWrappers::BlockSparseMatrix BlockSparseMatrix
Definition: global.h:272
void declare_parameters(ParameterHandler &prm)
bool heat_advection_by_melt
Definition: melt.h:405
std::vector< Tensor< 1, dim > > fluid_velocities
Definition: melt.h:65
std::vector< double > fluid_viscosities
Definition: melt.h:100
std::vector< double > compaction_viscosities
Definition: melt.h:95
std::vector< double > permeabilities
Definition: melt.h:105
Definition: compat.h:59
double melt_scaling_factor_threshold
Definition: melt.h:395
const std::unique_ptr< aspect::BoundaryFluidPressure::Interface< dim > > boundary_fluid_pressure
Definition: melt.h:549
std::vector< double > compaction_pressures
Definition: melt.h:59
void average(const AveragingOperation operation, const typename DoFHandler< dim >::active_cell_iterator &cell, const Quadrature< dim > &quadrature_formula, const Mapping< dim > &mapping, const MaterialProperties::Property &requested_properties, MaterialModelOutputs< dim > &values_out)
std::vector< double > fluid_densities
Definition: melt.h:110
bool use_discontinuous_p_c
Definition: melt.h:410
Definition: compat.h:42
AffineConstraints< double > current_constraints
Definition: melt.h:566