ASPECT
newton.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_newton_h
23 #define _aspect_newton_h
24 
27 #include <aspect/global.h>
29 
30 #include <deal.II/base/signaling_nan.h>
31 
32 namespace aspect
33 {
34  namespace MaterialModel
35  {
39  template <int dim>
41  {
42  public:
47  MaterialModelDerivatives (const unsigned int n_points);
48 
52  std::vector<double> viscosity_derivative_wrt_pressure;
53 
57  std::vector<SymmetricTensor<2,dim>> viscosity_derivative_wrt_strain_rate;
58 
64  };
65  }
66 
67 
68  namespace Newton
69  {
70  struct Parameters
71  {
72 
82  {
83  none = 0,
84  symmetric = 1,
85  PD = 2,
86  SPD = symmetric | PD
87  };
88 
89 
95  friend
98  const Stabilization b)
99  {
100  return static_cast<Stabilization>(
101  static_cast<int>(a) | static_cast<int>(b));
102  }
103 
104 
110  friend
112  operator& (const Stabilization a,
113  const Stabilization b)
114  {
115  return static_cast<Stabilization>(
116  static_cast<int>(a) & static_cast<int>(b));
117  }
118 
119 
124  static void declare_parameters (ParameterHandler &prm);
125 
130  void parse_parameters (ParameterHandler &prm);
131 
148 
151 
158 
165 
172  };
173 
174 
179  std::string
180  to_string(const Newton::Parameters::Stabilization preconditioner_stabilization);
181  }
182 
183 
184 
188  template <int dim>
189  class NewtonHandler : public SimulatorAccess<dim>
190  {
191  public:
197  void set_assemblers (Assemblers::Manager<dim> &assemblers) const;
198 
204  static void create_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &output);
205 
211  };
212 
213 
214  namespace Assemblers
215  {
220  template <int dim>
222  public SimulatorAccess<dim>
223  {
224  public:
230  void
231  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
232  };
233 
237  template <int dim>
239  {
240  public:
241  void
242  execute (internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
243  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
244 
249  void create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
250  };
251 
256  template <int dim>
258  {
259  public:
260  void
261  execute (internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
262  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
263 
269  void create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
270  };
271 
276  template <int dim>
278  public SimulatorAccess<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 
292  template <int dim>
294  public SimulatorAccess<dim>
295  {
296  public:
297  void
298  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
299  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
300  };
301 
310  template <int dim>
312  public SimulatorAccess<dim>
313  {
314  public:
315  void
316  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
317  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
318  };
319 
326  template <int dim>
328  public SimulatorAccess<dim>
329  {
330  public:
331  void
332  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
333  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
334  };
335 
344  template <int dim>
346  public SimulatorAccess<dim>
347  {
348  public:
349  void
352  };
353  }
354 }
355 
356 #endif
unsigned int max_pre_newton_nonlinear_iterations
Definition: newton.h:167
bool use_newton_residual_scaling_method
Definition: newton.h:169
double maximum_linear_stokes_solver_tolerance
Definition: newton.h:170
double nonlinear_switch_tolerance
Definition: newton.h:164
std::string to_string(const Newton::Parameters::Stabilization preconditioner_stabilization)
bool use_Eisenstat_Walker_method_for_Picard_iterations
Definition: newton.h:166
void declare_parameters(ParameterHandler &prm)
Stabilization preconditioner_stabilization
Definition: newton.h:149
std::vector< SymmetricTensor< 2, dim > > viscosity_derivative_wrt_strain_rate
Definition: newton.h:57
Stabilization velocity_block_stabilization
Definition: newton.h:150
double newton_derivative_scaling_factor
Definition: newton.h:147
Newton::Parameters parameters
Definition: newton.h:210
MaterialModelDerivatives(const unsigned int n_points)
Definition: compat.h:59
std::vector< double > viscosity_derivative_averaging_weights
Definition: newton.h:63
Dependence operator|(const Dependence d1, const Dependence d2)
Definition: interface.h:99
unsigned int max_newton_line_search_iterations
Definition: newton.h:168
std::vector< double > viscosity_derivative_wrt_pressure
Definition: newton.h:52