https://mooseframework.inl.gov
SolutionUserObjectBase.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 "SolutionUserObjectBase.h"
11 
12 // MOOSE includes
13 #include "ConsoleUtils.h"
14 #include "MooseError.h"
15 #include "MooseMesh.h"
16 #include "MooseUtils.h"
17 #include "MooseVariableFE.h"
18 #include "RotationMatrix.h"
19 #include "Function.h"
20 
21 // libMesh includes
22 #include "libmesh/equation_systems.h"
23 #include "libmesh/mesh_function.h"
24 #include "libmesh/numeric_vector.h"
25 #include "libmesh/nonlinear_implicit_system.h"
26 #include "libmesh/transient_system.h"
27 #include "libmesh/parallel_mesh.h"
28 #include "libmesh/serial_mesh.h"
29 #include "libmesh/exodusII_io.h"
30 #include "libmesh/exodusII_io_helper.h"
31 #include "libmesh/enum_xdr_mode.h"
32 #include "libmesh/string_to_enum.h"
33 
36 {
37  // Get the input parameters from the parent class
39 
40  // Add required parameters
41  params.addRequiredParam<MeshFileName>(
42  "mesh", "The name of the mesh file (must be xda/xdr or exodusII file).");
43  params.addParam<std::vector<std::string>>(
44  "system_variables",
45  std::vector<std::string>(),
46  "The name of the nodal and elemental variables from the file you want to use for values");
47 
48  // When using XDA/XDR files the following must be defined
49  params.addParam<FileName>(
50  "es",
51  "<not supplied>",
52  "The name of the file holding the equation system info in xda/xdr format (xda/xdr only).");
53  params.addParam<std::string>(
54  "system",
55  "nl0",
56  "The name of the system to pull values out of (xda/xdr only). The default name for the "
57  "nonlinear system is 'nl0', auxiliary system is 'aux0'");
58 
59  // When using ExodusII a specific time is extracted
60  params.addParam<std::string>("timestep",
61  "Index of the single timestep used or \"LATEST\" for "
62  "the last timestep (exodusII only). If not supplied, "
63  "time interpolation will occur.");
64 
65  // Choose the interpolation order for incoming nodal variables
66  params.addParam<MooseEnum>("nodal_variable_order",
67  MooseEnum("FIRST SECOND", "FIRST"),
68  "Specifies the order of the nodal solution data.");
69 
70  // Add ability to perform coordinate transformation: scale, factor
71  params.addParam<std::vector<Real>>(
72  "scale", std::vector<Real>(LIBMESH_DIM, 1), "Scale factor for points in the simulation");
73  params.addParam<std::vector<Real>>("scale_multiplier",
74  std::vector<Real>(LIBMESH_DIM, 1),
75  "Scale multiplying factor for points in the simulation");
76  params.addParam<std::vector<Real>>("translation",
77  std::vector<Real>(LIBMESH_DIM, 0),
78  "Translation factors for x,y,z coordinates of the simulation");
79  params.addParam<RealVectorValue>("rotation0_vector",
80  RealVectorValue(0, 0, 1),
81  "Vector about which to rotate points of the simulation.");
82  params.addParam<Real>(
83  "rotation0_angle",
84  0.0,
85  "Anticlockwise rotation angle (in degrees) to use for rotation about rotation0_vector.");
86  params.addParam<RealVectorValue>("rotation1_vector",
87  RealVectorValue(0, 0, 1),
88  "Vector about which to rotate points of the simulation.");
89  params.addParam<Real>(
90  "rotation1_angle",
91  0.0,
92  "Anticlockwise rotation angle (in degrees) to use for rotation about rotation1_vector.");
93 
94  // following lines build the default_transformation_order
95  MultiMooseEnum default_transformation_order(
96  "rotation0 translation scale rotation1 scale_multiplier", "translation scale");
97  params.addParam<MultiMooseEnum>(
98  "transformation_order",
99  default_transformation_order,
100  "The order to perform the operations in. Define R0 to be the rotation matrix encoded by "
101  "rotation0_vector and rotation0_angle. Similarly for R1. Denote the scale by s, the "
102  "scale_multiplier by m, and the translation by t. Then, given a point x in the simulation, "
103  "if transformation_order = 'rotation0 scale_multiplier translation scale rotation1' then "
104  "form p = R1*(R0*x*m - t)/s. Then the values provided by the SolutionUserObjectBase at "
105  "point x "
106  "in the simulation are the variable values at point p in the mesh.");
107  params.addParamNamesToGroup("scale scale_multiplier translation rotation0_vector rotation0_angle "
108  "rotation1_angle transformation_order",
109  "Coordinate system transformation");
110  params.addClassDescription("Reads a variable from a mesh in one simulation to another");
111  return params;
112 }
113 
114 // Static mutex definition
116 
118  : GeneralUserObject(parameters),
119  _file_type(MooseEnum("xda=0 exodusII=1 xdr=2")),
120  _mesh_file(getParam<MeshFileName>("mesh")),
121  _es_file(getParam<FileName>("es")),
122  _system_name(getParam<std::string>("system")),
123  _system_variables(getParam<std::vector<std::string>>("system_variables")),
124  _exodus_time_index(-1),
125  _interpolate_times(false),
126  _system(nullptr),
127  _system2(nullptr),
128  _interpolation_time(0.0),
129  _interpolation_factor(0.0),
130  _exodus_times(nullptr),
131  _exodus_index1(0),
132  _exodus_index2(0),
133  _nodal_variable_order(getParam<MooseEnum>("nodal_variable_order")),
134  _scale(getParam<std::vector<Real>>("scale")),
135  _scale_multiplier(getParam<std::vector<Real>>("scale_multiplier")),
136  _translation(getParam<std::vector<Real>>("translation")),
137  _rotation0_vector(getParam<RealVectorValue>("rotation0_vector")),
138  _rotation0_angle(getParam<Real>("rotation0_angle")),
139  _r0(RealTensorValue()),
140  _rotation1_vector(getParam<RealVectorValue>("rotation1_vector")),
141  _rotation1_angle(getParam<Real>("rotation1_angle")),
142  _r1(RealTensorValue()),
143  _transformation_order(getParam<MultiMooseEnum>("transformation_order")),
144  _initialized(false)
145 {
146  // form rotation matrices with the specified angles
147  Real halfPi = libMesh::pi / 2.0;
148  Real a;
149  Real b;
150 
151  a = std::cos(halfPi * -_rotation0_angle / 90.0);
152  b = std::sin(halfPi * -_rotation0_angle / 90.0);
153  // the following is an anticlockwise rotation about z
154  RealTensorValue rot0_z(a, -b, 0, b, a, 0, 0, 0, 1);
155  // form the rotation matrix that will take rotation0_vector to the z axis
157  // _r0 is then: rotate points so vec0 lies along z; then rotate about angle0; then rotate points
158  // back
159  _r0 = vec0_to_z.transpose() * (rot0_z * vec0_to_z);
160 
161  a = std::cos(halfPi * -_rotation1_angle / 90.0);
162  b = std::sin(halfPi * -_rotation1_angle / 90.0);
163  // the following is an anticlockwise rotation about z
164  RealTensorValue rot1_z(a, -b, 0, b, a, 0, 0, 0, 1);
165  // form the rotation matrix that will take rotation1_vector to the z axis
167  // _r1 is then: rotate points so vec1 lies along z; then rotate about angle1; then rotate points
168  // back
169  _r1 = vec1_to_z.transpose() * (rot1_z * vec1_to_z);
170 
171  if (isParamValid("timestep") && getParam<std::string>("timestep") == "-1")
172  mooseError("A \"timestep\" of -1 is no longer supported for interpolation. Instead simply "
173  "remove this parameter altogether for interpolation");
174 
175  // check for possible inconsistencies between mesh and input data
177  mooseInfo(
178  "The mesh has second-order elements, be sure to set 'nodal_variable_order' if needed.");
179 }
180 
181 void
183 {
184  if (!isParamSetByUser("es"))
185  paramError("es", "Equation system file (.xda or .xdr) should have been specified");
186 
187  // Check that the required files exist
190 
191  // Read the libmesh::mesh from the xda file
192  _mesh->read(_mesh_file);
193 
194  // Create the libmesh::EquationSystems
195  _es = std::make_unique<EquationSystems>(*_mesh);
196 
197  // Use new read syntax (binary)
198  if (_file_type == "xdr")
199  _es->read(_es_file,
201  EquationSystems::READ_HEADER | EquationSystems::READ_DATA |
202  EquationSystems::READ_ADDITIONAL_DATA);
203 
204  // Use new read syntax
205  else if (_file_type == "xda")
206  _es->read(_es_file,
208  EquationSystems::READ_HEADER | EquationSystems::READ_DATA |
209  EquationSystems::READ_ADDITIONAL_DATA);
210 
211  // This should never occur, just in case produce an error
212  else
213  mooseError("Failed to determine proper read method for XDA/XDR equation system file: ",
214  _es_file);
215 
216  // Update and store the EquationSystems name locally
217  _es->update();
218  _system = &_es->get_system(_system_name);
219 }
220 
221 void
223 {
224  // Define a default system name
225  if (_system_name == "")
226  _system_name = "SolutionUserObjectSystem";
227 
228  // Read the Exodus file
229  _exodusII_io = std::make_unique<libMesh::ExodusII_IO>(*_mesh);
230  _exodusII_io->read(_mesh_file);
232  _exodus_times = &_exodusII_io->get_time_steps();
233 
234  if (isParamValid("timestep"))
235  {
236  std::string s_timestep = getParam<std::string>("timestep");
237  int n_steps = _exodusII_io->get_num_time_steps();
238  if (s_timestep == "LATEST")
239  _exodus_time_index = n_steps;
240  else
241  {
242  std::istringstream ss(s_timestep);
243  if (!((ss >> _exodus_time_index) && ss.eof()) || _exodus_time_index > n_steps)
244  mooseError("Invalid value passed as \"timestep\". Expected \"LATEST\" or a valid integer "
245  "less than ",
246  n_steps,
247  ", received ",
248  s_timestep);
249  }
250  }
251  else
252  // Interpolate between times rather than using values from a set timestep
253  _interpolate_times = true;
254 
255  // Check that the number of time steps is valid
256  int num_exo_times = _exodus_times->size();
257  if (num_exo_times == 0)
258  mooseError("In SolutionUserObjectBase, exodus file contains no timesteps.");
259 
260  // Account for parallel mesh
261  if (dynamic_cast<DistributedMesh *>(_mesh.get()))
262  {
263  _mesh->allow_renumbering(true);
264  _mesh->prepare_for_use();
265  }
266  else
267  {
268  _mesh->allow_renumbering(false);
269  _mesh->prepare_for_use();
270  }
271 
272  // Create EquationSystems object for solution
273  _es = std::make_unique<EquationSystems>(*_mesh);
275  _system = &_es->get_system(_system_name);
276 
277  // Get the variable name lists as set; these need to be sets to perform set_intersection
278  const std::vector<std::string> & all_nodal(_exodusII_io->get_nodal_var_names());
279  const std::vector<std::string> & all_elemental(_exodusII_io->get_elem_var_names());
280  const std::vector<std::string> & all_scalar(_exodusII_io->get_global_var_names());
281 
282  // Build nodal/elemental variable lists, limit to variables listed in 'system_variables', if
283  // provided
284  // This function could be called more than once, so clear the member variables so we don't keep
285  // adding to the vectors
286  _nodal_variables.clear();
287  _elemental_variables.clear();
288  _scalar_variables.clear();
289  if (!_system_variables.empty())
290  {
291  for (const auto & var_name : _system_variables)
292  {
293  if (std::find(all_nodal.begin(), all_nodal.end(), var_name) != all_nodal.end())
294  _nodal_variables.push_back(var_name);
295  else if (std::find(all_elemental.begin(), all_elemental.end(), var_name) !=
296  all_elemental.end())
297  _elemental_variables.push_back(var_name);
298  else if (std::find(all_scalar.begin(), all_scalar.end(), var_name) != all_scalar.end())
299  // We checked other variables types first, so if a postprocessor has the same name as a
300  // field variable, it won't get loaded as a scalar variable
301  _scalar_variables.push_back(var_name);
302  else
303  paramError("system_variables", "Variable '" + var_name + "' was not found in Exodus file");
304  }
305  }
306  else
307  {
308  _nodal_variables = all_nodal;
309  _elemental_variables = all_elemental;
310 
311  for (auto var_name : all_scalar)
312  // Check if the scalar matches any field variables, and ignore the var if it does. This means
313  // its a Postprocessor.
314  if (std::find(begin(_nodal_variables), end(_nodal_variables), var_name) ==
315  _nodal_variables.end() &&
316  std::find(begin(_elemental_variables), end(_elemental_variables), var_name) ==
317  _elemental_variables.end())
318  _scalar_variables.push_back(var_name);
319  }
320 
321  // Add the variables to the system
322  for (const auto & var_name : _nodal_variables)
323  _system->add_variable(var_name, Utility::string_to_enum<Order>(_nodal_variable_order));
324 
325  for (const auto & var_name : _elemental_variables)
326  _system->add_variable(var_name, CONSTANT, MONOMIAL);
327 
328  for (const auto & var_name : _scalar_variables)
329  _system->add_variable(var_name, FIRST, SCALAR);
330 
331  // Initialize the equations systems
332  _es->init();
333 
334  // Interpolate times
335  if (_interpolate_times)
336  {
337  // Create a second equation system
338  _es2 = std::make_unique<EquationSystems>(*_mesh);
340  _system2 = &_es2->get_system(_system_name);
341 
342  // Add the variables to the system
343  for (const auto & var_name : _nodal_variables)
344  _system2->add_variable(var_name, Utility::string_to_enum<Order>(_nodal_variable_order));
345 
346  for (const auto & var_name : _elemental_variables)
347  _system2->add_variable(var_name, CONSTANT, MONOMIAL);
348 
349  for (const auto & var_name : _scalar_variables)
350  _system2->add_variable(var_name, FIRST, SCALAR);
351 
352  // Initialize
353  _es2->init();
354 
355  // Update the times for interpolation (initially start at 0)
357 
358  // Copy the solutions from the first system
359  for (const auto & var_name : _nodal_variables)
360  {
361  _exodusII_io->copy_nodal_solution(*_system, var_name, var_name, _exodus_index1 + 1);
362  _exodusII_io->copy_nodal_solution(*_system2, var_name, var_name, _exodus_index2 + 1);
363  }
364 
365  for (const auto & var_name : _elemental_variables)
366  {
367  _exodusII_io->copy_elemental_solution(*_system, var_name, var_name, _exodus_index1 + 1);
368  _exodusII_io->copy_elemental_solution(*_system2, var_name, var_name, _exodus_index2 + 1);
369  }
370 
371  if (_scalar_variables.size() > 0)
372  {
373  _exodusII_io->copy_scalar_solution(
375  _exodusII_io->copy_scalar_solution(
377  }
378 
379  // Update the systems
380  _system->update();
381  _es->update();
382  _system2->update();
383  _es2->update();
384  }
385 
386  // Non-interpolated times
387  else
388  {
389  if (_exodus_time_index > num_exo_times)
390  mooseError("In SolutionUserObjectBase, timestep = ",
392  ", but there are only ",
393  num_exo_times,
394  " time steps.");
395 
396  // Copy the values from the ExodusII file
397  for (const auto & var_name : _nodal_variables)
398  _exodusII_io->copy_nodal_solution(*_system, var_name, var_name, _exodus_time_index);
399 
400  for (const auto & var_name : _elemental_variables)
401  _exodusII_io->copy_elemental_solution(*_system, var_name, var_name, _exodus_time_index);
402 
403  if (!_scalar_variables.empty())
404  _exodusII_io->copy_scalar_solution(
406 
407  // Update the equations systems
408  _system->update();
409  _es->update();
410  }
411 }
412 
413 Real
414 SolutionUserObjectBase::directValue(const Node * node, const std::string & var_name) const
415 {
416  // Get the libmesh variable and system numbers
417  unsigned int var_num = _system->variable_number(var_name);
418  unsigned int sys_num = _system->number();
419 
420  // Get the node id and associated dof
421  dof_id_type node_id = node->id();
422  const Node & sys_node = _system->get_mesh().node_ref(node_id);
423  mooseAssert(sys_node.n_dofs(sys_num, var_num) > 0,
424  "Variable " << var_name << " has no DoFs on node " << sys_node.id());
425  dof_id_type dof_id = sys_node.dof_number(sys_num, var_num, 0);
426 
427  // Return the desired value for the dof
428  return directValue(dof_id);
429 }
430 
431 Real
432 SolutionUserObjectBase::directValue(const Elem * elem, const std::string & var_name) const
433 {
434  // Get the libmesh variable and system numbers
435  unsigned int var_num = _system->variable_number(var_name);
436  unsigned int sys_num = _system->number();
437 
438  // Get the element id and associated dof
439  dof_id_type elem_id = elem->id();
440  const Elem & sys_elem = _system->get_mesh().elem_ref(elem_id);
441  mooseAssert(sys_elem.n_dofs(sys_num, var_num) > 0,
442  "Variable " << var_name << " has no DoFs on element " << sys_elem.id());
443  dof_id_type dof_id = sys_elem.dof_number(sys_num, var_num, 0);
444 
445  // Return the desired value
446  return directValue(dof_id);
447 }
448 
449 void
451 {
452 }
453 
454 void
456 {
457 }
458 
459 void
461 {
462  // Update time interpolation for ExodusII solution
463  if (_file_type == 1 && _interpolate_times)
465 
466  // Clear the caches
469 }
470 
471 void
473 {
474 }
475 
476 void
478 {
479 
480  // Make sure this only happens once
481  if (_initialized)
482  return;
483 
484  // Create a libmesh::Mesh object for storing the loaded data.
485  // Several aspects of SolutionUserObjectBase won't work with a DistributedMesh:
486  // .) ExodusII_IO::copy_nodal_solution() doesn't work in parallel.
487  // .) We don't know if directValue will be used, which may request
488  // a value on a Node we don't have.
489  // We force the Mesh used here to be a ReplicatedMesh.
490  _mesh = std::make_unique<ReplicatedMesh>(_communicator);
491 
492  // ExodusII mesh file supplied
493  if (MooseUtils::hasExtension(_mesh_file, "e", /*strip_exodus_ext =*/true) ||
494  MooseUtils::hasExtension(_mesh_file, "exo", true))
495  {
496  _file_type = "exodusII";
497  readExodusII();
498  }
499 
500  // XDA mesh file supplied
501  else if (MooseUtils::hasExtension(_mesh_file, "xda"))
502  {
503  _file_type = "xda";
504  readXda();
505  }
506 
507  // XDR mesh file supplied
508  else if (MooseUtils::hasExtension(_mesh_file, "xdr"))
509  {
510  _file_type = "xdr";
511  readXda();
512  }
513 
514  // Produce an error for an unknown file type
515  else
516  mooseError(
517  "In SolutionUserObjectBase, invalid file type: only .xda, .xdr, .exo and .e supported");
518 
519  // Initialize the serial solution vector
522 
523  // Pull down a full copy of this vector on every processor so we can get values in parallel
525 
526  // Vector of variable numbers to apply the MeshFunction to
527  std::vector<unsigned int> var_nums;
528 
529  // If no variables were given, use all of them
530  if (_system_variables.empty())
531  {
533  for (const auto & var_num : var_nums)
534  _system_variables.push_back(_system->variable_name(var_num));
535  }
536 
537  // Otherwise, gather the numbers for the variables given
538  else
539  {
540  for (const auto & var_name : _system_variables)
541  var_nums.push_back(_system->variable_number(var_name));
542  }
543 
544  // Create the MeshFunction for working with the solution data
545  _mesh_function = std::make_unique<libMesh::MeshFunction>(
546  *_es, *_serialized_solution, _system->get_dof_map(), var_nums);
547  _mesh_function->init();
548 
549  // Tell the MeshFunctions that we might be querying them outside the
550  // mesh, so we can handle any errors at the MOOSE rather than at the
551  // libMesh level.
552  DenseVector<Number> default_values;
553  _mesh_function->enable_out_of_mesh_mode(default_values);
554 
555  // Build second MeshFunction for interpolation
556  if (_interpolate_times)
557  {
558  // Need to pull down a full copy of this vector on every processor so we can get values in
559  // parallel
563 
564  // Create the MeshFunction for the second copy of the data
565  _mesh_function2 = std::make_unique<libMesh::MeshFunction>(
567  _mesh_function2->init();
568  _mesh_function2->enable_out_of_mesh_mode(default_values);
569  }
570 
571  // Populate the MeshFunction variable index
572  for (unsigned int i = 0; i < _system_variables.size(); ++i)
573  {
574  std::string name = _system_variables[i];
576  }
577 
578  // If the start time is not the same as in the exodus file, we may need this on INITIAL
579  if (_file_type == 1 && _interpolate_times)
581 
582  // Set initialization flag
583  _initialized = true;
584 }
585 
586 MooseEnum
588 {
589  return _file_type;
590 }
591 
592 void
594 {
595  if (_t != _interpolation_time)
596  {
598  {
599 
600  for (const auto & var_name : _nodal_variables)
601  _exodusII_io->copy_nodal_solution(*_system, var_name, var_name, _exodus_index1 + 1);
602  for (const auto & var_name : _elemental_variables)
603  _exodusII_io->copy_elemental_solution(*_system, var_name, var_name, _exodus_index1 + 1);
604  if (_scalar_variables.size() > 0)
605  _exodusII_io->copy_scalar_solution(
607 
608  _system->update();
609  _es->update();
611 
612  for (const auto & var_name : _nodal_variables)
613  _exodusII_io->copy_nodal_solution(*_system2, var_name, var_name, _exodus_index2 + 1);
614  for (const auto & var_name : _elemental_variables)
615  _exodusII_io->copy_elemental_solution(*_system2, var_name, var_name, _exodus_index2 + 1);
616  if (_scalar_variables.size() > 0)
617  _exodusII_io->copy_scalar_solution(
619 
620  _system2->update();
621  _es2->update();
623  }
625  }
626 }
627 
628 bool
630 {
631  if (_file_type != 1)
632  mooseError("In SolutionUserObjectBase, getTimeInterpolationData only applicable for exodusII "
633  "file type");
634 
635  int old_index1 = _exodus_index1;
636  int old_index2 = _exodus_index2;
637 
638  int num_exo_times = _exodus_times->size();
639 
640  const auto solution_time = solutionSampleTime();
641 
642  if (solution_time < (*_exodus_times)[0])
643  {
644  _exodus_index1 = 0;
645  _exodus_index2 = 0;
646  _interpolation_factor = 0.0;
647  }
648  else
649  {
650  for (int i = 0; i < num_exo_times - 1; ++i)
651  {
652  if (solution_time <= (*_exodus_times)[i + 1])
653  {
654  _exodus_index1 = i;
655  _exodus_index2 = i + 1;
657  (solution_time - (*_exodus_times)[i]) / ((*_exodus_times)[i + 1] - (*_exodus_times)[i]);
658  break;
659  }
660  else if (i == num_exo_times - 2)
661  {
662  _exodus_index1 = num_exo_times - 1;
663  _exodus_index2 = num_exo_times - 1;
664  _interpolation_factor = 1.0;
665  break;
666  }
667  }
668  }
669 
670  bool indices_modified(false);
671 
672  if (_exodus_index1 != old_index1 || _exodus_index2 != old_index2)
673  indices_modified = true;
674 
675  return indices_modified;
676 }
677 
678 unsigned int
679 SolutionUserObjectBase::getLocalVarIndex(const std::string & var_name) const
680 {
681  // Extract the variable index for the MeshFunction(s)
682  std::map<std::string, unsigned int>::const_iterator it = _local_variable_index.find(var_name);
683  if (it == _local_variable_index.end())
684  mooseError("Value requested for nonexistent variable '",
685  var_name,
686  "' in the '",
687  name(),
688  "' SolutionUserObjectBase.\nSystem selected: ",
689  _system_name,
690  "\nAvailable variables:\n",
692  return it->second;
693 }
694 
695 Real
697  const Point & p,
698  const std::string & var_name,
699  const MooseEnum & weighting_type,
700  const std::set<subdomain_id_type> * subdomain_ids) const
701 {
702  // first check if the FE type is continuous because in that case the value is
703  // unique and we can take a short cut, the default weighting_type found_first also
704  // shortcuts out
705  auto family =
707  .getVariable(
709  .feType()
710  .family;
711 
712  if (weighting_type == 1 ||
713  (family != L2_LAGRANGE && family != MONOMIAL && family != L2_HIERARCHIC))
714  return pointValue(t, p, var_name, subdomain_ids);
715 
716  // the shape function is discontinuous so we need to compute a suitable unique value
717  std::map<const Elem *, Real> values = discontinuousPointValue(t, p, var_name);
718  switch (weighting_type)
719  {
720  case 2:
721  {
722  Real average = 0.0;
723  for (auto & v : values)
724  average += v.second;
725  return average / Real(values.size());
726  }
727  case 4:
728  {
729  Real smallest_elem_id_value = std::numeric_limits<Real>::max();
730  dof_id_type smallest_elem_id = _fe_problem.mesh().maxElemId();
731  for (auto & v : values)
732  if (v.first->id() < smallest_elem_id)
733  {
734  smallest_elem_id = v.first->id();
735  smallest_elem_id_value = v.second;
736  }
737  return smallest_elem_id_value;
738  }
739  case 8:
740  {
741  Real largest_elem_id_value = std::numeric_limits<Real>::lowest();
742  dof_id_type largest_elem_id = 0;
743  for (auto & v : values)
744  if (v.first->id() > largest_elem_id)
745  {
746  largest_elem_id = v.first->id();
747  largest_elem_id_value = v.second;
748  }
749  return largest_elem_id_value;
750  }
751  }
752 
753  mooseError("SolutionUserObjectBase::pointValueWrapper reaches line that it should not be able to "
754  "reach.");
755  return 0.0;
756 }
757 
758 Real
760  const Point & p,
761  const std::string & var_name,
762  const std::set<subdomain_id_type> * subdomain_ids) const
763 {
764  const unsigned int local_var_index = getLocalVarIndex(var_name);
765  return pointValue(t, p, local_var_index, subdomain_ids);
766 }
767 
768 Real
769 SolutionUserObjectBase::pointValue(Real libmesh_dbg_var(t),
770  const Point & p,
771  const unsigned int local_var_index,
772  const std::set<subdomain_id_type> * subdomain_ids) const
773 {
774  // Create copy of point
775  Point pt(p);
776 
777  // do the transformations
778  for (unsigned int trans_num = 0; trans_num < _transformation_order.size(); ++trans_num)
779  {
780  if (_transformation_order[trans_num] == "rotation0")
781  pt = _r0 * pt;
782  else if (_transformation_order[trans_num] == "translation")
783  for (const auto i : make_range(Moose::dim))
784  pt(i) -= _translation[i];
785  else if (_transformation_order[trans_num] == "scale")
786  for (const auto i : make_range(Moose::dim))
787  pt(i) /= _scale[i];
788  else if (_transformation_order[trans_num] == "scale_multiplier")
789  for (const auto i : make_range(Moose::dim))
790  pt(i) *= _scale_multiplier[i];
791  else if (_transformation_order[trans_num] == "rotation1")
792  pt = _r1 * pt;
793  }
794 
795  // Extract the value at the current point
796  Real val = evalMeshFunction(pt, local_var_index, 1, subdomain_ids);
797 
798  // Interpolate
799  if (_file_type == 1 && _interpolate_times)
800  {
801  mooseAssert(t == _interpolation_time,
802  "Time passed into value() must match time at last call to timestepSetup()");
803  Real val2 = evalMeshFunction(pt, local_var_index, 2, subdomain_ids);
804  val = val + (val2 - val) * _interpolation_factor;
805  }
806 
807  return val;
808 }
809 
810 std::map<const Elem *, Real>
812  Real t,
813  const Point & p,
814  const std::string & var_name,
815  const std::set<subdomain_id_type> * subdomain_ids) const
816 {
817  const unsigned int local_var_index = getLocalVarIndex(var_name);
818  return discontinuousPointValue(t, p, local_var_index, subdomain_ids);
819 }
820 
821 std::map<const Elem *, Real>
823  Real libmesh_dbg_var(t),
824  Point pt,
825  const unsigned int local_var_index,
826  const std::set<subdomain_id_type> * subdomain_ids) const
827 {
828  // do the transformations
829  for (unsigned int trans_num = 0; trans_num < _transformation_order.size(); ++trans_num)
830  {
831  if (_transformation_order[trans_num] == "rotation0")
832  pt = _r0 * pt;
833  else if (_transformation_order[trans_num] == "translation")
834  for (const auto i : make_range(Moose::dim))
835  pt(i) -= _translation[i];
836  else if (_transformation_order[trans_num] == "scale")
837  for (const auto i : make_range(Moose::dim))
838  pt(i) /= _scale[i];
839  else if (_transformation_order[trans_num] == "scale_multiplier")
840  for (const auto i : make_range(Moose::dim))
841  pt(i) *= _scale_multiplier[i];
842  else if (_transformation_order[trans_num] == "rotation1")
843  pt = _r1 * pt;
844  }
845 
846  // Extract the value at the current point
847  std::map<const Elem *, Real> map =
848  evalMultiValuedMeshFunction(pt, local_var_index, 1, subdomain_ids);
849 
850  // Interpolate
851  if (_file_type == 1 && _interpolate_times)
852  {
853  mooseAssert(t == _interpolation_time,
854  "Time passed into value() must match time at last call to timestepSetup()");
855  std::map<const Elem *, Real> map2 = evalMultiValuedMeshFunction(pt, local_var_index, 2);
856 
857  if (map.size() != map2.size())
858  mooseError(
859  "In SolutionUserObjectBase::discontinuousPointValue map and map2 have different size");
860 
861  // construct the interpolated map
862  for (auto & k : map)
863  {
864  if (map2.find(k.first) == map2.end())
865  mooseError(
866  "In SolutionUserObjectBase::discontinuousPointValue map and map2 have differing keys");
867  Real val = k.second;
868  Real val2 = map2[k.first];
869  map[k.first] = val + (val2 - val) * _interpolation_factor;
870  }
871  }
872 
873  return map;
874 }
875 
878  Real t,
879  const Point & p,
880  const std::string & var_name,
881  const MooseEnum & weighting_type,
882  const std::set<subdomain_id_type> * subdomain_ids) const
883 {
884  // the default weighting_type found_first shortcuts out
885  if (weighting_type == 1)
886  return pointValueGradient(t, p, var_name, subdomain_ids);
887 
888  // the shape function is discontinuous so we need to compute a suitable unique value
889  std::map<const Elem *, RealGradient> values =
890  discontinuousPointValueGradient(t, p, var_name, subdomain_ids);
891  switch (weighting_type)
892  {
893  case 2:
894  {
895  RealGradient average = RealGradient(0.0, 0.0, 0.0);
896  for (auto & v : values)
897  average += v.second;
898  return average / Real(values.size());
899  }
900  case 4:
901  {
902  RealGradient smallest_elem_id_value;
903  dof_id_type smallest_elem_id = _fe_problem.mesh().maxElemId();
904  for (auto & v : values)
905  if (v.first->id() < smallest_elem_id)
906  {
907  smallest_elem_id = v.first->id();
908  smallest_elem_id_value = v.second;
909  }
910  return smallest_elem_id_value;
911  }
912  case 8:
913  {
914  RealGradient largest_elem_id_value;
915  dof_id_type largest_elem_id = 0;
916  for (auto & v : values)
917  if (v.first->id() > largest_elem_id)
918  {
919  largest_elem_id = v.first->id();
920  largest_elem_id_value = v.second;
921  }
922  return largest_elem_id_value;
923  }
924  }
925 
926  mooseError("SolutionUserObjectBase::pointValueGradientWrapper reaches line that it should not be "
927  "able to reach.");
928  return RealGradient(0.0, 0.0, 0.0);
929 }
930 
933  const Point & p,
934  const std::string & var_name,
935  const std::set<subdomain_id_type> * subdomain_ids) const
936 {
937  const unsigned int local_var_index = getLocalVarIndex(var_name);
938  return pointValueGradient(t, p, local_var_index, subdomain_ids);
939 }
940 
942 SolutionUserObjectBase::pointValueGradient(Real libmesh_dbg_var(t),
943  Point pt,
944  const unsigned int local_var_index,
945  const std::set<subdomain_id_type> * subdomain_ids) const
946 {
947  // do the transformations
948  for (unsigned int trans_num = 0; trans_num < _transformation_order.size(); ++trans_num)
949  {
950  if (_transformation_order[trans_num] == "rotation0")
951  pt = _r0 * pt;
952  else if (_transformation_order[trans_num] == "translation")
953  for (const auto i : make_range(Moose::dim))
954  pt(i) -= _translation[i];
955  else if (_transformation_order[trans_num] == "scale")
956  for (const auto i : make_range(Moose::dim))
957  pt(i) /= _scale[i];
958  else if (_transformation_order[trans_num] == "scale_multiplier")
959  for (const auto i : make_range(Moose::dim))
960  pt(i) *= _scale_multiplier[i];
961  else if (_transformation_order[trans_num] == "rotation1")
962  pt = _r1 * pt;
963  }
964 
965  // Extract the value at the current point
966  RealGradient val = evalMeshFunctionGradient(pt, local_var_index, 1, subdomain_ids);
967 
968  // Interpolate
969  if (_file_type == 1 && _interpolate_times)
970  {
971  mooseAssert(t == _interpolation_time,
972  "Time passed into value() must match time at last call to timestepSetup()");
973  RealGradient val2 = evalMeshFunctionGradient(pt, local_var_index, 2, subdomain_ids);
974  val = val + (val2 - val) * _interpolation_factor;
975  }
976 
977  return val;
978 }
979 
980 std::map<const Elem *, RealGradient>
982  Real t,
983  const Point & p,
984  const std::string & var_name,
985  const std::set<subdomain_id_type> * subdomain_ids) const
986 {
987  const unsigned int local_var_index = getLocalVarIndex(var_name);
988  return discontinuousPointValueGradient(t, p, local_var_index, subdomain_ids);
989 }
990 
991 std::map<const Elem *, RealGradient>
993  Real libmesh_dbg_var(t),
994  Point pt,
995  const unsigned int local_var_index,
996  const std::set<subdomain_id_type> * subdomain_ids) const
997 {
998  // do the transformations
999  for (unsigned int trans_num = 0; trans_num < _transformation_order.size(); ++trans_num)
1000  {
1001  if (_transformation_order[trans_num] == "rotation0")
1002  pt = _r0 * pt;
1003  else if (_transformation_order[trans_num] == "translation")
1004  for (const auto i : make_range(Moose::dim))
1005  pt(i) -= _translation[i];
1006  else if (_transformation_order[trans_num] == "scale")
1007  for (const auto i : make_range(Moose::dim))
1008  pt(i) /= _scale[i];
1009  else if (_transformation_order[trans_num] == "scale_multiplier")
1010  for (const auto i : make_range(Moose::dim))
1011  pt(i) *= _scale_multiplier[i];
1012  else if (_transformation_order[trans_num] == "rotation1")
1013  pt = _r1 * pt;
1014  }
1015 
1016  // Extract the value at the current point
1017  std::map<const Elem *, RealGradient> map =
1018  evalMultiValuedMeshFunctionGradient(pt, local_var_index, 1, subdomain_ids);
1019 
1020  // Interpolate
1021  if (_file_type == 1 && _interpolate_times)
1022  {
1023  mooseAssert(t == _interpolation_time,
1024  "Time passed into value() must match time at last call to timestepSetup()");
1025  std::map<const Elem *, RealGradient> map2 =
1026  evalMultiValuedMeshFunctionGradient(pt, local_var_index, 1, subdomain_ids);
1027 
1028  if (map.size() != map2.size())
1029  mooseError(
1030  "In SolutionUserObjectBase::discontinuousPointValue map and map2 have different size");
1031 
1032  // construct the interpolated map
1033  for (auto & k : map)
1034  {
1035  if (map2.find(k.first) == map2.end())
1036  mooseError(
1037  "In SolutionUserObjectBase::discontinuousPointValue map and map2 have differing keys");
1038  RealGradient val = k.second;
1039  RealGradient val2 = map2[k.first];
1040  map[k.first] = val + (val2 - val) * _interpolation_factor;
1041  }
1042  }
1043 
1044  return map;
1045 }
1046 
1047 Real
1049 {
1050  Real val = (*_serialized_solution)(dof_index);
1051  if (_file_type == 1 && _interpolate_times)
1052  {
1053  Real val2 = (*_serialized_solution2)(dof_index);
1054  val = val + (val2 - val) * _interpolation_factor;
1055  }
1056  return val;
1057 }
1058 
1059 Real
1061  const unsigned int local_var_index,
1062  unsigned int func_num,
1063  const std::set<subdomain_id_type> * subdomain_ids) const
1064 {
1065  // Storage for mesh function output
1066  DenseVector<Number> output;
1067 
1068  // Extract a value from the _mesh_function
1069  {
1070  Threads::spin_mutex::scoped_lock lock(_solution_user_object_mutex);
1071 
1072  if (func_num == 1)
1073  {
1074  // Check the cache
1075  if (p == _cached_p && (!subdomain_ids || (*subdomain_ids == _cached_subdomain_ids)))
1076  return _cached_values(local_var_index);
1077 
1078  // else get a new value
1079  (*_mesh_function)(p, 0.0, output, subdomain_ids);
1080  // and cache it
1081  _cached_p = p;
1082  if (subdomain_ids)
1083  _cached_subdomain_ids = *subdomain_ids;
1084  _cached_values = output;
1085  }
1086 
1087  // Extract a value from _mesh_function2
1088  else if (func_num == 2)
1089  {
1090  // Check the cache
1091  if (p == _cached_p2 && (!subdomain_ids || (*subdomain_ids == _cached_subdomain_ids2)))
1092  return _cached_values2(local_var_index);
1093 
1094  // else get a new value
1095  (*_mesh_function2)(p, 0.0, output, subdomain_ids);
1096  // and cache it
1097  _cached_p2 = p;
1098  if (subdomain_ids)
1099  _cached_subdomain_ids2 = *subdomain_ids;
1100  _cached_values2 = output;
1101  }
1102  else
1103  mooseError("The func_num must be 1 or 2");
1104  }
1105 
1106  // Error if the data is out-of-range, which will be the case if the mesh functions are evaluated
1107  // outside the domain
1108  if (output.size() == 0)
1109  {
1110  std::ostringstream oss;
1111  p.print(oss);
1112  mooseError("Failed to access the data for variable '",
1113  _system_variables[local_var_index],
1114  "' at point ",
1115  oss.str(),
1116  " in the '",
1117  name(),
1118  "' SolutionUserObjectBase");
1119  }
1120  return output(local_var_index);
1121 }
1122 
1123 std::map<const Elem *, Real>
1125  const Point & p,
1126  const unsigned int local_var_index,
1127  unsigned int func_num,
1128  const std::set<subdomain_id_type> * subdomain_ids) const
1129 {
1130  // Storage for mesh function output
1131  std::map<const Elem *, DenseVector<Number>> temporary_output;
1132 
1133  // Extract a value from the _mesh_function
1134  {
1135  Threads::spin_mutex::scoped_lock lock(_solution_user_object_mutex);
1136  if (func_num == 1)
1137  _mesh_function->discontinuous_value(p, 0.0, temporary_output, subdomain_ids);
1138 
1139  // Extract a value from _mesh_function2
1140  else if (func_num == 2)
1141  _mesh_function2->discontinuous_value(p, 0.0, temporary_output, subdomain_ids);
1142 
1143  else
1144  mooseError("The func_num must be 1 or 2");
1145  }
1146 
1147  // Error if the data is out-of-range, which will be the case if the mesh functions are evaluated
1148  // outside the domain
1149  if (temporary_output.size() == 0)
1150  {
1151  std::ostringstream oss;
1152  p.print(oss);
1153  mooseError("Failed to access the data for variable '",
1154  _system_variables[local_var_index],
1155  "' at point ",
1156  oss.str(),
1157  " in the '",
1158  name(),
1159  "' SolutionUserObjectBase");
1160  }
1161 
1162  // Fill the actual map that is returned
1163  std::map<const Elem *, Real> output;
1164  for (auto & k : temporary_output)
1165  {
1166  mooseAssert(
1167  k.second.size() > local_var_index,
1168  "In SolutionUserObjectBase::evalMultiValuedMeshFunction variable with local_var_index "
1169  << local_var_index << " does not exist");
1170  output[k.first] = k.second(local_var_index);
1171  }
1172 
1173  return output;
1174 }
1175 
1178  const Point & p,
1179  const unsigned int local_var_index,
1180  unsigned int func_num,
1181  const std::set<subdomain_id_type> * subdomain_ids) const
1182 {
1183  // Storage for mesh function output
1184  std::vector<Gradient> output;
1185 
1186  // Extract a value from the _mesh_function
1187  {
1188  Threads::spin_mutex::scoped_lock lock(_solution_user_object_mutex);
1189  if (func_num == 1)
1190  _mesh_function->gradient(p, 0.0, output, subdomain_ids);
1191 
1192  // Extract a value from _mesh_function2
1193  else if (func_num == 2)
1194  _mesh_function2->gradient(p, 0.0, output, subdomain_ids);
1195 
1196  else
1197  mooseError("The func_num must be 1 or 2");
1198  }
1199 
1200  // Error if the data is out-of-range, which will be the case if the mesh functions are evaluated
1201  // outside the domain
1202  if (output.size() == 0)
1203  {
1204  std::ostringstream oss;
1205  p.print(oss);
1206  mooseError("Failed to access the data for variable '",
1207  _system_variables[local_var_index],
1208  "' at point ",
1209  oss.str(),
1210  " in the '",
1211  name(),
1212  "' SolutionUserObjectBase");
1213  }
1214  return output[local_var_index];
1215 }
1216 
1217 std::map<const Elem *, RealGradient>
1219  const Point & p,
1220  const unsigned int local_var_index,
1221  unsigned int func_num,
1222  const std::set<subdomain_id_type> * subdomain_ids) const
1223 {
1224  // Storage for mesh function output
1225  std::map<const Elem *, std::vector<Gradient>> temporary_output;
1226 
1227  // Extract a value from the _mesh_function
1228  {
1229  Threads::spin_mutex::scoped_lock lock(_solution_user_object_mutex);
1230  if (func_num == 1)
1231  _mesh_function->discontinuous_gradient(p, 0.0, temporary_output, subdomain_ids);
1232 
1233  // Extract a value from _mesh_function2
1234  else if (func_num == 2)
1235  _mesh_function2->discontinuous_gradient(p, 0.0, temporary_output, subdomain_ids);
1236 
1237  else
1238  mooseError("The func_num must be 1 or 2");
1239  }
1240 
1241  // Error if the data is out-of-range, which will be the case if the mesh functions are evaluated
1242  // outside the domain
1243  if (temporary_output.size() == 0)
1244  {
1245  std::ostringstream oss;
1246  p.print(oss);
1247  mooseError("Failed to access the data for variable '",
1248  _system_variables[local_var_index],
1249  "' at point ",
1250  oss.str(),
1251  " in the '",
1252  name(),
1253  "' SolutionUserObjectBase");
1254  }
1255 
1256  // Fill the actual map that is returned
1257  std::map<const Elem *, RealGradient> output;
1258  for (auto & k : temporary_output)
1259  {
1260  mooseAssert(
1261  k.second.size() > local_var_index,
1262  "In SolutionUserObjectBase::evalMultiValuedMeshFunction variable with local_var_index "
1263  << local_var_index << " does not exist");
1264  output[k.first] = k.second[local_var_index];
1265  }
1266 
1267  return output;
1268 }
1269 
1270 const std::vector<std::string> &
1272 {
1273  return _system_variables;
1274 }
1275 
1276 bool
1277 SolutionUserObjectBase::isVariableNodal(const std::string & var_name) const
1278 {
1279  return std::find(_nodal_variables.begin(), _nodal_variables.end(), var_name) !=
1280  _nodal_variables.end();
1281 }
1282 
1283 Real
1284 SolutionUserObjectBase::scalarValue(Real /*t*/, const std::string & var_name) const
1285 {
1286  unsigned int var_num = _system->variable_number(var_name);
1287  const DofMap & dof_map = _system->get_dof_map();
1288  std::vector<dof_id_type> dofs;
1289  dof_map.SCALAR_dof_indices(dofs, var_num);
1290  // We can handle only FIRST order scalar variables
1291  return directValue(dofs[0]);
1292 }
1293 
1294 void
1296 {
1297 #ifdef LIBMESH_HAVE_EXODUS_API
1298  libMesh::ExodusII_IO_Helper & exio_helper = _exodusII_io->get_exio_helper();
1299  const auto & id_to_block = exio_helper.id_to_block_names;
1300  _block_name_to_id.clear();
1301  _block_id_to_name.clear();
1302  for (const auto & it : id_to_block)
1303  {
1304  _block_name_to_id[it.second] = it.first;
1305  _block_id_to_name[it.first] = it.second;
1306  }
1307 #endif
1308 }
L2_HIERARCHIC
MultiMooseEnum _transformation_order
transformations (rotations, translation, scales) are performed in this order
void mooseInfo(Args &&... args) const
Definition: MooseBase.h:321
std::vector< std::string > _system_variables
A list of variables to extract from the read system.
static InputParameters validParams()
SolutionUserObjectBase(const InputParameters &parameters)
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
KOKKOS_INLINE_FUNCTION const T * find(const T &target, const T *const begin, const T *const end)
Find a value in an array.
Definition: KokkosUtils.h:42
std::unique_ptr< libMesh::MeshFunction > _mesh_function2
Pointer to second libMesh::MeshFuntion, used for interpolation.
const libMesh::FEType & feType() const
Get the type of finite element object.
std::unique_ptr< libMesh::EquationSystems > _es
Pointer to the libMesh::EquationSystems object.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
RealVectorValue _rotation0_vector
vector about which to rotate
const std::vector< Real > * _exodus_times
The times available in the ExodusII file.
virtual void initialSetup() override
Initialize the System and Mesh objects for the solution being read.
virtual Real solutionSampleTime()=0
Get the time at which to sample the solution.
SCALAR
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:439
static InputParameters validParams()
void readXda()
Method for reading XDA mesh and equation systems file(s) This method is called by the constructor whe...
static Threads::spin_mutex _solution_user_object_mutex
std::map< const Elem *, libMesh::RealGradient > discontinuousPointValueGradient(Real t, const Point &p, const std::string &var_name, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns the gradient at a specific location and variable for cases where the gradient is multivalued ...
FIRST
TypeTensor< T > transpose() const
std::unique_ptr< NumericVector< Number > > _serialized_solution2
Pointer to second serial solution, used for interpolation.
std::vector< Real > _scale
Scale parameter.
std::vector< Real > _translation
Translation.
Real _rotation1_angle
angle (in degrees) which to rotate through about vector _rotation1_vector
unsigned int size() const
Return the number of active items in the MultiMooseEnum.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual dof_id_type maxElemId() const
Definition: MooseMesh.C:3173
Real _rotation0_angle
angle (in degrees) which to rotate through about vector _rotation0_vector
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:162
MooseEnum getSolutionFileType() const
Get the type of file that was read.
void get_all_variable_numbers(std::vector< unsigned int > &all_variable_numbers) const
int _exodus_time_index
Current ExodusII time index.
std::unique_ptr< libMesh::ExodusII_IO > _exodusII_io
Pointer to the libMesh::ExodusII used to read the files.
Real _interpolation_time
Time in the current simulation at which the solution interpolation was last updated.
std::unique_ptr< libMesh::EquationSystems > _es2
Pointer to second libMesh::EquationSystems object, used for interpolation.
const Parallel::Communicator & _communicator
std::map< const Elem *, Real > evalMultiValuedMeshFunction(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method for calling the various MeshFunctions that calls the mesh function functionality for...
libMesh::RealGradient evalMeshFunctionGradient(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method interfacing with the libMesh mesh function for evaluating the gradient.
std::unique_ptr< NumericVector< Number > > _serialized_solution
Pointer to the serial solution vector.
Point _cached_p
Cached points.
bool _initialized
True if initial_setup has executed.
const MeshBase & get_mesh() const
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
virtual const MooseVariableFieldBase & getVariable(const THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) const override
Returns the variable reference for requested variable which must be of the expected_var_type (Nonline...
dof_id_type n_dofs() const
std::map< SubdomainName, SubdomainID > _block_name_to_id
Map from block ID to block names. Read from the ExodusII file.
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
auto max(const L &left, const R &right)
unsigned int variable_number(std::string_view var) const
virtual void timestepSetup() override
When reading ExodusII files, this will update the interpolation times.
DenseVector< Number > _cached_values2
virtual void execute() override
Execute method.
DenseVector< Number > _cached_values
Cached values.
std::unique_ptr< libMesh::MeshBase > _mesh
Pointer the libMesh::mesh object.
libMesh::System * _system
Pointer libMesh::System class storing the read solution.
bool updateExodusBracketingTimeIndices()
Updates the time indices to interpolate between for ExodusII data.
CONSTANT
unsigned int number() const
RealVectorValue _rotation1_vector
vector about which to rotate
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:103
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cos(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cos
std::string formatString(std::string message, const std::string &prefix)
Add new lines and prefixes to a string for pretty display in output NOTE: This makes a copy of the st...
Definition: ConsoleUtils.C:581
std::string _mesh_file
The XDA or ExodusII file that is being read.
std::unique_ptr< NumericVector< Number > > solution
std::vector< std::string > _scalar_variables
Stores names of scalar variables.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
RealTensorValue _r1
Rotation matrix that performs the "_rotation1_angle about rotation1_vector".
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
const std::string & variable_name(const unsigned int i) const
void readExodusII()
Method for reading an ExodusII file, which is called when &#39;file_type = exodusII is set in the input f...
std::string stringify(const T &t)
conversion to string
Definition: Conversion.h:64
bool hasSecondOrderElements()
check if the mesh has SECOND order elements
Definition: MooseMesh.C:3791
virtual void initialize() override
Called before execute() is ever called so that data can be cleared.
Real evalMeshFunction(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method for calling the various MeshFunctions used for reading the data. ...
MONOMIAL
void updateExodusTimeInterpolation()
Updates the times for interpolating ExodusII data.
Real scalarValue(Real t, const std::string &var_name) const
Returns a value of a global variable.
const MooseEnum _nodal_variable_order
Nodal variable order, used when reading in solution data.
libMesh::System * _system2
Pointer to a second libMesh::System object, used for interpolation.
GenericRealTensorValue< is_ad > rotVecToZ(GenericRealVectorValue< is_ad > vec)
provides a rotation matrix that will rotate the vector vec to the z axis (the "2" direction) ...
virtual void update()
bool hasExtension(const std::string &filename, std::string ext, bool strip_exodus_ext)
Definition: MooseUtils.C:387
Real directValue(const Node *node, const std::string &var_name) const
Return a value directly from a Node.
void readBlockIdMapFromExodusII()
Read block ID map from the ExodusII file.
const std::vector< std::string > & variableNames() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
FEProblemBase & _fe_problem
Reference to the FEProblemBase for this user object.
Definition: UserObject.h:211
bool isVariableNodal(const std::string &var_name) const
std::map< std::string, unsigned int > _local_variable_index
Stores the local index need by MeshFunction.
std::map< const Elem *, libMesh::RealGradient > evalMultiValuedMeshFunctionGradient(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method interfacing with the libMesh mesh function that calls the gradient functionality for...
Real pointValueWrapper(Real t, const Point &p, const std::string &var_name, const MooseEnum &weighting_type=weightingType(), const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns a value at a specific location and variable checking for multiple values and weighting these ...
int _exodus_index1
Time index 1, used for interpolation.
virtual const Elem & elem_ref(const dof_id_type i) const
std::vector< std::string > _nodal_variables
Stores names of nodal variables.
std::map< int, std::string > id_to_block_names
std::string _es_file
The XDA file that contians the EquationSystems data (xda only)
const THREAD_ID _tid
Thread ID of this postprocessor.
Definition: UserObject.h:218
IntRange< T > make_range(T beg, T end)
virtual MooseMesh & mesh() override
static std::unique_ptr< NumericVector< Number > > build(const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC)
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:271
RealTensorValue _r0
Rotation matrix that performs the "_rotation0_angle about rotation0_vector".
virtual unsigned int size() const override final
L2_LAGRANGE
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
std::vector< std::string > _elemental_variables
Stores names of elemental variables.
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
virtual const Node & node_ref(const dof_id_type i) const
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:199
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...
MooseEnum _file_type
File type to read (0 = xda; 1 = ExodusII)
libMesh::RealGradient pointValueGradientWrapper(Real t, const Point &p, const std::string &var_name, const MooseEnum &weighting_type=weightingType(), const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns the gradient at a specific location and variable checking for multiple values and weighting t...
std::vector< Real > _scale_multiplier
scale_multiplier parameter
std::unique_ptr< libMesh::MeshFunction > _mesh_function
Pointer the libMesh::MeshFunction object that the read data is stored.
bool _interpolate_times
Flag for triggering interpolation of ExodusII data.
bool checkFileReadable(const std::string &filename, bool check_line_endings, bool throw_on_unreadable, bool check_for_git_lfs_pointer)
Definition: MooseUtils.C:254
std::string _system_name
The system name to extract from the XDA file (xda only)
const DofMap & get_dof_map() const
bool isParamSetByUser(const std::string &name) const
Test if the supplied parameter is set by a user, as opposed to not set or set to default.
Definition: MooseBase.h:205
std::map< const Elem *, Real > discontinuousPointValue(Real t, Point pt, const unsigned int local_var_index, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns a value at a specific location and variable for cases where the solution is multivalued at el...
Real pointValue(Real t, const Point &p, const unsigned int local_var_index, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns a value at a specific location and variable (see SolutionFunction)
unsigned int getLocalVarIndex(const std::string &var_name) const
Returns the local index for a given variable name.
int _exodus_index2
Time index 2, used for interpolation.
Real _interpolation_factor
Interpolation weight factor.
std::map< SubdomainID, SubdomainName > _block_id_to_name
Map from block names to block IDs. Read from the ExodusII file.
virtual void finalize() override
Finalize.
uint8_t dof_id_type
const Real pi
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...
std::set< subdomain_id_type > _cached_subdomain_ids2
libMesh::RealGradient pointValueGradient(Real t, const Point &p, const std::string &var_name, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns the gradient at a specific location and variable (see SolutionFunction)
std::set< subdomain_id_type > _cached_subdomain_ids
Cached subdomain ids.