https://mooseframework.inl.gov
ExplicitDirichletBCBase.C
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 #include "ExplicitMixedOrder.h"
12 #include "MooseError.h"
13 #include "NonlinearSystemBase.h"
14 #include "libmesh/numeric_vector.h"
15 #include <iostream>
16 #include <ostream>
17 
20 {
22  return params;
23 }
24 
26  : NodalBC(parameters),
27  _mass_diag(initLumpedMass()),
28  _u_old(_var.nodalValueOld()),
29  _u_dot_old(_var.nodalValueDotOld()),
30  _explicit_integrator(
31  dynamic_cast<const ExplicitMixedOrder *>(&_sys.getTimeIntegrator(_var.number())))
32 {
34  mooseError("Time integrator for the variable is not of the right type.");
35 }
36 
37 Real
39 {
40  // Get dof for current var
41  const auto dofnum = _variable->nodalDofIndex();
42  Real resid = 0;
43  // Compute residual to enforce BC based on time order
44  switch (_var_time_order)
45  {
47  resid = (computeQpValue() - _u_old) / _dt;
48  resid /= -_mass_diag(dofnum);
49  break;
50 
52  Real avg_dt = (_dt + _dt_old) / 2;
53  resid = (computeQpValue() - _u_old) / (avg_dt * _dt) - (_u_dot_old) / avg_dt;
54  resid /= -_mass_diag(dofnum);
55  break;
56  }
57  return resid;
58 }
59 
60 void
62 {
63  // Now is the point that the time integrator has the variable time orders setup
65 }
66 
69 {
70  const auto & nl = _fe_problem.getNonlinearSystemBase(_sys.number());
71  if (nl.hasVector("mass_matrix_diag_inverted"))
72  return nl.getVector("mass_matrix_diag_inverted");
73 
74  mooseError("Lumped mass matrix is missing. Make sure ExplicitMixedOrder is being used as the "
75  "time integrator.");
76 }
static InputParameters validParams()
ExplicitMixedOrder::TimeOrder _var_time_order
variable time order need for computing BC
unsigned int number() const
Real & _dt_old
virtual Real computeQpValue()=0
Compute the value of the DirichletBC at the current quadrature point.
ExplicitDirichletBCBase(const InputParameters &parameters)
TimeOrder findVariableTimeOrder(unsigned int var_num) const
Retrieve the order of the highest time derivative of a variable.
Real & _dt
Implements a form of the central difference time integrator that calculates acceleration directly fro...
virtual Real computeQpResidual() override
SystemBase & _sys
const NumericVector< Number > & initLumpedMass()
Initialize the lumped mass matrix.
MooseVariable & _var
FEProblemBase & _fe_problem
const NumericVector< Number > & _mass_diag
The diagonal of the mass matrix.
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
unsigned int number() const
MooseVariableFE< Real > * _variable
virtual void timestepSetup() override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const dof_id_type & nodalDofIndex() const override
void mooseError(Args &&... args) const
static InputParameters validParams()
const ExplicitMixedOrder * _explicit_integrator
explicit time integrator
virtual NumericVector< Number > & getVector(const std::string &name)