Line data Source code
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 : 21 : /** 22 : * Implements a form of the central difference time integrator that calculates acceleration directly 23 : * from the residual forces. 24 : */ 25 : class ExplicitMixedOrder : public ExplicitTimeIntegrator 26 : { 27 : public: 28 : static InputParameters validParams(); 29 : 30 : ExplicitMixedOrder(const InputParameters & parameters); 31 : 32 0 : 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 738 : virtual bool overridesSolve() const override { return true; } 38 : 39 369 : 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. 42 369 : _fe_problem.advanceState(); 43 369 : } 44 369 : virtual bool advancesProblemState() const override { return true; } 45 : 46 : virtual bool performExplicitSolve(SparseMatrix<Number> & mass_matrix) override; 47 : 48 0 : void computeADTimeDerivatives(ADReal &, const dof_id_type &, ADReal &) const override 49 : { 50 0 : mooseError("NOT SUPPORTED"); 51 : } 52 : virtual void init() override; 53 : 54 : virtual void meshChanged() override; 55 : 56 : enum TimeOrder 57 : { 58 : FIRST, 59 : SECOND 60 : }; 61 : 62 : /** 63 : * Retrieve the order of the highest time derivative of a variable. 64 : * @return Returns the time order enum of this variable. 65 : */ 66 : TimeOrder findVariableTimeOrder(unsigned int var_num) const; 67 : 68 : protected: 69 : /// compile the dof indices for first and second order in time variables 70 : void updateDOFIndices(); 71 : 72 : virtual TagID massMatrixTagID() const override; 73 : 74 : /// Evaluate the RHS residual 75 : virtual void evaluateRHSResidual(); 76 : 77 : /// Whether we are reusing the mass matrix 78 : const bool & _constant_mass; 79 : 80 : /** 81 : * Must be set to true to use adaptivity with a constant mass 82 : * matrix. This will recompute the mass matrix when the mesh changes. The user must make sure that 83 : * the underlying material density stays constant, otherwise simulation results will depend on 84 : * adaptivity. 85 : */ 86 : const bool & _recompute_mass_matrix_on_mesh_change; 87 : 88 : /// Whether the mesh changed just before the current solve 89 : bool _mesh_changed; 90 : 91 : /// Mass matrix name 92 : const TagName & _mass_matrix_name; 93 : 94 : /// Lumped mass matrix 95 : NumericVector<Real> * _mass_matrix_lumped; 96 : 97 : /// The older solution 98 : const NumericVector<Number> & _solution_older; 99 : 100 : // Variables that forward Euler time integration will be used for 101 : std::unordered_set<unsigned int> & _vars_first; 102 : 103 : // local dofs that will have forward euler time integration 104 : std::vector<dof_id_type> & _local_first_order_indices; 105 : 106 : // Variables that central difference time integration will be used for 107 : std::unordered_set<unsigned int> & _vars_second; 108 : 109 : // local dofs that will have central difference time integration 110 : std::vector<dof_id_type> & _local_second_order_indices; 111 : 112 : /** 113 : * Helper function that actually does the math for computing the time derivative 114 : */ 115 : template <typename T, typename T2, typename T3, typename T4> 116 : void 117 : computeTimeDerivativeHelper(T & u_dot, T2 & u_dotdot, const T3 & u_old, const T4 & u_older) const; 118 : 119 : void computeICs(); 120 : }; 121 : 122 : template <typename T, typename T2, typename T3, typename T4> 123 : void 124 : ExplicitMixedOrder::computeTimeDerivativeHelper(T & u_dot, 125 : T2 & u_dotdot, 126 : const T3 & u_old, 127 : const T4 & u_older) const 128 : { 129 : // computing first derivative 130 : // using the Central Difference method 131 : // u_dot_old = (first_term - second_term) / 2 / dt 132 : // first_term = u 133 : // second_term = u_older 134 : u_dot -= u_older; // 'older than older' solution 135 : u_dot *= 1.0 / (2.0 * _dt); 136 : 137 : // computing second derivative 138 : // using the Central Difference method 139 : // u_dotdot_old = (first_term - second_term + third_term) / dt / dt 140 : // first_term = u 141 : // second_term = 2 * u_old 142 : // third_term = u_older 143 : u_dotdot -= u_old; 144 : u_dotdot -= u_old; 145 : u_dotdot += u_older; 146 : u_dotdot *= 1.0 / (_dt * _dt); 147 : }