ASPECT
newton.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2016 - 2022 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  using namespace dealii;
35 
36  namespace MaterialModel
37  {
41  template <int dim>
43  {
44  public:
49  MaterialModelDerivatives (const unsigned int n_points);
50 
54  std::vector<double> viscosity_derivative_wrt_pressure;
55 
59  std::vector<SymmetricTensor<2,dim>> viscosity_derivative_wrt_strain_rate;
60 
66  };
67  }
68 
69 
70  namespace Newton
71  {
72  struct Parameters
73  {
74 
84  {
85  none = 0,
86  symmetric = 1,
87  PD = 2,
88  SPD = symmetric | PD
89  };
90 
91 
97  friend
100  const Stabilization b)
101  {
102  return static_cast<Stabilization>(
103  static_cast<int>(a) | static_cast<int>(b));
104  }
105 
106 
112  friend
114  operator& (const Stabilization a,
115  const Stabilization b)
116  {
117  return static_cast<Stabilization>(
118  static_cast<int>(a) & static_cast<int>(b));
119  }
120 
121 
126  static void declare_parameters (ParameterHandler &prm);
127 
132  void parse_parameters (ParameterHandler &prm);
133 
150 
153 
160 
167 
174  };
175 
176 
181  std::string
182  to_string(const Newton::Parameters::Stabilization preconditioner_stabilization);
183  }
184 
185 
186 
190  template <int dim>
191  class NewtonHandler : public SimulatorAccess<dim>
192  {
193  public:
199  void set_assemblers (Assemblers::Manager<dim> &assemblers) const;
200 
206  static void create_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &output);
207 
213  };
214 
215 
216  namespace Assemblers
217  {
222  template <int dim>
224  public SimulatorAccess<dim>
225  {
226  public:
232  void
233  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
234  };
235 
239  template <int dim>
241  {
242  public:
243  void
244  execute (internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
245  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
246 
251  void create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
252  };
253 
258  template <int dim>
260  {
261  public:
262  void
263  execute (internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
264  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
265 
271  void create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const override;
272  };
273 
278  template <int dim>
280  public SimulatorAccess<dim>
281  {
282  public:
283  void
284  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
285  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
286  };
287 
294  template <int dim>
296  public SimulatorAccess<dim>
297  {
298  public:
299  void
300  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
301  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
302  };
303 
312  template <int dim>
314  public SimulatorAccess<dim>
315  {
316  public:
317  void
318  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
319  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
320  };
321 
328  template <int dim>
330  public SimulatorAccess<dim>
331  {
332  public:
333  void
334  execute(internal::Assembly::Scratch::ScratchBase<dim> &scratch_base,
335  internal::Assembly::CopyData::CopyDataBase<dim> &data_base) const override;
336  };
337 
346  template <int dim>
348  public SimulatorAccess<dim>
349  {
350  public:
351  void
354  };
355  }
356 }
357 
358 #endif
unsigned int max_pre_newton_nonlinear_iterations
Definition: newton.h:169
bool use_newton_residual_scaling_method
Definition: newton.h:171
double maximum_linear_stokes_solver_tolerance
Definition: newton.h:172
double nonlinear_switch_tolerance
Definition: newton.h:166
std::string to_string(const Newton::Parameters::Stabilization preconditioner_stabilization)
bool use_Eisenstat_Walker_method_for_Picard_iterations
Definition: newton.h:168
Stabilization preconditioner_stabilization
Definition: newton.h:151
std::vector< SymmetricTensor< 2, dim > > viscosity_derivative_wrt_strain_rate
Definition: newton.h:59
void declare_parameters(ParameterHandler &prm)
Stabilization velocity_block_stabilization
Definition: newton.h:152
double newton_derivative_scaling_factor
Definition: newton.h:149
Newton::Parameters parameters
Definition: newton.h:212
std::vector< double > viscosity_derivative_averaging_weights
Definition: newton.h:65
Dependence operator|(const Dependence d1, const Dependence d2)
Definition: interface.h:101
Definition: compat.h:42
unsigned int max_newton_line_search_iterations
Definition: newton.h:170
std::vector< double > viscosity_derivative_wrt_pressure
Definition: newton.h:54