Skip to content

Commit

Permalink
Add time transformation (idaholab#29639)
Browse files Browse the repository at this point in the history
  • Loading branch information
dschwen committed Jan 3, 2025
1 parent cf0cd9d commit 5e1aea7
Show file tree
Hide file tree
Showing 2 changed files with 11 additions and 5 deletions.
3 changes: 3 additions & 0 deletions framework/include/userobjects/SolutionUserObject.h
Original file line number Diff line number Diff line change
Expand Up @@ -498,6 +498,9 @@ class SolutionUserObject : public GeneralUserObject
/// transformations (rotations, translation, scales) are performed in this order
MultiMooseEnum _transformation_order;

/// Transform from current simulation time to time at which to sample the solution
const Function & _time_transformation;

/// True if initial_setup has executed
bool _initialized;

Expand Down
13 changes: 8 additions & 5 deletions framework/src/userobjects/SolutionUserObject.C
Original file line number Diff line number Diff line change
Expand Up @@ -86,6 +86,8 @@ SolutionUserObject::validParams()
0.0,
"Anticlockwise rotation angle (in degrees) to use for rotation about rotation1_vector.");

params.addParam<FunctionName>("time_transformation", "t", "Transfrom between simulation time and solution time.");

// following lines build the default_transformation_order
MultiMooseEnum default_transformation_order(
"rotation0 translation scale rotation1 scale_multiplier", "translation scale");
Expand Down Expand Up @@ -134,6 +136,7 @@ SolutionUserObject::SolutionUserObject(const InputParameters & parameters)
_rotation1_angle(getParam<Real>("rotation1_angle")),
_r1(RealTensorValue()),
_transformation_order(getParam<MultiMooseEnum>("transformation_order")),
_time_transformation(getFunction("time_transformation")),
_initialized(false)
{
// form rotation matrices with the specified angles
Expand Down Expand Up @@ -452,7 +455,7 @@ SolutionUserObject::timestepSetup()
{
// Update time interpolation for ExodusII solution
if (_file_type == 1 && _interpolate_times)
updateExodusTimeInterpolation(_t);
updateExodusTimeInterpolation(_time_transformation.value(_t, {}));
}

void
Expand Down Expand Up @@ -778,7 +781,7 @@ SolutionUserObject::pointValue(Real libmesh_dbg_var(t),
// Interpolate
if (_file_type == 1 && _interpolate_times)
{
mooseAssert(t == _interpolation_time,
mooseAssert(_time_transformation.value(t, {}) == _interpolation_time,
"Time passed into value() must match time at last call to timestepSetup()");
Real val2 = evalMeshFunction(pt, local_var_index, 2, subdomain_ids);
val = val + (val2 - val) * _interpolation_factor;
Expand Down Expand Up @@ -828,7 +831,7 @@ SolutionUserObject::discontinuousPointValue(Real libmesh_dbg_var(t),
// Interpolate
if (_file_type == 1 && _interpolate_times)
{
mooseAssert(t == _interpolation_time,
mooseAssert(_time_transformation.value(t, {}) == _interpolation_time,
"Time passed into value() must match time at last call to timestepSetup()");
std::map<const Elem *, Real> map2 = evalMultiValuedMeshFunction(pt, local_var_index, 2);

Expand Down Expand Up @@ -945,7 +948,7 @@ SolutionUserObject::pointValueGradient(Real libmesh_dbg_var(t),
// Interpolate
if (_file_type == 1 && _interpolate_times)
{
mooseAssert(t == _interpolation_time,
mooseAssert(_time_transformation.value(t, {}) == _interpolation_time,
"Time passed into value() must match time at last call to timestepSetup()");
RealGradient val2 = evalMeshFunctionGradient(pt, local_var_index, 2, subdomain_ids);
val = val + (val2 - val) * _interpolation_factor;
Expand Down Expand Up @@ -997,7 +1000,7 @@ SolutionUserObject::discontinuousPointValueGradient(
// Interpolate
if (_file_type == 1 && _interpolate_times)
{
mooseAssert(t == _interpolation_time,
mooseAssert(_time_transformation.value(t, {}) == _interpolation_time,
"Time passed into value() must match time at last call to timestepSetup()");
std::map<const Elem *, RealGradient> map2 =
evalMultiValuedMeshFunctionGradient(pt, local_var_index, 1, subdomain_ids);
Expand Down

0 comments on commit 5e1aea7

Please sign in to comment.