ASPECT
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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  namespace MaterialModel
34  {
40  template <int dim>
42  {
43  public:
52  MeltInputs (const unsigned int n_points);
53 
57  std::vector<double> compaction_pressures;
58 
63  std::vector<Tensor<1,dim>> fluid_velocities;
64 
68  void fill (const LinearAlgebra::BlockVector &solution,
69  const FEValuesBase<dim> &fe_values,
70  const Introspection<dim> &introspection) override;
71  };
72 
73  template <int dim>
75  {
76  public:
85  MeltOutputs (const unsigned int n_points,
86  const unsigned int n_comp);
87 
93  std::vector<double> compaction_viscosities;
94 
98  std::vector<double> fluid_viscosities;
99 
103  std::vector<double> permeabilities;
104 
108  std::vector<double> fluid_densities;
109 
116  std::vector<Tensor<1,dim>> fluid_density_gradients;
117 
123  void average (const MaterialAveraging::AveragingOperation operation,
124  const FullMatrix<double> &projection_matrix,
125  const FullMatrix<double> &expansion_matrix) override;
126  };
127 
139  template <int dim>
141  {
142  public:
147  virtual ~MeltFractionModel () = default;
148 
157  virtual void melt_fractions (const MaterialModel::MaterialModelInputs<dim> &in,
158  std::vector<double> &melt_fractions) const = 0;
159 
167  template <typename ModelType>
168  static
169  bool is_melt_fraction_model (const ModelType &model_object);
170 
176  template <typename ModelType>
177  static
178  const MeltFractionModel<dim> &
179  as_melt_fraction_model (const ModelType &model_object);
180  };
181 
182 
183 
184  template <int dim>
185  template <typename ModelType>
186  inline
187  bool
188  MeltFractionModel<dim>::is_melt_fraction_model (const ModelType &model_object)
189  {
190  return (dynamic_cast<const MeltFractionModel<dim>*>(&model_object)
191  != nullptr);
192  }
193 
194 
195  template <int dim>
196  template <typename ModelType>
197  inline
198  const MeltFractionModel<dim> &
199  MeltFractionModel<dim>::as_melt_fraction_model (const ModelType &model_object)
200  {
201  Assert (is_melt_fraction_model(model_object) == true,
202  ExcMessage ("This function can only be called for model objects "
203  "whose types are derived from MeltFractionModel."));
204 
205  return dynamic_cast<const MeltFractionModel<dim>&>(model_object);
206  }
207 
208 
209 
213  template <int dim>
215  {
216  public:
221  virtual double reference_darcy_coefficient () const = 0;
222 
240  double p_c_scale (const MaterialModel::MaterialModelInputs<dim> &inputs,
242  const MeltHandler<dim> &melt_handler,
243  const bool consider_is_melt_cell) const;
244  };
245 
246 
247  }
248 
249 
250 
251  namespace Assemblers
252  {
258  template <int dim>
260  public SimulatorAccess<dim>
261  {
262  public:
268  void
269  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
270  };
271 
276  template <int dim>
278  {
279  public:
280  void
281  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
282  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
283  };
284 
289  template <int dim>
290  class MeltStokesSystem : public MeltInterface<dim>
291  {
292  public:
293  void
294  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
295  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
296  };
297 
298 
304  template <int dim>
306  {
307  public:
308  void
309  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
310  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
311  };
312 
317  template <int dim>
319  {
320  public:
321  void
322  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
323  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
324 
329  std::vector<double>
330  compute_residual(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base) const override;
331  };
332 
338  template <int dim>
340  {
341  public:
342  void
343  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
344  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
345  };
346 
350  template <int dim>
352  {
353  public:
354  void
355  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
356  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
357  };
358  }
359 
360 
361  namespace Melt
362  {
363  template <int dim>
364  struct Parameters
365  {
370  static void declare_parameters (ParameterHandler &prm);
371 
381  void parse_parameters (ParameterHandler &prm);
382 
394 
404 
409 
419  };
420  }
421 
422 
429  template <int dim>
430  class MeltHandler: public SimulatorAccess<dim>
431  {
432  public:
433  MeltHandler(ParameterHandler &prm);
434 
441  static void create_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &output);
442 
448  void edit_finite_element_variables(const Parameters<dim> &parameters,
449  std::vector<VariableDeclaration<dim>> &variables);
450 
456  void set_assemblers (Assemblers::Manager<dim> &assemblers) const;
457 
458 
465  void initialize() const;
466 
470  void initialize_simulator (const Simulator<dim> &simulator_object) override;
471 
486  void compute_melt_variables(LinearAlgebra::BlockSparseMatrix &system_matrix,
487  LinearAlgebra::BlockVector &solution,
488  LinearAlgebra::BlockVector &system_rhs) const;
489 
493  bool is_porosity (const typename Simulator<dim>::AdvectionField &advection_field) const;
494 
499  void apply_free_surface_stabilization_with_melt (const double free_surface_theta,
500  const typename DoFHandler<dim>::active_cell_iterator &cell,
503 
510  void add_current_constraints(AffineConstraints<double> &constraints);
511 
517  bool is_melt_cell(const typename DoFHandler<dim>::active_cell_iterator &cell) const;
518 
526  double limited_darcy_coefficient(const double K_D,
527  const bool is_melt_cell) const;
528 
533  get_boundary_fluid_pressure () const;
534 
540 
541  private:
547  const std::unique_ptr<aspect::BoundaryFluidPressure::Interface<dim>> boundary_fluid_pressure;
548 
556  std::vector<bool> is_melt_cell_vector;
557 
564  AffineConstraints<double> current_constraints;
565 
566  };
567 
568 }
569 
570 #endif
bool average_melt_velocity
Definition: melt.h:418
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:268
std::vector< bool > is_melt_cell_vector
Definition: melt.h:556
Melt::Parameters< dim > melt_parameters
Definition: melt.h:539
void declare_parameters(ParameterHandler &prm)
static const MeltFractionModel< dim > & as_melt_fraction_model(const ModelType &model_object)
Definition: melt.h:199
std::vector< Tensor< 1, dim > > fluid_density_gradients
Definition: melt.h:116
::TrilinosWrappers::BlockSparseMatrix BlockSparseMatrix
Definition: global.h:279
bool heat_advection_by_melt
Definition: melt.h:403
std::vector< Tensor< 1, dim > > fluid_velocities
Definition: melt.h:63
std::vector< double > fluid_viscosities
Definition: melt.h:98
std::vector< double > compaction_viscosities
Definition: melt.h:93
std::vector< double > permeabilities
Definition: melt.h:103
Definition: compat.h:59
double melt_scaling_factor_threshold
Definition: melt.h:393
const std::unique_ptr< aspect::BoundaryFluidPressure::Interface< dim > > boundary_fluid_pressure
Definition: melt.h:547
MeltInputs(const unsigned int n_points)
void fill(const LinearAlgebra::BlockVector &solution, const FEValuesBase< dim > &fe_values, const Introspection< dim > &introspection) override
std::vector< double > compaction_pressures
Definition: melt.h:57
static bool is_melt_fraction_model(const ModelType &model_object)
Definition: melt.h:188
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:108
bool use_discontinuous_p_c
Definition: melt.h:408
AffineConstraints< double > current_constraints
Definition: melt.h:564