Line data Source code
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 : #pragma once
11 :
12 : // MOOSE includes
13 : #include "GeneralUserObject.h"
14 :
15 : // Forward declarations
16 : namespace libMesh
17 : {
18 : class ExodusII_IO;
19 : class Nemesis_IO;
20 : class EquationSystems;
21 : class System;
22 : class MeshFunction;
23 : template <class T>
24 : class NumericVector;
25 : }
26 :
27 : /**
28 : * User object that reads an existing solution from an input file and
29 : * uses it in the current simulation.
30 : */
31 : class SolutionUserObjectBase : public GeneralUserObject
32 : {
33 : public:
34 : static InputParameters validParams();
35 :
36 : SolutionUserObjectBase(const InputParameters & parameters);
37 :
38 : /**
39 : * Get the time at which to sample the solution
40 : */
41 : virtual Real solutionSampleTime() = 0;
42 :
43 : /**
44 : * When reading ExodusII files, this will update the interpolation times
45 : */
46 : virtual void timestepSetup() override;
47 :
48 : /**
49 : * Returns the local index for a given variable name
50 : * @param var_name The name of the variable for which the index is located
51 : * @return The local index of the variable
52 : */
53 : unsigned int getLocalVarIndex(const std::string & var_name) const;
54 :
55 : /**
56 : * Returns a value at a specific location and variable checking for multiple values and weighting
57 : * these values to
58 : * obtain a single unique value (see SolutionFunction)
59 : * @param t The time at which to extract (not used, it is handled automatically when reading the
60 : * data)
61 : * @param p The location at which to return a value
62 : * @param var_name The variable to be evaluated
63 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
64 : * @return The desired value for the given variable at a location
65 : */
66 : Real pointValueWrapper(Real t,
67 : const Point & p,
68 : const std::string & var_name,
69 : const MooseEnum & weighting_type = weightingType(),
70 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
71 :
72 : /**
73 : * Returns a value at a specific location and variable (see SolutionFunction)
74 : * @param t The time at which to extract (not used, it is handled automatically when reading the
75 : * data)
76 : * @param p The location at which to return a value
77 : * @param local_var_index The local index of the variable to be evaluated
78 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
79 : * @return The desired value for the given variable at a location
80 : */
81 : Real pointValue(Real t,
82 : const Point & p,
83 : const unsigned int local_var_index,
84 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
85 :
86 : /**
87 : * Returns a value at a specific location and variable (see SolutionFunction)
88 : * @param t The time at which to extract (not used, it is handled automatically when reading the
89 : * data)
90 : * @param p The location at which to return a value
91 : * @param var_name The variable to be evaluated
92 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
93 : * @return The desired value for the given variable at a location
94 : */
95 : Real pointValue(Real t,
96 : const Point & p,
97 : const std::string & var_name,
98 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
99 :
100 : /**
101 : * Returns a value at a specific location and variable for cases where the solution is
102 : * multivalued at element faces
103 : * Use pointValue for continuous shape functions or if you are sure your point is within an
104 : * element!
105 : * @param t The time at which to extract (not used, it is handled automatically when reading the
106 : * data)
107 : * @param p The location at which to return a value
108 : * @param local_var_index The local index of the variable to be evaluated
109 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
110 : * @return The desired value for the given variable at a location
111 : */
112 : std::map<const Elem *, Real>
113 : discontinuousPointValue(Real t,
114 : Point pt,
115 : const unsigned int local_var_index,
116 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
117 :
118 : /**
119 : * Returns a value at a specific location and variable for cases where the solution is
120 : * multivalued at element faces
121 : * Use pointValue for continuous shape functions or if you are sure your point is within an
122 : * element!
123 : * @param t The time at which to extract (not used, it is handled automatically when reading the
124 : * data)
125 : * @param p The location at which to return a value
126 : * @param var_name The variable to be evaluated
127 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
128 : * @return The desired value for the given variable at a location
129 : */
130 : std::map<const Elem *, Real>
131 : discontinuousPointValue(Real t,
132 : const Point & p,
133 : const std::string & var_name,
134 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
135 :
136 : /**
137 : * Returns the gradient at a specific location and variable checking for multiple values and
138 : * weighting these values to
139 : * obtain a single unique value (see SolutionFunction)
140 : * @param t The time at which to extract (not used, it is handled automatically when reading the
141 : * data)
142 : * @param p The location at which to return a value
143 : * @param var_name The variable to be evaluated
144 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
145 : * @return The desired value for the given variable at a location
146 : */
147 : libMesh::RealGradient
148 : pointValueGradientWrapper(Real t,
149 : const Point & p,
150 : const std::string & var_name,
151 : const MooseEnum & weighting_type = weightingType(),
152 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
153 :
154 : /**
155 : * Returns the gradient at a specific location and variable (see SolutionFunction)
156 : * @param t The time at which to extract (not used, it is handled automatically when reading the
157 : * data)
158 : * @param p The location at which to return a value
159 : * @param var_name The variable to be evaluated
160 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
161 : * @return The desired value for the given variable at a location
162 : */
163 : libMesh::RealGradient
164 : pointValueGradient(Real t,
165 : const Point & p,
166 : const std::string & var_name,
167 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
168 :
169 : /**
170 : * Returns the gradient at a specific location and variable (see SolutionFunction)
171 : * @param t The time at which to extract (not used, it is handled automatically when reading the
172 : * data)
173 : * @param p The location at which to return a value
174 : * @param local_var_index The local index of the variable to be evaluated
175 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
176 : * @return The desired value for the given variable at a location
177 : */
178 : libMesh::RealGradient
179 : pointValueGradient(Real t,
180 : Point pt,
181 : const unsigned int local_var_index,
182 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
183 :
184 : /**
185 : * Returns the gradient at a specific location and variable for cases where the gradient is
186 : * multivalued (e.g. at element faces)
187 : * Use pointValueGradient for continuous gradients or if you are sure your point is within an
188 : * element!
189 : * @param t The time at which to extract (not used, it is handled automatically when reading the
190 : * data)
191 : * @param p The location at which to return a value
192 : * @param var_name The variable to be evaluated
193 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
194 : * @return The desired value for the given variable at a location
195 : */
196 : std::map<const Elem *, libMesh::RealGradient> discontinuousPointValueGradient(
197 : Real t,
198 : const Point & p,
199 : const std::string & var_name,
200 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
201 :
202 : /**
203 : * Returns the gradient at a specific location and variable for cases where the gradient is
204 : * multivalued (e.g. at element faces)
205 : * Use pointValueGradient for continuous gradients or if you are sure your point is within an
206 : * element!
207 : * @param t The time at which to extract (not used, it is handled automatically when reading the
208 : * data)
209 : * @param p The location at which to return a value
210 : * @param local_var_index The local index of the variable to be evaluated
211 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
212 : * @return The desired value for the given variable at a location
213 : */
214 : std::map<const Elem *, libMesh::RealGradient> discontinuousPointValueGradient(
215 : Real t,
216 : Point pt,
217 : const unsigned int local_var_index,
218 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
219 :
220 : /**
221 : * Return a value directly from a Node
222 : * @param node A pointer to the node at which a value is desired
223 : * @param var_name The variable from which to extract a value
224 : * @return The desired value for the given node and variable name
225 : */
226 : Real directValue(const Node * node, const std::string & var_name) const;
227 :
228 : /**
229 : * Return a value from the centroid of an element
230 : * @param elem A pointer to the element at which a value is desired
231 : * @param var_name The variable from which to extract a value
232 : * @return The desired value for the given element and variable name
233 : */
234 : Real directValue(const Elem * elem, const std::string & var_name) const;
235 :
236 : /**
237 : * Returns a value of a global variable
238 : * @param t The time at which to extract (not used, it is handled automatically when reading the
239 : * data)
240 : * @param var_name The variable from which to extract a value
241 : * @return The desired value for the given variable
242 : */
243 : Real scalarValue(Real t, const std::string & var_name) const;
244 :
245 : // Required pure virtual function (not used)
246 : virtual void initialize() override;
247 :
248 : // Required pure virtual function (not used)
249 : virtual void finalize() override;
250 :
251 : // Required pure virtual function (not used)
252 : virtual void execute() override;
253 :
254 : /// Initialize the System and Mesh objects for the solution being read
255 : virtual void initialSetup() override;
256 :
257 : const std::vector<std::string> & variableNames() const;
258 :
259 : bool isVariableNodal(const std::string & var_name) const;
260 :
261 : static MooseEnum weightingType()
262 : {
263 : return MooseEnum("found_first=1 average=2 smallest_element_id=4 largest_element_id=8",
264 : "found_first");
265 : }
266 :
267 : /**
268 : * Return the spatial dimension of the mesh file
269 : * @return The spatial dimension of the mesh file
270 : */
271 : unsigned int getMeshFileDimension() const { return _mesh->spatial_dimension(); }
272 :
273 : /**
274 : * Return the name of the mesh file this object read the solution from
275 : */
276 3 : const std::string getMeshFileName() const { return _mesh_file; }
277 :
278 : /**
279 : * Get the map from block name to block ID. Only works for ExodusII files.
280 : *
281 : * @return Map from block name to block ID
282 : */
283 99 : const std::map<SubdomainName, SubdomainID> & getBlockNamesToIds() const
284 : {
285 99 : return _block_name_to_id;
286 : }
287 :
288 : /**
289 : * Get the map from block id to block name. Only works for ExodusII files.
290 : *
291 : * @return Map from block ID to block name
292 : */
293 63 : const std::map<SubdomainID, SubdomainName> & getBlockIdsToNames() const
294 : {
295 63 : return _block_id_to_name;
296 : }
297 :
298 : /**
299 : * Get the type of file that was read
300 : * @return Returns a MooseEnum that specifies the type of file read
301 : */
302 : MooseEnum getSolutionFileType() const;
303 :
304 : protected:
305 : /**
306 : * Method for reading XDA mesh and equation systems file(s)
307 : * This method is called by the constructor when 'file_type = xda' is set
308 : * in the input file.
309 : */
310 : void readXda();
311 :
312 : /**
313 : * Method for reading an ExodusII or Nemesis file, which is called when
314 : * a mesh file with a .e, .exo, .n or .nem extension is provided in the input file.
315 : */
316 : void readExodusIIOrNemesis();
317 :
318 : /**
319 : * Method for extracting value of solution based on the DOF,
320 : * this is called by the public overloaded function that accept
321 : * a node or element pointer.
322 : * @param dof_index The DOF of the solution desired
323 : * @return The solution at the given DOF
324 : */
325 : virtual Real directValue(dof_id_type dof_index) const;
326 :
327 : /**
328 : * Updates the times for interpolating ExodusII or Nemesis data
329 : */
330 : void updateTimeInterpolationFromFile();
331 :
332 : /**
333 : * Updates the time indices to interpolate between for ExodusII or Nemesis data
334 : */
335 : bool updateInterpolationBracketingTimeIndices();
336 :
337 : /**
338 : * A wrapper method for calling the various MeshFunctions used for reading the data
339 : * @param p The location at which data is desired
340 : * @param local_var_index The local index of the variable to extract data from
341 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
342 : * @param func_num The MeshFunction index to use (1 = _mesh_function; 2 = _mesh_function2)
343 : */
344 : Real evalMeshFunction(const Point & p,
345 : const unsigned int local_var_index,
346 : unsigned int func_num,
347 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
348 :
349 : /**
350 : * A wrapper method for calling the various MeshFunctions that calls the mesh function
351 : * functionality for evaluating discontinuous shape functions near a face (where it's multivalued)
352 : * @param p The location at which data is desired
353 : * @param local_var_index The local index of the variable to extract data from
354 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
355 : * @param func_num The MeshFunction index to use (1 = _mesh_function; 2 = _mesh_function2)
356 : */
357 : std::map<const Elem *, Real>
358 : evalMultiValuedMeshFunction(const Point & p,
359 : const unsigned int local_var_index,
360 : unsigned int func_num,
361 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
362 :
363 : /**
364 : * A wrapper method interfacing with the libMesh mesh function for evaluating the gradient
365 : * @param p The location at which data is desired
366 : * @param local_var_index The local index of the variable to extract data from
367 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
368 : * @param func_num The MeshFunction index to use (1 = _mesh_function; 2 = _mesh_function2)
369 : */
370 : libMesh::RealGradient
371 : evalMeshFunctionGradient(const Point & p,
372 : const unsigned int local_var_index,
373 : unsigned int func_num,
374 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
375 :
376 : /**
377 : * A wrapper method interfacing with the libMesh mesh function that calls the gradient
378 : * functionality for evaluating potentially discontinuous gradients at element's faces (where it's
379 : * multivalued)
380 : * @param p The location at which data is desired
381 : * @param local_var_index The local index of the variable to extract data from
382 : * @param subdomain_ids Subdomains IDs where to look for the value, if nullptr look everywhere
383 : * @param func_num The MeshFunction index to use (1 = _mesh_function; 2 = _mesh_function2)
384 : */
385 : std::map<const Elem *, libMesh::RealGradient> evalMultiValuedMeshFunctionGradient(
386 : const Point & p,
387 : const unsigned int local_var_index,
388 : unsigned int func_num,
389 : const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
390 :
391 : /**
392 : * Read block ID map from the ExodusII file
393 : */
394 : void readBlockIdMapFromExodusII();
395 :
396 : /// File type to read (0 = xda; 1 = ExodusII, 2 = xdr, 3 = Nemesis)
397 : MooseEnum _file_type;
398 :
399 : /// The XDA/ExodusII/XDR/Nemesis file that is being read
400 : std::string _mesh_file;
401 :
402 : /// The XDA/XDR file that contians the EquationSystems data (xda/xdr only)
403 : std::string _es_file;
404 :
405 : /// The system name to extract from the XDA/XDR file (xda/xdr only)
406 : std::string _system_name;
407 :
408 : /// A list of variables to extract from the read system
409 : std::vector<std::string> _system_variables;
410 :
411 : /// Stores the local index need by MeshFunction
412 : std::map<std::string, unsigned int> _local_variable_index;
413 :
414 : /// Stores names of nodal variables
415 : std::vector<std::string> _nodal_variables;
416 :
417 : /// Stores names of elemental variables
418 : std::vector<std::string> _elemental_variables;
419 :
420 : /// Stores names of scalar variables
421 : std::vector<std::string> _scalar_variables;
422 :
423 : /// Current ExodusII time index
424 : int _exodus_time_index;
425 :
426 : /// Flag for triggering interpolation of ExodusII data
427 : bool _interpolate_times;
428 :
429 : /// Pointer the libMesh::mesh object
430 : std::unique_ptr<libMesh::MeshBase> _mesh;
431 :
432 : /// Pointer to the libMesh::EquationSystems object
433 : std::unique_ptr<libMesh::EquationSystems> _es;
434 :
435 : /// Pointer libMesh::System class storing the read solution
436 : libMesh::System * _system;
437 :
438 : /// Pointer the libMesh::MeshFunction object that the read data is stored
439 : std::unique_ptr<libMesh::MeshFunction> _mesh_function;
440 :
441 : /// Pointer to the libMesh::ExodusII used to read the files
442 : std::unique_ptr<libMesh::ExodusII_IO> _exodusII_io;
443 :
444 : /// Pointer to the libMesh::Nemesis_IO used to read the files
445 : std::unique_ptr<libMesh::Nemesis_IO> _nemesis_io;
446 :
447 : /// Pointer to the serial solution vector
448 : std::unique_ptr<NumericVector<Number>> _serialized_solution;
449 :
450 : /// Pointer to second libMesh::EquationSystems object, used for interpolation
451 : std::unique_ptr<libMesh::EquationSystems> _es2;
452 :
453 : /// Pointer to a second libMesh::System object, used for interpolation
454 : libMesh::System * _system2;
455 :
456 : /// Pointer to second libMesh::MeshFuntion, used for interpolation
457 : std::unique_ptr<libMesh::MeshFunction> _mesh_function2;
458 :
459 : /// Pointer to second serial solution, used for interpolation
460 : std::unique_ptr<NumericVector<Number>> _serialized_solution2;
461 :
462 : /// Time in the current simulation at which the solution interpolation was last updated
463 : Real _interpolation_time;
464 :
465 : /// Interpolation weight factor
466 : Real _interpolation_factor;
467 :
468 : /// The times available in the ExodusII file
469 : const std::vector<Real> * _exodus_times;
470 :
471 : /// Time index 1, used for interpolation
472 : int _exodus_index1;
473 :
474 : /// Time index 2, used for interpolation
475 : int _exodus_index2;
476 :
477 : /// Nodal variable order, used when reading in solution data
478 : const MooseEnum _nodal_variable_order;
479 :
480 : /// Scale parameter
481 : std::vector<Real> _scale;
482 :
483 : /// scale_multiplier parameter
484 : std::vector<Real> _scale_multiplier;
485 :
486 : /// Translation
487 : std::vector<Real> _translation;
488 :
489 : /// vector about which to rotate
490 : RealVectorValue _rotation0_vector;
491 :
492 : /// angle (in degrees) which to rotate through about vector _rotation0_vector
493 : Real _rotation0_angle;
494 :
495 : /// Rotation matrix that performs the "_rotation0_angle about rotation0_vector"
496 : RealTensorValue _r0;
497 :
498 : /// vector about which to rotate
499 : RealVectorValue _rotation1_vector;
500 :
501 : /// angle (in degrees) which to rotate through about vector _rotation1_vector
502 : Real _rotation1_angle;
503 :
504 : /// Rotation matrix that performs the "_rotation1_angle about rotation1_vector"
505 : RealTensorValue _r1;
506 :
507 : /// transformations (rotations, translation, scales) are performed in this order
508 : MultiMooseEnum _transformation_order;
509 :
510 : /// Whether to force the source mesh to be replicated
511 : const bool _force_replicated_source;
512 :
513 : /// True if initial_setup has executed
514 : bool _initialized;
515 :
516 : /// Map from block ID to block names. Read from the ExodusII file
517 : std::map<SubdomainName, SubdomainID> _block_name_to_id;
518 :
519 : /// Map from block names to block IDs. Read from the ExodusII file
520 : std::map<SubdomainID, SubdomainName> _block_id_to_name;
521 :
522 : // All the caches are duplicated between the before-timestep and after-timestep mesh functions
523 : /// Cached points
524 : mutable Point _cached_p = Point(std::numeric_limits<Real>::max(), 0., 0.);
525 : mutable Point _cached_p2 = Point(std::numeric_limits<Real>::max(), 0., 0.);
526 : /// Cached subdomain ids
527 : mutable std::set<subdomain_id_type> _cached_subdomain_ids;
528 : mutable std::set<subdomain_id_type> _cached_subdomain_ids2;
529 : /// Cached values
530 : mutable DenseVector<Number> _cached_values;
531 : mutable DenseVector<Number> _cached_values2;
532 :
533 : private:
534 : static Threads::spin_mutex _solution_user_object_mutex;
535 : };
|