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