https://mooseframework.inl.gov
ExplicitMixedOrder.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "ExplicitTimeIntegrator.h"
13 
14 // Forward declarations
15 namespace libMesh
16 {
17 template <typename T>
18 class SparseMatrix;
19 }
20 
26 {
27 public:
29 
31 
32  virtual int order() override { return 1; }
33  virtual void computeTimeDerivatives() override;
34 
35  virtual void solve() override;
36  virtual void postResidual(NumericVector<Number> & residual) override;
37  virtual bool overridesSolve() const override { return true; }
38 
39  virtual void postSolve() override
40  { // Once we have the new solution, we want to adanceState to make sure the
41  // coupling between the solution and the computed material properties is kept correctly.
43  }
44  virtual bool advancesProblemState() const override { return true; }
45 
46  virtual bool performExplicitSolve(SparseMatrix<Number> & mass_matrix) override;
47 
48  void computeADTimeDerivatives(ADReal &, const dof_id_type &, ADReal &) const override
49  {
50  mooseError("NOT SUPPORTED");
51  }
52  virtual void init() override;
53 
54  enum TimeOrder
55  {
58  };
59 
64  TimeOrder findVariableTimeOrder(unsigned int var_num) const;
65 
66 protected:
67  virtual TagID massMatrixTagID() const override;
68 
70  const bool & _constant_mass;
71 
73  const TagName & _mass_matrix;
74 
77 
78  // Variables that forward Euler time integration will be used for
79  std::unordered_set<unsigned int> & _vars_first;
80 
81  // local dofs that will have forward euler time integration
82  std::vector<dof_id_type> & _local_first_order_indices;
83 
84  // Variables that central difference time integration will be used for
85  std::unordered_set<unsigned int> & _vars_second;
86 
87  // local dofs that will have central difference time integration
88  std::vector<dof_id_type> & _local_second_order_indices;
89 
93  template <typename T, typename T2, typename T3, typename T4>
94  void
95  computeTimeDerivativeHelper(T & u_dot, T2 & u_dotdot, const T3 & u_old, const T4 & u_older) const;
96 
97  void computeICs();
98 };
99 
100 template <typename T, typename T2, typename T3, typename T4>
101 void
103  T2 & u_dotdot,
104  const T3 & u_old,
105  const T4 & u_older) const
106 {
107  // computing first derivative
108  // using the Central Difference method
109  // u_dot_old = (first_term - second_term) / 2 / dt
110  // first_term = u
111  // second_term = u_older
112  u_dot -= u_older; // 'older than older' solution
113  u_dot *= 1.0 / (2.0 * _dt);
114 
115  // computing second derivative
116  // using the Central Difference method
117  // u_dotdot_old = (first_term - second_term + third_term) / dt / dt
118  // first_term = u
119  // second_term = 2 * u_old
120  // third_term = u_older
121  u_dotdot -= u_old;
122  u_dotdot -= u_old;
123  u_dotdot += u_older;
124  u_dotdot *= 1.0 / (_dt * _dt);
125 }
virtual TagID massMatrixTagID() const override
unsigned int TagID
std::unordered_set< unsigned int > & _vars_second
FEProblemBase & _fe_problem
const NumericVector< Number > & _solution_older
The older solution.
const bool & _constant_mass
Whether we are reusing the mass matrix.
virtual void postSolve() override
const TagName & _mass_matrix
Mass matrix name.
TimeOrder findVariableTimeOrder(unsigned int var_num) const
Retrieve the order of the highest time derivative of a variable.
virtual void postResidual(NumericVector< Number > &residual) override
Implements a form of the central difference time integrator that calculates acceleration directly fro...
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
virtual void advanceState()
virtual bool overridesSolve() const override
std::unordered_set< unsigned int > & _vars_first
virtual bool advancesProblemState() const override
void computeTimeDerivativeHelper(T &u_dot, T2 &u_dotdot, const T3 &u_old, const T4 &u_older) const
Helper function that actually does the math for computing the time derivative.
std::vector< dof_id_type > & _local_second_order_indices
static InputParameters validParams()
virtual void computeTimeDerivatives() override
virtual int order() override
void computeADTimeDerivatives(ADReal &, const dof_id_type &, ADReal &) const override
virtual void init() override
virtual bool performExplicitSolve(SparseMatrix< Number > &mass_matrix) override
void mooseError(Args &&... args) const
const InputParameters & parameters() const
virtual void solve() override
ExplicitMixedOrder(const InputParameters &parameters)
std::vector< dof_id_type > & _local_first_order_indices
uint8_t dof_id_type