https://mooseframework.inl.gov
CFLTimeStepSize.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 // Navier-Stokes includes
11 #include "CFLTimeStepSize.h"
12 
13 // MOOSE includes
14 #include "libmesh/quadrature.h"
15 #include "metaphysicl/raw_type.h"
16 
17 #include <algorithm>
18 #include <limits>
19 
20 registerMooseObject("NavierStokesApp", CFLTimeStepSize);
21 registerMooseObject("NavierStokesApp", ADCFLTimeStepSize);
22 
23 template <bool is_ad>
26 {
28 
29  params.addParam<Real>("CFL", 0.5, "The CFL number to use in computing time step size");
30  params.addRequiredParam<std::vector<MaterialPropertyName>>("vel_names",
31  "Velocity material property name(s)");
32  params.addRequiredParam<std::vector<MaterialPropertyName>>(
33  "c_names", "Sound speed material property name(s)");
34 
35  // Because this post-processor is meant to be used with PostprocessorDT, it
36  // should be executed on initial (not included by default) and timestep end.
37  params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
38  params.suppressParameter<ExecFlagEnum>("execute_on");
39 
40  params.addClassDescription("Computes a time step size based on a user-specified CFL number");
41 
42  return params;
43 }
44 
45 template <bool is_ad>
47  : ElementPostprocessor(parameters),
48  _CFL(getParam<Real>("CFL")),
49  _vel_names(getParam<std::vector<MaterialPropertyName>>("vel_names")),
50  _c_names(getParam<std::vector<MaterialPropertyName>>("c_names")),
51  _n_phases(_vel_names.size()),
52  _dt(std::numeric_limits<Real>::max())
53 {
54  if (_vel_names.size() != _c_names.size())
55  mooseError("The number of elements in the parameters 'vel_names' and 'c_names' must be equal.");
56 
57  for (unsigned int k = 0; k < _n_phases; ++k)
58  {
59  _vel.push_back(&getGenericMaterialPropertyByName<Real, is_ad>(_vel_names[k]));
60  _c.push_back(&getGenericMaterialPropertyByName<Real, is_ad>(_c_names[k]));
61  }
62 }
63 
64 template <bool is_ad>
65 void
67 {
68  // start with the max
69  _dt = std::numeric_limits<Real>::max();
70 }
71 
72 template <bool is_ad>
73 Real
75 {
76  return _dt;
77 }
78 
79 template <bool is_ad>
80 void
82 {
83  gatherMin(_dt);
84 }
85 
86 template <bool is_ad>
87 void
89 {
90  const auto & pps = static_cast<const CFLTimeStepSizeTempl<is_ad> &>(y);
91 
92  _dt = std::min(_dt, pps._dt);
93 }
94 
95 template <bool is_ad>
96 void
98 {
99  // get minimum element diameter for element
100  const Real h_min_element = _current_elem->hmin();
101 
102  // loop over quadrature points
103  for (unsigned int qp = 0; qp < _qrule->n_points(); qp++)
104  {
105  // determine minimum time step size over all phases
106  for (unsigned int k = 0; k < _n_phases; ++k)
107  {
108  const Real dt_phase = _CFL * h_min_element /
109  (std::fabs(MetaPhysicL::raw_value((*_vel[k])[qp])) +
110  MetaPhysicL::raw_value((*_c[k])[qp]));
111  _dt = std::min(_dt, dt_phase);
112  }
113  }
114 }
virtual void initialize() override
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
T & set(const std::string &name, bool quiet_mode=false)
static InputParameters validParams()
auto raw_value(const Eigen::Map< T > &in)
const std::vector< double > y
const ExecFlagType EXEC_TIMESTEP_END
registerMooseObject("NavierStokesApp", CFLTimeStepSize)
void addRequiredParam(const std::string &name, const std::string &doc_string)
auto max(const L &left, const R &right)
void suppressParameter(const std::string &name)
const std::vector< MaterialPropertyName > & _vel_names
Velocity material property name(s)
const std::vector< MaterialPropertyName > & _c_names
Sound speed material property name(s)
static InputParameters validParams()
CFLTimeStepSizeTempl(const InputParameters &parameters)
virtual void threadJoin(const UserObject &y) override
Real _dt
Time step size.
Computes a time step size based on user-specified CFL number.
virtual Real getValue() const override
std::vector< const GenericMaterialProperty< Real, is_ad > * > _vel
Velocity material properties.
virtual void finalize() override
const unsigned int _n_phases
Number of phases.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void execute() override
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
std::vector< const GenericMaterialProperty< Real, is_ad > * > _c
Sound speed material properties.
static const std::string k
Definition: NS.h:130
const ExecFlagType EXEC_INITIAL