20 #include "libmesh/equation_systems.h"
21 #include "libmesh/mesh_function.h"
22 #include "libmesh/numeric_vector.h"
23 #include "libmesh/nonlinear_implicit_system.h"
24 #include "libmesh/transient_system.h"
25 #include "libmesh/parallel_mesh.h"
26 #include "libmesh/serial_mesh.h"
27 #include "libmesh/exodusII_io.h"
28 #include "libmesh/enum_xdr_mode.h"
42 "mesh",
"The name of the mesh file (must be xda or exodusII file).");
43 params.
addParam<std::vector<std::string>>(
45 std::vector<std::string>(),
46 "The name of the nodal and elemental variables from the file you want to use for values");
52 "The name of the file holding the equation system info in xda format (xda only).");
54 "system",
"nl0",
"The name of the system to pull values out of (xda only).");
57 params.
addParam<std::string>(
"timestep",
58 "Index of the single timestep used or \"LATEST\" for "
59 "the last timestep (exodusII only). If not supplied, "
60 "time interpolation will occur.");
64 "scale", std::vector<Real>(LIBMESH_DIM, 1),
"Scale factor for points in the simulation");
65 params.
addParam<std::vector<Real>>(
"scale_multiplier",
66 std::vector<Real>(LIBMESH_DIM, 1),
67 "Scale multiplying factor for points in the simulation");
68 params.
addParam<std::vector<Real>>(
"translation",
69 std::vector<Real>(LIBMESH_DIM, 0),
70 "Translation factors for x,y,z coordinates of the simulation");
73 "Vector about which to rotate points of the simulation.");
77 "Anticlockwise rotation angle (in degrees) to use for rotation about rotation0_vector.");
80 "Vector about which to rotate points of the simulation.");
84 "Anticlockwise rotation angle (in degrees) to use for rotation about rotation1_vector.");
88 "rotation0 translation scale rotation1 scale_multiplier",
"translation scale");
90 "transformation_order",
91 default_transformation_order,
92 "The order to perform the operations in. Define R0 to be the rotation matrix encoded by "
93 "rotation0_vector and rotation0_angle. Similarly for R1. Denote the scale by s, the "
94 "scale_multiplier by m, and the translation by t. Then, given a point x in the simulation, "
95 "if transformation_order = 'rotation0 scale_multiplier translation scale rotation1' then "
96 "form p = R1*(R0*x*m - t)/s. Then the values provided by the SolutionUserObject at point x "
97 "in the simulation are the variable values at point p in the mesh.");
108 _file_type(
MooseEnum(
"xda=0 exodusII=1 xdr=2")),
109 _mesh_file(getParam<MeshFileName>(
"mesh")),
110 _es_file(getParam<FileName>(
"es")),
111 _system_name(getParam<
std::string>(
"system")),
112 _system_variables(getParam<
std::vector<
std::string>>(
"system_variables")),
113 _exodus_time_index(-1),
114 _interpolate_times(false),
117 _interpolation_time(0.0),
118 _interpolation_factor(0.0),
119 _exodus_times(nullptr),
122 _scale(getParam<
std::vector<Real>>(
"scale")),
123 _scale_multiplier(getParam<
std::vector<Real>>(
"scale_multiplier")),
124 _translation(getParam<
std::vector<Real>>(
"translation")),
126 _rotation0_angle(getParam<Real>(
"rotation0_angle")),
129 _rotation1_angle(getParam<Real>(
"rotation1_angle")),
131 _transformation_order(getParam<
MultiMooseEnum>(
"transformation_order")),
135 Real halfPi = std::acos(0.0);
147 _r0 = vec0_to_z.transpose() * (rot0_z * vec0_to_z);
157 _r1 = vec1_to_z.transpose() * (rot1_z * vec1_to_z);
159 if (
isParamValid(
"timestep") && getParam<std::string>(
"timestep") ==
"-1")
160 mooseError(
"A \"timestep\" of -1 is no longer supported for interpolation. Instead simply "
161 "remove this parameter altogether for interpolation");
177 _es = libmesh_make_unique<EquationSystems>(*
_mesh);
183 EquationSystems::READ_HEADER | EquationSystems::READ_DATA |
184 EquationSystems::READ_ADDITIONAL_DATA);
190 EquationSystems::READ_HEADER | EquationSystems::READ_DATA |
191 EquationSystems::READ_ADDITIONAL_DATA);
195 mooseError(
"Failed to determine proper read method for XDA/XDR equation system file: ",
217 std::string s_timestep = getParam<std::string>(
"timestep");
219 if (s_timestep ==
"LATEST")
223 std::istringstream ss(s_timestep);
225 mooseError(
"Invalid value passed as \"timestep\". Expected \"LATEST\" or a valid integer "
238 if (num_exo_times == 0)
239 mooseError(
"In SolutionUserObject, exodus file contains no timesteps.");
242 if (dynamic_cast<DistributedMesh *>(
_mesh.get()))
244 _mesh->allow_renumbering(
true);
245 _mesh->prepare_for_use();
249 _mesh->allow_renumbering(
false);
250 _mesh->prepare_for_use();
254 _es = libmesh_make_unique<EquationSystems>(*
_mesh);
259 const std::vector<std::string> & all_nodal(
_exodusII_io->get_nodal_var_names());
260 const std::vector<std::string> & all_elemental(
_exodusII_io->get_elem_var_names());
261 const std::vector<std::string> & all_scalar(
_exodusII_io->get_global_var_names());
269 if (std::find(all_nodal.begin(), all_nodal.end(), var_name) != all_nodal.end())
271 if (std::find(all_elemental.begin(), all_elemental.end(), var_name) != all_elemental.end())
273 if (std::find(all_scalar.begin(), all_scalar.end(), var_name) != all_scalar.end())
286 _system->add_variable(var_name, FIRST);
289 _system->add_variable(var_name, CONSTANT, MONOMIAL);
292 _system->add_variable(var_name, FIRST, SCALAR);
301 _es2 = libmesh_make_unique<EquationSystems>(*
_mesh);
307 _system2->add_variable(var_name, FIRST);
310 _system2->add_variable(var_name, CONSTANT, MONOMIAL);
313 _system2->add_variable(var_name, FIRST, SCALAR);
353 mooseError(
"In SolutionUserObject, timestep = ",
355 ", but there are only ",
380 unsigned int var_num =
_system->variable_number(var_name);
381 unsigned int sys_num =
_system->number();
384 dof_id_type node_id = node->id();
385 const Node & sys_node =
_system->get_mesh().node_ref(node_id);
386 mooseAssert(sys_node.n_dofs(sys_num, var_num) > 0,
387 "Variable " << var_name <<
" has no DoFs on node " << sys_node.id());
388 dof_id_type dof_id = sys_node.dof_number(sys_num, var_num, 0);
398 unsigned int var_num =
_system->variable_number(var_name);
399 unsigned int sys_num =
_system->number();
402 dof_id_type elem_id = elem->id();
403 const Elem & sys_elem =
_system->get_mesh().elem_ref(elem_id);
404 mooseAssert(sys_elem.n_dofs(sys_num, var_num) > 0,
405 "Variable " << var_name <<
" has no DoFs on element " << sys_elem.id());
406 dof_id_type dof_id = sys_elem.dof_number(sys_num, var_num, 0);
449 _mesh = libmesh_make_unique<ReplicatedMesh>(_communicator);
473 mooseError(
"In SolutionUserObject, invalid file type (only .xda, .xdr, and .e supported)");
483 std::vector<unsigned int> var_nums;
488 _system->get_all_variable_numbers(var_nums);
489 for (
const auto & var_num : var_nums)
497 var_nums.push_back(
_system->variable_number(var_name));
508 DenseVector<Number> default_values;
585 "In SolutionUserObject, getTimeInterpolationData only applicable for exodusII file type");
600 for (
int i = 0; i < num_exo_times - 1; ++i)
607 (time - (*_exodus_times)[i]) / ((*
_exodus_times)[i + 1] - (*_exodus_times)[i]);
610 else if (i == num_exo_times - 2)
620 bool indices_modified(
false);
623 indices_modified =
true;
625 return indices_modified;
634 mooseError(
"Value requested for nonexistent variable '",
638 "' SolutionUserObject");
645 const std::string & var_name,
647 const std::set<subdomain_id_type> * subdomain_ids)
const
659 if (weighting_type == 1 ||
660 (family != L2_LAGRANGE && family != MONOMIAL && family != L2_HIERARCHIC))
661 return pointValue(t, p, var_name, subdomain_ids);
665 switch (weighting_type)
670 for (
auto & v : values)
672 return average / Real(values.size());
676 Real smallest_elem_id_value = std::numeric_limits<Real>::max();
678 for (
auto & v : values)
679 if (v.first->id() < smallest_elem_id)
681 smallest_elem_id = v.first->id();
682 smallest_elem_id_value = v.second;
684 return smallest_elem_id_value;
688 Real largest_elem_id_value = std::numeric_limits<Real>::lowest();
689 dof_id_type largest_elem_id = 0;
690 for (
auto & v : values)
691 if (v.first->id() > largest_elem_id)
693 largest_elem_id = v.first->id();
694 largest_elem_id_value = v.second;
696 return largest_elem_id_value;
701 "SolutionUserObject::pointValueWrapper reaches line that it should not be able to reach.");
708 const std::string & var_name,
709 const std::set<subdomain_id_type> * subdomain_ids)
const
712 return pointValue(t, p, local_var_index, subdomain_ids);
718 const unsigned int local_var_index,
719 const std::set<subdomain_id_type> * subdomain_ids)
const
730 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
733 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
736 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
749 "Time passed into value() must match time at last call to timestepSetup()");
757 std::map<const Elem *, Real>
760 const std::string & var_name,
761 const std::set<subdomain_id_type> * subdomain_ids)
const
767 std::map<const Elem *, Real>
770 const unsigned int local_var_index,
771 const std::set<subdomain_id_type> * subdomain_ids)
const
779 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
782 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
785 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
792 std::map<const Elem *, Real> map =
799 "Time passed into value() must match time at last call to timestepSetup()");
802 if (map.size() != map2.size())
803 mooseError(
"In SolutionUserObject::discontinuousPointValue map and map2 have different size");
808 if (map2.find(k.first) == map2.end())
810 "In SolutionUserObject::discontinuousPointValue map and map2 have differing keys");
812 Real val2 = map2[k.first];
824 const std::string & var_name,
826 const std::set<subdomain_id_type> * subdomain_ids)
const
829 if (weighting_type == 1)
833 std::map<const Elem *, RealGradient> values =
835 switch (weighting_type)
839 RealGradient average = RealGradient(0.0, 0.0, 0.0);
840 for (
auto & v : values)
842 return average / Real(values.size());
846 RealGradient smallest_elem_id_value;
848 for (
auto & v : values)
849 if (v.first->id() < smallest_elem_id)
851 smallest_elem_id = v.first->id();
852 smallest_elem_id_value = v.second;
854 return smallest_elem_id_value;
858 RealGradient largest_elem_id_value;
859 dof_id_type largest_elem_id = 0;
860 for (
auto & v : values)
861 if (v.first->id() > largest_elem_id)
863 largest_elem_id = v.first->id();
864 largest_elem_id_value = v.second;
866 return largest_elem_id_value;
870 mooseError(
"SolutionUserObject::pointValueGradientWrapper reaches line that it should not be "
872 return RealGradient(0.0, 0.0, 0.0);
878 const std::string & var_name,
879 const std::set<subdomain_id_type> * subdomain_ids)
const
888 const unsigned int local_var_index,
889 const std::set<subdomain_id_type> * subdomain_ids)
const
897 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
900 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
903 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
916 "Time passed into value() must match time at last call to timestepSetup()");
924 std::map<const Elem *, RealGradient>
928 const std::string & var_name,
929 const std::set<subdomain_id_type> * subdomain_ids)
const
935 std::map<const Elem *, RealGradient>
937 Real libmesh_dbg_var(t),
939 const unsigned int local_var_index,
940 const std::set<subdomain_id_type> * subdomain_ids)
const
948 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
951 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
954 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
961 std::map<const Elem *, RealGradient> map =
968 "Time passed into value() must match time at last call to timestepSetup()");
969 std::map<const Elem *, RealGradient> map2 =
972 if (map.size() != map2.size())
973 mooseError(
"In SolutionUserObject::discontinuousPointValue map and map2 have different size");
978 if (map2.find(k.first) == map2.end())
980 "In SolutionUserObject::discontinuousPointValue map and map2 have differing keys");
981 RealGradient val = k.second;
982 RealGradient val2 = map2[k.first];
993 Real val = (*_serialized_solution)(dof_index);
996 Real val2 = (*_serialized_solution2)(dof_index);
1004 const unsigned int local_var_index,
1005 unsigned int func_num,
1006 const std::set<subdomain_id_type> * subdomain_ids)
const
1009 DenseVector<Number> output;
1015 (*_mesh_function)(p, 0.0, output, subdomain_ids);
1018 else if (func_num == 2)
1019 (*_mesh_function2)(p, 0.0, output, subdomain_ids);
1027 if (output.size() == 0)
1029 std::ostringstream oss;
1031 mooseError(
"Failed to access the data for variable '",
1037 "' SolutionUserObject");
1039 return output(local_var_index);
1042 std::map<const Elem *, Real>
1045 const unsigned int local_var_index,
1046 unsigned int func_num,
1047 const std::set<subdomain_id_type> * subdomain_ids)
const
1050 std::map<const Elem *, DenseVector<Number>> temporary_output;
1056 _mesh_function->discontinuous_value(p, 0.0, temporary_output, subdomain_ids);
1059 else if (func_num == 2)
1060 _mesh_function2->discontinuous_value(p, 0.0, temporary_output, subdomain_ids);
1068 if (temporary_output.size() == 0)
1070 std::ostringstream oss;
1072 mooseError(
"Failed to access the data for variable '",
1078 "' SolutionUserObject");
1082 std::map<const Elem *, Real> output;
1083 for (
auto & k : temporary_output)
1085 mooseAssert(k.second.size() > local_var_index,
1086 "In SolutionUserObject::evalMultiValuedMeshFunction variable with local_var_index "
1087 << local_var_index <<
" does not exist");
1088 output[k.first] = k.second(local_var_index);
1097 const unsigned int local_var_index,
1098 unsigned int func_num,
1099 const std::set<subdomain_id_type> * subdomain_ids)
const
1102 std::vector<Gradient> output;
1111 else if (func_num == 2)
1120 if (output.size() == 0)
1122 std::ostringstream oss;
1124 mooseError(
"Failed to access the data for variable '",
1130 "' SolutionUserObject");
1132 return output[local_var_index];
1135 std::map<const Elem *, RealGradient>
1138 const unsigned int local_var_index,
1139 unsigned int func_num,
1140 const std::set<subdomain_id_type> * subdomain_ids)
const
1143 std::map<const Elem *, std::vector<Gradient>> temporary_output;
1149 _mesh_function->discontinuous_gradient(p, 0.0, temporary_output, subdomain_ids);
1152 else if (func_num == 2)
1153 _mesh_function2->discontinuous_gradient(p, 0.0, temporary_output, subdomain_ids);
1161 if (temporary_output.size() == 0)
1163 std::ostringstream oss;
1165 mooseError(
"Failed to access the data for variable '",
1171 "' SolutionUserObject");
1175 std::map<const Elem *, RealGradient> output;
1176 for (
auto & k : temporary_output)
1178 mooseAssert(k.second.size() > local_var_index,
1179 "In SolutionUserObject::evalMultiValuedMeshFunction variable with local_var_index "
1180 << local_var_index <<
" does not exist");
1181 output[k.first] = k.second[local_var_index];
1187 const std::vector<std::string> &
1203 unsigned int var_num =
_system->variable_number(var_name);
1205 std::vector<dof_id_type> dofs;
1206 dof_map.SCALAR_dof_indices(dofs, var_num);