Line data Source code
1 : /********************************************************************/
2 : /* SOFTWARE COPYRIGHT NOTIFICATION */
3 : /* Cardinal */
4 : /* */
5 : /* (c) 2021 UChicago Argonne, LLC */
6 : /* ALL RIGHTS RESERVED */
7 : /* */
8 : /* Prepared by UChicago Argonne, LLC */
9 : /* Under Contract No. DE-AC02-06CH11357 */
10 : /* With the U. S. Department of Energy */
11 : /* */
12 : /* Prepared by Battelle Energy Alliance, LLC */
13 : /* Under Contract No. DE-AC07-05ID14517 */
14 : /* With the U. S. Department of Energy */
15 : /* */
16 : /* See LICENSE for full restrictions */
17 : /********************************************************************/
18 :
19 : #pragma once
20 :
21 : #include "OpenMCProblemBase.h"
22 : #include "SymmetryPointGenerator.h"
23 :
24 : /// Tally/filter includes.
25 : #include "TallyBase.h"
26 : #include "FilterBase.h"
27 :
28 : #ifdef ENABLE_DAGMC
29 : #include "MoabSkinner.h"
30 : #include "DagMC.hpp"
31 : #endif
32 :
33 : /// Forward declarations to avoid cyclic dependencies.
34 : class OpenMCVolumeCalculation;
35 :
36 : /**
37 : * Mapping of OpenMC to a collection of MOOSE elements, with temperature and/or
38 : * density feedback. The mappind is established authomatically by looping over
39 : * all the MOOSE elements and finding the OpenMC cell at each element's centroid.
40 : *
41 : * All feedback into OpenMC is performed via element averages. The
42 : * 'temperature_blocks' parameter is used to indicate which MOOSE blocks should
43 : * provide temperature feedback, while the 'density_blocks' parameter is used to
44 : * indicate which MOOSE blocks should provide density feedback. Tallies are
45 : * automatically added to OpenMC using either cell or mesh tallies.
46 : *
47 : * Each OpenMC cell shall not have ambiguous data transfers. That is, a cell
48 : * should map to a set of elements that are ALL/NONE providing temperature
49 : * feedback, ALL/NONE providing density feedback, ALL/NONE providing cell
50 : * tallies, and ALL/NONE being uncoupled altogether.
51 : *
52 : * TODO: If this is too restrictive in the future, we could implement some type
53 : * of weighted averaging process. Also, if a cell maps to a phase and an
54 : * unmapped region, perhaps we want to allow that.
55 : *
56 : * Other considerations you should be aware of:
57 : * - The density being transferred into OpenMC from MOOSE is in units of kg/m3; this
58 : * is the unit employed by the MOOSE fluid properties module.
59 : * - The temperature being transferred into OpenMC from MOOSE is in units of K; this
60 : * is the unit employed by the MOOSE fluid and solid properties modules.
61 : * - If your geometry is highly fine-scale (such as TRISOs), you might be able to get a
62 : * huge speedup in your runtime if you (i) build your OpenMC model by repeating the same
63 : * TRISO universe in each of your repeatable-units (e.g. pebbles, compacts, plates)
64 : * AND (ii) leverage the 'identical_cell_fills' option.
65 : */
66 : class OpenMCCellAverageProblem : public OpenMCProblemBase
67 : {
68 : public:
69 : OpenMCCellAverageProblem(const InputParameters & params);
70 : static InputParameters validParams();
71 :
72 : virtual void initialSetup() override;
73 : virtual void externalSolve() override;
74 : virtual void syncSolutions(ExternalProblem::Direction direction) override;
75 2290 : virtual bool converged(unsigned int) override { return true; }
76 :
77 : /**
78 : * Read a 2d vector of subdomain names, and check that there are no duplications
79 : * and that all provided values exist on the mesh.
80 : * @param[in] name string name for the 2d vector parameter
81 : * @param[out] names subdomain names
82 : * @param[out] flattened_ids flattened 1d vector of subdomain IDs
83 : */
84 : void read2DBlockParameters(const std::string name,
85 : std::vector<std::vector<SubdomainName>> & names,
86 : std::vector<SubdomainID> & flattened_ids);
87 :
88 : /**
89 : * Check that the specified blocks are in the mesh
90 : * @param[in] name name for throwing an error
91 : * @param[in] ids block IDs to check
92 : * @param[in] names block subdomain names for throwing an error
93 : */
94 : void checkBlocksInMesh(const std::string name,
95 : const std::vector<SubdomainID> & ids,
96 : const std::vector<SubdomainName> & names) const;
97 :
98 : /// Initialize the mapping of OpenMC to the MooseMesh and perform additional setup actions
99 : void setupProblem();
100 :
101 : /**
102 : * Add the tally variable(s) (to receive OpenMC tally values), temperature variable(s)
103 : * (to write into OpenMC cells), and density variable(s) (to write into OpenMC materials)
104 : */
105 : virtual void addExternalVariables() override;
106 :
107 : /**
108 : * Get the cell volume from a stochastic calculation
109 : * @param[in] cell_info cell index, instance pair
110 : * @return stochastically-computed OpenMC cell volume
111 : */
112 : virtual Real cellVolume(const cellInfo & cell_info) const;
113 :
114 : /**
115 : * Reference to stochastic volume calculation
116 : * @return reference to stochastic volume calculation
117 : */
118 609314 : virtual const OpenMCVolumeCalculation * volumeCalculation() const { return _volume_calc; }
119 :
120 : /**
121 : * Get the mapping of cells to MOOSE elements
122 : * @return mapping of cells to MOOSE elements
123 : */
124 9054 : virtual const std::map<cellInfo, std::vector<unsigned int>> & cellToElem() const
125 : {
126 9054 : return _cell_to_elem;
127 : }
128 :
129 : /**
130 : * Get the MOOSE subdomains associated with an OpenMC cell
131 : * @param info the cell info
132 : * @return MOOSE subdomains associated with an OpenMC cell
133 : */
134 74804 : virtual std::unordered_set<SubdomainID> getCellToElementSub(const cellInfo & info)
135 : {
136 74804 : return _cell_to_elem_subdomain.at(info);
137 : }
138 :
139 : /**
140 : * Whether transformations are applied to the [Mesh] points when mapping to OpenMC
141 : * @return whether transformations are applied
142 : */
143 8728156 : virtual bool hasPointTransformations() const { return _symmetry != nullptr; }
144 :
145 : /**
146 : * Get all the scores added to the tally
147 : * @return scores
148 : */
149 336 : virtual const std::vector<std::string> & getTallyScores() const { return _all_tally_scores; }
150 :
151 : /**
152 : * Check to see if this problem contains a specific tally score.
153 : * @param[in] score the tally score
154 : * @return whether this problem contains the tally score in a tally object
155 : */
156 298 : bool hasScore(const std::string & score)
157 : {
158 298 : return std::find(_all_tally_scores.begin(), _all_tally_scores.end(), score) !=
159 298 : _all_tally_scores.end();
160 : }
161 :
162 : /**
163 : * Fetch the tally that contains the requested score.
164 : * @param[in] score the tally score
165 : * @return a vector of Cardinal tally objects containing the scores
166 : */
167 : std::vector<const TallyBase *> getTalliesByScore(const std::string & score);
168 :
169 : /**
170 : * Get the variable(s) associated with an OpenMC tally score.
171 : * @param[in] score the OpenMC score
172 : * @param[in] tid the thread ID associated with the current MOOSE object
173 : * @param[in] output the output variable (relative error, standard deviation, etc.) to fetch
174 : * @param[in] skip_func_exp whether functional expansion filter bins should be skipped or not when
175 : * fetching variable values
176 : * @return a vector of variable values associated with score
177 : */
178 : std::vector<const MooseVariableFE<Real> *> getTallyScoreVariables(const std::string & score,
179 : THREAD_ID tid,
180 : const std::string & output = "",
181 : bool skip_func_exp = false);
182 :
183 : /**
184 : * Get the variable value(s) associated with an OpenMC tally score.
185 : * @param[in] score the OpenMC score
186 : * @param[in] tid the thread ID associated with the current MOOSE object
187 : * @param[in] output the output variable (relative error, standard deviation, etc.) to fetch
188 : * @param[in] skip_func_exp whether functional expansion filter bins should be skipped or not when
189 : * fetching variable values
190 : * @return a vector of variable values associated with score
191 : */
192 : std::vector<const VariableValue *> getTallyScoreVariableValues(const std::string & score,
193 : THREAD_ID tid,
194 : const std::string & output = "",
195 : bool skip_func_exp = false);
196 :
197 : /**
198 : * Get the variable value(s) associated with an OpenMC tally score.
199 : * @param[in] score the OpenMC score
200 : * @param[in] tid the thread ID associated with the current MOOSE object
201 : * @param[in] output the output variable (relative error, standard deviation, etc.) to fetch
202 : * @param[in] skip_func_exp whether functional expansion filter bins should be skipped or not when
203 : * fetching variable values
204 : * @return a vector of variable values associated with score
205 : */
206 : std::vector<const VariableValue *>
207 : getTallyScoreNeighborVariableValues(const std::string & score,
208 : THREAD_ID tid,
209 : const std::string & output = "",
210 : bool skip_func_exp = false);
211 :
212 : /**
213 : * Whether a tally contains a specified output or not.
214 : * @param[in] score the tally score to check
215 : * @param[in] output the additional output (unrelaxed standard deviation, relative error, or
216 : * tally)
217 : * @return whether an added tally has the output or not
218 : */
219 : bool hasOutput(const std::string & score, const std::string & output) const;
220 :
221 : /**
222 : * Apply transformations to point
223 : * @param[in] pt point
224 : * @return transformed point
225 : */
226 8679002 : virtual Point transformPoint(const Point & pt) const
227 : {
228 8679002 : return this->hasPointTransformations() ? _symmetry->transformPoint(pt) : pt;
229 : }
230 :
231 : /**
232 : * This class uses elem->volume() in order to normalize the tally values. However,
233 : * elem->volume() is expensive, so whenever MOOSE does integration, they set
234 : * _current_elem_volume to the volume as set by the sum of the quadrature weights.
235 : * The quadrature rule that MOOSE provides when you only have CONSTANT MONOMIALS is
236 : * insufficient for exactly integrating the element Jacobian mapping type (which
237 : * is FIRST LAGRANGE for a first order element), so you get an error relative to
238 : * the libmesh volume computation.
239 : *
240 : * So, we need to make sure that a minimum order quadrature rule is used
241 : * so that the total tally as computed by an
242 : * ElementIntegralVariablePostprocessor actually matches the specified total
243 : * (for low quadrature orders, there can be an error up to about 5% or so in total
244 : * power). This override simply forces the volume quadrature order to be 2 or higher
245 : * when using Gauss (default), monomial, or Gauss-Lobatto quadrature.
246 : *
247 : * For other quadrature rules, the approximations made in elem->volume() are never
248 : * going to match the volume integrations in MOOSE (unless the quadrature order is
249 : * very very high). For these orders, we print an error message informing the user
250 : * that they should switch to a different order.
251 : */
252 : virtual void createQRules(QuadratureType type,
253 : Order order,
254 : Order volume_order,
255 : Order face_order,
256 : SubdomainID block,
257 : bool allow_negative_weights = true) override;
258 :
259 : /**
260 : * Type definition for cells contained within a parent cell; the first value
261 : * is the cell index, while the second is the set of cell instances
262 : */
263 : typedef std::unordered_map<int32_t, std::vector<int32_t>> containedCells;
264 :
265 : /**
266 : * Get the cell index from the element ID; will return UNMAPPED for unmapped elements
267 : * @param[in] elem_id element ID
268 : * @return cell index
269 : */
270 : int32_t elemToCellIndex(const int & elem_id) const { return elemToCellInfo(elem_id).first; }
271 :
272 : /**
273 : * Get the cell ID from the element ID. Note that this function requires that the elem_id
274 : * maps to an OpenMC cell, or else an error will be raised from OpenMC in cellID.
275 : * @param[in] elem_id element ID
276 : * @return cell ID
277 : */
278 70422272 : int32_t elemToCellID(const int & elem_id) const { return cellID(elemToCellIndex(elem_id)); }
279 :
280 : /**
281 : * Get the cell instance from the element ID; will return UNMAPPED for unmapped elements
282 : * @param[in] elem_id element ID
283 : * @return cell instance
284 : */
285 : int32_t elemToCellInstance(const int & elem_id) const { return elemToCellInfo(elem_id).second; }
286 :
287 : /**
288 : * Get the cell index, instance pair from element ID; if the element doesn't map to an OpenMC
289 : * cell, the index and instance are both set to UNMAPPED
290 : * @param[in] elem_id element ID
291 : * @return cell index, instance pair
292 : */
293 90375296 : cellInfo elemToCellInfo(const int & elem_id) const { return _elem_to_cell[elem_id]; }
294 :
295 : /**
296 : * Get the cell material index based on index, instance pair. Note that this function requires
297 : * a valid instance, index pair for cellInfo - you cannot pass in an unmapped cell, i.e.
298 : * (UNMAPPED, UNMAPPED)
299 : * @param[in] cell_info cell index, instance pair
300 : * @return material index
301 : */
302 : int32_t cellToMaterialIndex(const cellInfo & cell_info) const
303 : {
304 2838912 : return _cell_to_material.at(cell_info);
305 : }
306 :
307 : /**
308 : * Get the fields coupled for each cell; because we require that each cell maps to a consistent
309 : * set, we simply look up the coupled fields of the first element that this cell maps to. Note
310 : * that this function requires a valid instance, index pair for cellInfo - you cannot pass in an
311 : * unmapped cell, i.e. (UNMAPPED, UNMAPPED)
312 : * @param[in] cell_info cell index, instance pair
313 : * @return coupling fields
314 : */
315 : coupling::CouplingFields cellFeedback(const cellInfo & cell_info) const;
316 :
317 : /**
318 : * Whether a cell has density feedback
319 : * @param[in] cell_info cell index,instance pair
320 : * @return if cell has density feedback
321 : */
322 6959629 : bool hasDensityFeedback(const cellInfo & cell_info) const
323 : {
324 : std::vector<coupling::CouplingFields> phase = {coupling::density,
325 6959629 : coupling::density_and_temperature};
326 6959629 : return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
327 : }
328 :
329 : /**
330 : * Whether a cell has temperature feedback
331 : * @param[in] cell_info cell index,instance pair
332 : * @return if cell has temperature feedback
333 : */
334 8106434 : bool hasTemperatureFeedback(const cellInfo & cell_info) const
335 : {
336 : std::vector<coupling::CouplingFields> phase = {coupling::temperature,
337 8106434 : coupling::density_and_temperature};
338 8106434 : return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
339 : }
340 :
341 : /**
342 : * Checks if the [Problem/Filters] block contains a specific filter.
343 : * @param[in] filter_name the MOOSE object name of the filter
344 : * @return whether the problem contains the specified filter
345 : */
346 : bool hasFilter(const std::string & filter_name) const { return _filters.count(filter_name) > 0; }
347 :
348 : /**
349 : * Get a filter added by the [Problem/Filters] block by it's MOOSE object name.
350 : * @param[in] filter_name the MOOSE object name of the filter
351 : * @return the filter object
352 : */
353 : std::shared_ptr<FilterBase> & getFilter(const std::string & filter_name)
354 : {
355 1120 : return _filters.at(filter_name);
356 : }
357 :
358 : /**
359 : * Get the local tally
360 : * @return local tally
361 : */
362 : const std::vector<std::shared_ptr<TallyBase>> & getLocalTally() const { return _local_tallies; }
363 :
364 : /**
365 : * Get the temperature of a cell; for cells not filled with materials, this will return
366 : * the temperature of the first material-type cell
367 : * @param[in] cell_info cell index, instance pair
368 : */
369 : double cellTemperature(const cellInfo & cell_info) const;
370 :
371 : /**
372 : * Get the volume that each OpenMC cell mapped to
373 : * @param[in] cell_info cell index, instance pair
374 : */
375 : double cellMappedVolume(const cellInfo & cell_info) const;
376 :
377 : /// Reconstruct the DAGMC geometry after skinning
378 : void reloadDAGMC();
379 :
380 : /**
381 : * Add a Filter object using the filter system.
382 : * @param[in] type the new tally type
383 : * @param[in] name the name of the new tally
384 : * @param[in] moose_object_pars the input parameters of the new tally
385 : */
386 : void addFilter(const std::string & type,
387 : const std::string & name,
388 : InputParameters & moose_object_pars);
389 :
390 : /**
391 : * Add a Tally object using the tally system.
392 : * @param[in] type the new tally type
393 : * @param[in] name the name of the new tally
394 : * @param[in] moose_object_pars the input parameters of the new tally
395 : */
396 : void
397 : addTally(const std::string & type, const std::string & name, InputParameters & moose_object_pars);
398 :
399 : /**
400 : * Multiplier on the normalized tally results; for fixed source runs,
401 : * we multiply the tally (which has units of eV/source)
402 : * by the source strength and the eV to joule conversion, while for k-eigenvalue runs, we
403 : * multiply the normalized tally (which is unitless and has an integral
404 : * value of 1.0) by the power.
405 : * @param[in] global_score tally score
406 : */
407 : Real tallyMultiplier(unsigned int global_score) const;
408 :
409 : /**
410 : * Check whether a vector extracted with getParam is empty
411 : * @param[in] vector vector
412 : * @param[in] name name to use for printing error if empty
413 : */
414 : template <typename T>
415 4520 : void checkEmptyVector(const std::vector<T> & vector, const std::string & name) const
416 : {
417 4404 : if (vector.empty())
418 8 : mooseError(name + " cannot be empty!");
419 4512 : }
420 :
421 3548 : int fixedPointIteration() const { return _fixed_point_iteration; }
422 :
423 : /**
424 : * Checks if the problem uses adaptivity or not.
425 : * @return if the problem uses adaptivity.
426 : */
427 2540 : bool hasAdaptivity() const { return _has_adaptivity; }
428 :
429 : /// Constant flag to indicate that a cell/element was unmapped
430 : static constexpr int32_t UNMAPPED{-1};
431 :
432 : /// Spatial dimension of the Monte Carlo problem
433 : static constexpr int DIMENSION{3};
434 :
435 : /// Get a modifyable non-const reference to the Moose mesh
436 : virtual MooseMesh & getMooseMesh();
437 :
438 : /// Get a modifyable const reference to the Moose mesh
439 : virtual const MooseMesh & getMooseMesh() const;
440 :
441 : /**
442 : * Whether a moving mesh is used
443 : * @return whether the [Mesh] is moving
444 : */
445 : const bool & useDisplaced() const { return _use_displaced; }
446 :
447 : protected:
448 : /**
449 : * Get the cell level in OpenMC to use for coupling
450 : * @param[in] c point
451 : * @return cell level
452 : */
453 : unsigned int getCellLevel(const Point & c) const;
454 :
455 : /**
456 : * Read the names of the MOOSE variables used for sending feedback into OpenMC
457 : * @param[in] param feedback term to read
458 : * @param[in] default_name default name to use for MOOSE variables holding this field
459 : * @param[out] vars_to_specified_blocks map from MOOSE variable names to the blocks on which they
460 : * are defined
461 : * @param[out] specified_blocks user-specified blocks for feedback
462 : */
463 : void
464 : readBlockVariables(const std::string & param,
465 : const std::string & default_name,
466 : std::map<std::string, std::vector<SubdomainName>> & vars_to_specified_blocks,
467 : std::vector<SubdomainID> & specified_blocks);
468 :
469 : /**
470 : * Whether this cell has an identical fill
471 : * @param[in] cell_info cell index, instance pair
472 : * @return whether this cell has an identical fill
473 : */
474 : bool cellHasIdenticalFill(const cellInfo & cell_info) const;
475 :
476 : /**
477 : * When using the 'identical_cell_fills' feature, this is used to determine the
478 : * contained material cells in each parent cell by applying a uniform shift
479 : * @param[in] cell_info cell index, instance pair
480 : * @return material cells contained within the given cell
481 : */
482 : containedCells shiftCellInstances(const cellInfo & cell_info) const;
483 :
484 : /**
485 : * Whether this cell overlaps with ANY value in the given subdomain set
486 : * @param[in] cell_info cell index, instance pair
487 : * @param[in] id subdomain IDs
488 : * @return whether the cell overlaps with the subdomain
489 : */
490 : bool cellMapsToSubdomain(const cellInfo & cell_info,
491 : const std::unordered_set<SubdomainID> & id) const;
492 :
493 : /**
494 : * Get the first material cell contained in the given cell
495 : * @param[in] cell_info cell index, instance pair
496 : * @return material cell index, instance pair
497 : */
498 : cellInfo firstContainedMaterialCell(const cellInfo & cell_info) const;
499 :
500 : /**
501 : * Get all of the material cells contained within this cell
502 : * @param[in] cell_info cell index, instance pair
503 : * @return all material cells contained in the given cell
504 : */
505 : containedCells containedMaterialCells(const cellInfo & cell_info) const;
506 :
507 : /**
508 : * Delete the OpenMC DAGMC geometry and re-generate the CSG geometry data structures in-place.
509 : */
510 : void updateOpenMCGeometry();
511 :
512 : /**
513 : * Re-generate the OpenMC materials in-place, needed for skinning operation where
514 : * we create new OpenMC materials on-the-fly in order to receive density feedback.
515 : */
516 : virtual void updateMaterials();
517 :
518 : /**
519 : * Get a list of each material in the problem, sorted by subdomain. This function also checks
520 : * that there is just one OpenMC material in each subdomain, necessary for the DAGMC skinning.
521 : * @return material in each subdomain
522 : */
523 : std::vector<std::string> getMaterialInEachSubdomain() const;
524 :
525 : /**
526 : * Apply transformations and scale point from MOOSE into the OpenMC domain
527 : * @param[in] pt point
528 : * @return transformed point
529 : */
530 : Point transformPointToOpenMC(const Point & pt) const;
531 :
532 : /**
533 : * Check that the tally normalization gives a total tally sum of 1.0 (when normalized
534 : * against the total tally value).
535 : * @param[in] sum sum of the tally
536 : * @param[in] score tally score
537 : */
538 : void checkNormalization(const Real & sum, unsigned int global_score) const;
539 :
540 : /**
541 : * For geometries with fine-scale details (e.g. TRISO), Cardinal's default settings can
542 : * take a very long time to initialize the problem (but we can't change those defaults
543 : * because they are not 100% applicable all the time). So, we print out a message to
544 : * the user to point them in the right direction if their initialization is taking a
545 : * long time.
546 : * @param[in] start time to use for evaluating whether we've exceeded our limit for printing the
547 : * message
548 : */
549 : void
550 : printTrisoHelp(const std::chrono::time_point<std::chrono::high_resolution_clock> & start) const;
551 :
552 : /**
553 : * Print to the console the names of the auxvariables used for I/O with OpenMC.
554 : * We only print these tables once, upon initialization, because this data does
555 : * not change if re-initializing the spatial mapping for moving-mesh problems,
556 : * adaptive refinement, skinning, etc.
557 : */
558 : void printAuxVariableIO();
559 :
560 : /**
561 : * Get all the material indices within the set of cells
562 : * @param[in] contained_cells set of cells
563 : * @return contained materials
564 : */
565 : std::vector<int32_t> materialsInCells(const containedCells & contained_cells) const;
566 :
567 : /// Loop over the mapped cells, and build a map between subdomains to OpenMC materials
568 : void subdomainsToMaterials();
569 :
570 : /**
571 : * Get a set of all subdomains that have at least 1 element coupled to an OpenMC cell
572 : * @return subdomains with at least 1 element coupled to OpenMC
573 : */
574 : std::set<SubdomainID> coupledSubdomains() const;
575 :
576 : /**
577 : * Gather a vector of values to be summed for each cell
578 : * @param[in] local local values to be summed for the cells
579 : * @param[out] global global mapping of the summed values to the cells
580 : */
581 : template <typename T>
582 : void gatherCellSum(std::vector<T> & local, std::map<cellInfo, T> & global) const;
583 :
584 : /**
585 : * Gather a vector of values to be pushed back to for each cell
586 : * @param[in] local local values to be pushed back for the cells
587 : * @param[in] n_local number of local values contributed to each cell
588 : * @param[out] global global mapping of the pushed back values to the cells
589 : */
590 : template <typename T>
591 : void gatherCellVector(std::vector<T> & local, std::vector<unsigned int> & n_local, std::map<cellInfo, std::vector<T>> & global);
592 :
593 : /**
594 : * Get the feedback which this element provides to OpenMC
595 : * @param[in] elem
596 : * @return coupling phase
597 : */
598 : coupling::CouplingFields elemFeedback(const Elem * elem) const;
599 :
600 : /**
601 : * Read the parameters needed for triggers
602 : * @param[in] params input parameters
603 : */
604 : void getTallyTriggerParameters(const InputParameters & params);
605 :
606 : /**
607 : * Read the block parameters based on user settings
608 : * @param[in] name name of input parameter representing a vector of subdomain names
609 : * @param[in] blocks list of block ids to write
610 : */
611 : void readBlockParameters(const std::string name, std::unordered_set<SubdomainID> & blocks);
612 :
613 : /**
614 : * Cache the material cells contained within each coupling cell;
615 : * depending on user settings, this may attempt to take shortcuts
616 : * by assuming each cell has the same fills
617 : */
618 : void cacheContainedCells();
619 :
620 : /**
621 : * Fill the cached contained cells data structure for a given cell
622 : * @param[in] cell_info cell index, instance pair
623 : * @param[in] hint location hint used to accelerate the search
624 : * @param[out] map contained cell map
625 : */
626 : void setContainedCells(const cellInfo & cell_info,
627 : const Point & hint,
628 : std::map<cellInfo, containedCells> & map);
629 :
630 : /**
631 : * Check that the structure of the contained material cells for two cell matches;
632 : * i.e. this checks that the keys are the same and that the *number* of instances
633 : * of each filling material cell match.
634 : * @param[in] cell_info cell index, instance pair
635 : * @param[in] reference map we want to check against
636 : * @param[in] compare map we want to check
637 : */
638 : void checkContainedCellsStructure(const cellInfo & cell_info,
639 : containedCells & reference,
640 : containedCells & compare) const;
641 :
642 : /**
643 : * Set a minimum order for a volume quadrature rule
644 : * @param[in] volume_order order of the volume quadrature rule
645 : * @param[in] type string type of quadrature rule for printing a console message
646 : */
647 : void setMinimumVolumeQRules(Order & volume_order, const std::string & type);
648 :
649 : /// For keeping the output neat when using verbose
650 : std::string printNewline() const
651 : {
652 : if (_verbose)
653 : return "\n";
654 : else
655 : return "";
656 : }
657 :
658 : /// Loop over the elements in the MOOSE mesh and store the type of feedback applied by each.
659 : void storeElementPhase();
660 :
661 : /**
662 : * Relax the tally and normalize it according to some "global" tally. If you set
663 : * 'normalize_by_global_tally = true', you will be normalizing by a tally over the ENTIRE
664 : * OpenMC geometry. Otherwise, you normalize only by the sum of the tally bins themselves.
665 : *
666 : * NOTE: This function relaxes the tally _distribution_, and not the actual magnitude of the sum.
667 : * That is, we relax the shape distribution and then multiply it by the power
668 : * (for k-eigenvalue) or source strength (for fixed source) of the current step before
669 : * applying it to MOOSE. If the magnitude of the power or source strength is constant in time,
670 : * there is zero error in this. But if the magnitude of the tally varies in time, we are basically
671 : * relaxing the distribution of the tally, but then multiplying it by the _current_ mean tally
672 : * magnitude.
673 : *
674 : * There will be very small errors in these approximations unless the power/source strength
675 : * change dramatically with iteration. But because relaxation is itself a numerical approximation,
676 : * this is still inconsequential at the end of the day as long as your problem has converged
677 : * the relaxed tally to the raw (unrelaxed) tally.
678 : * @param[in] global_score the global index of the tally score
679 : * @param[in] local_score the local index of the tally score
680 : * @param[in] local_tally the tally to relax and normalize
681 : */
682 : void relaxAndNormalizeTally(unsigned int global_score,
683 : unsigned int local_score,
684 : std::shared_ptr<TallyBase> local_tally);
685 :
686 : /**
687 : * Loop over all the OpenMC cells and count the number of MOOSE elements to which the cell
688 : * is mapped based on phase.
689 : */
690 : void getCellMappedPhase();
691 :
692 : /// This function is used to ensure that each OpenMC cell only maps to a single phase
693 : void checkCellMappedPhase();
694 :
695 : /// Loop over all the OpenMC cells and get the element subdomain IDs that map to each cell
696 : void getCellMappedSubdomains();
697 :
698 : /**
699 : * Loop over all the OpenMC cells and compute the volume of the MOOSE elements that each
700 : * cell maps to
701 : */
702 : void computeCellMappedVolumes();
703 :
704 : /// Set up the mapping from MOOSE elements to OpenMC cells
705 : void initializeElementToCellMapping();
706 :
707 : /// Populate maps of MOOSE elements to OpenMC cells
708 : void mapElemsToCells();
709 :
710 : /**
711 : * A function which validates local tallies. This is done to ensure that at least one of the
712 : * tallies contains a heating score when running in eigenvalue mode. This must be done outside
713 : * of the constructor as tallies are added from an external system.
714 : */
715 : void validateLocalTallies();
716 :
717 : /// Add OpenMC tallies to facilitate the coupling
718 : void initializeTallies();
719 :
720 : /**
721 : * Reset any tallies previously added by Cardinal, by deleting them from OpenMC.
722 : * Also delete any mesh filters and meshes added to OpenMC for mesh filters.
723 : */
724 : void resetTallies();
725 :
726 : /// Find the material filling each cell which receives density feedback
727 : void getMaterialFills();
728 :
729 : /**
730 : * Get one point inside each cell, for accelerating the particle search routine.
731 : * This function will get the centroid of the first global element in the lowest
732 : * rank in the cell.
733 : */
734 : void getPointInCell();
735 :
736 : /**
737 : * Compute the product of volume with a field across ranks and sum into a global map
738 : * @param[in] var_num variable to weight with volume, mapped by subdomain ID
739 : * @param[in] phase phases to compute the operation for
740 : * @return volume-weighted field for each cell, in a global sense
741 : */
742 : std::map<cellInfo, Real> computeVolumeWeightedCellInput(
743 : const std::map<SubdomainID, std::pair<unsigned int, std::string>> & var_num,
744 : const std::vector<coupling::CouplingFields> * phase) const;
745 :
746 : /**
747 : * Send temperature from MOOSE to OpenMC by computing a volume average
748 : * and applying a single temperature per OpenMC cell
749 : */
750 : void sendTemperatureToOpenMC() const;
751 :
752 : /**
753 : * Send density from MOOSE to OpenMC by computing a volume average
754 : * and applying a single density per OpenMC cell.
755 : */
756 : void sendDensityToOpenMC() const;
757 :
758 : /**
759 : * Factor by which to normalize a tally
760 : * @param[in] global_score global index for the tally score
761 : * @return value to divide tally sum by for normalization
762 : */
763 : Real tallyNormalization(unsigned int global_score) const;
764 :
765 : /**
766 : * Check the sum of the tallies against the global tally
767 : * @param[in] score tally score
768 : */
769 : void checkTallySum(const unsigned int & score) const;
770 :
771 : /**
772 : * Check if a mapped location is in the outer universe of a lattice
773 : * @param[in] level lattice level
774 : * @return whether the location is in the outer universe
775 : */
776 : void latticeOuterCheck(const Point & c, int level) const;
777 :
778 : /**
779 : * Report an error for a mapped location in an outer universe of a lattice
780 : * @param[in] c Mapped location
781 : * @param[in] level level of the mapped cell
782 : */
783 : void latticeOuterError(const Point & c, int level) const;
784 :
785 : /**
786 : * Find the OpenMC cell at a given point in space
787 : * @param[in] point point
788 : * @return whether OpenMC reported an error
789 : */
790 : bool findCell(const Point & point);
791 :
792 : /**
793 : * Checks that the contained material cells exactly match between a reference obtained
794 : * by calling openmc::Cell::get_contained_cells for each cell and a shortcut
795 : * approach that assumes all identical cells (which aren't simply just material fills)
796 : * has exactly the same contained material cells.
797 : * @param[in] reference reference map to compare against
798 : * @param[in] compare shortcut map to compare
799 : */
800 : void compareContainedCells(std::map<cellInfo, containedCells> & reference,
801 : std::map<cellInfo, containedCells> & compare) const;
802 :
803 : NumericVector<Number> & _serialized_solution;
804 :
805 : /**
806 : * Return all IDs of all Cardinal-mapped Tallies
807 : * @return all Cardinal-mapped Tally IDs
808 : */
809 : virtual std::vector<int32_t> getMappedTallyIDs() const override;
810 :
811 : /**
812 : * Whether to automatically compute the mapping of OpenMC cell IDs and
813 : * instances to the [Mesh].
814 : */
815 : const bool & _output_cell_mapping;
816 :
817 : /**
818 : * Where to get the initial OpenMC temperatures and densities from;
819 : * can be either hdf5 (from a properties.h5 file), xml (whatever is already
820 : * set in the XML files), or moose (meaning whatever ICs are set on the 'temperature_variables'
821 : * and 'density_variables'
822 : */
823 : const coupling::OpenMCInitialCondition _initial_condition;
824 :
825 : /// Type of relaxation to apply to the OpenMC tallies
826 : const relaxation::RelaxationEnum _relaxation;
827 :
828 : /**
829 : * Type of trigger to apply to k eigenvalue to indicate when
830 : * the simulation is complete. These can be used to on-the-fly adjust the number
831 : * of active batches in order to reach some desired criteria (which is specified
832 : * by this parameter).
833 : */
834 : const trigger::TallyTriggerTypeEnum _k_trigger;
835 :
836 : /**
837 : * Coordinate level in the OpenMC domain to use for mapping cells to mesh.
838 : * When using 'lowest_cell_level', this parameter indicates that the lowest
839 : * cell level is used, up until _cell_level.
840 : */
841 : unsigned int _cell_level;
842 :
843 : /**
844 : * Whether OpenMC properties (temperature and density) should be exported
845 : * after being updated in syncSolutions.
846 : */
847 : const bool & _export_properties;
848 :
849 : /**
850 : * How to normalize the OpenMC tally into units of W/volume. If 'true',
851 : * normalization is performed by dividing each local tally against a problem-global
852 : * tally. The advantage of this approach is that some non-zero tally regions of the
853 : * OpenMC domain can be excluded from multiphysics feedback (without us having to guess
854 : * what the power of the *included* part of the domain is). This can let us do
855 : * "zooming" type calculations, where perhaps we only want to send T/H feedback to
856 : * one bundle in a full core.
857 : *
858 : * If 'false', normalization is performed by dividing each local tally by the sum
859 : * of the local tally itself. The advantage of this approach becomes evident when
860 : * using mesh tallies. If a mesh tally does not perfectly align with an OpenMC cell -
861 : * for instance, a first-order sphere mesh will not perfectly match the volume of a
862 : * TRISO pebble - then not all of the power actually produced in the pebble is
863 : * tallies on the mesh approximation to that pebble. Therefore, if you set a core
864 : * power of 1 MW and you normalized based on a global tally, you'd always
865 : * miss some of that power when sending to MOOSE. So, in this case, it is better to
866 : * normalize against the local tally itself so that the correct power is preserved.
867 : */
868 : const bool _normalize_by_global;
869 :
870 : /// Whether or not the problem uses a skinner to regenerate the OpenMC geometry.
871 : const bool _using_skinner;
872 :
873 : /**
874 : * When the mesh changes during the simulation (either from adaptive mesh refinement
875 : * or deformation), the mapping from OpenMC cells to the [Mesh] must be re-established
876 : * after each OpenMC run.
877 : */
878 : bool _need_to_reinit_coupling;
879 :
880 : /**
881 : * Whether to check the tallies against the global tally;
882 : * if set to true, and the tallies added for the 'tally_blocks' do not
883 : * sum to the global tally, an error is thrown. If you are
884 : * only performing multiphysics feedback for, say, a single assembly in a
885 : * full-core OpenMC model, you must set this check to false, because there
886 : * are known fission sources outside the domain of interest.
887 : *
888 : * If not specified, then this is set to 'true' if normalizing by a global
889 : * tally, and to 'false' if normalizing by the local tally (because when we choose
890 : * to normalize by the local tally, we're probably using mesh tallies). But you can
891 : * of course still set a value for this parameter to override the default.
892 : */
893 : const bool _check_tally_sum;
894 :
895 : /// Constant relaxation factor
896 : const Real & _relaxation_factor;
897 :
898 : /**
899 : * If known a priori by the user, whether the tally cells (which are not simply material
900 : * fills) have EXACTLY the same contained material cells. This is a big optimization for
901 : * TRISO problems in setting up homogenized temperature/density feedback to OpenMC.
902 : *
903 : * The concept can best be explained with a pebble bed reactor.
904 : * If every pebble is filled with an identical TRISO universe, then the material fills
905 : * in each pebble are identical to one another except for a constant offset. This idea
906 : * can be used to then skip all but the first two openmc::Cell::get_contained_cells
907 : * calls (which are required in order to figure out the pattern by which pebble N is
908 : * incremented relative to pebble 1).
909 : *
910 : * When using this parameter, we HIGHLY recommend setting 'check_identical_cell_fills =
911 : * true' the first time you run your model. This will figure out the material cell fills using a
912 : * method that calls openmc::Cell::get_contained_cells for every tally cell, i.e. without assuming
913 : * anything about repeated structure in your OpenMC model. Setting 'identical_cell_fills'
914 : * without also setting 'check_identical_cell_fills = true' may result in SILENT
915 : * errors!!! So it is essential to be sure you've removed any error sources before you turn the
916 : * error check off to actually leverage the speedup.
917 : *
918 : * Note: for any tally cells that are just filled with a material, we use the approach
919 : * where openmc::Cell::get_contained_cells is called in full.
920 : *
921 : * This optimization will not work (and 'check_identical_cell_fills = true' *will*
922 : * catch these) for:
923 : * - any situation where tallied, non-material-fill pebbles have different fills
924 : * (such as if you have different TRISO lattices in each pebble)
925 : * - any situation where there is a "gap" in the incrementing of the material fill
926 : * instances (such as if pebble 89 does not map to 'tally_blocks', then the instance
927 : * shift for pebble 90 relative to pebble 1 is 89, when it should have been 90).
928 : */
929 : const bool _has_identical_cell_fills;
930 :
931 : /**
932 : * Whether we should rigorously check that each tally cell has identical fills;
933 : * this is SLOW for large TRISO problems, but is essential to ensure the accuracy of
934 : * 'identical_cell_fills'. Please set 'check_identical_cell_fills' to 'true' at least
935 : * once before running production cases to be sure the optimization can be applied.
936 : */
937 : const bool & _check_identical_cell_fills;
938 :
939 : /**
940 : * Whether it can be assumed that all of the tallies (both those set by the user
941 : * in the XML file, as well as those created automatically by Cardinal) are
942 : * spatially separate. This means that once a particle scores to one tally bin, it wouldn't
943 : * score to ANY other tally bins. This can dramatically increase tracking rates
944 : * for problems with many tallies.
945 : */
946 : const bool & _assume_separate_tallies;
947 :
948 : /**
949 : * Whether to map density according to each individual OpenMC cell (in which case an
950 : * error is thrown if you don't have a unique material in each cell) or by material.
951 : */
952 : bool _map_density_by_cell;
953 :
954 : /**
955 : * Whether the problem has density feedback blocks specified; note that this is NOT necessarily
956 : * indicative that the mapping was successful in finding any cells corresponding to those blocks
957 : */
958 : const bool _specified_density_feedback;
959 :
960 : /**
961 : * Whether the problem has temperature feedback blocks specified; note that this is NOT
962 : * necessarily indicative that the mapping was successful in finding any cells corresponding to
963 : * those blocks
964 : */
965 : const bool _specified_temperature_feedback;
966 :
967 : /// Whether any cell tallies exist.
968 : bool _has_cell_tallies = false;
969 :
970 : /// Whether any spatial mapping from OpenMC's cells to the mesh is needed
971 : bool _needs_to_map_cells;
972 :
973 : /**
974 : * Whether a global tally is required for the sake of normalization and/or checking
975 : * the tally sum
976 : */
977 : const bool _needs_global_tally;
978 :
979 : /**
980 : * A map of the filter objects created by the [Problem/Filters] block. The key for each filter is
981 : * it's corresponding MOOSE name to allow tallies to look up filters.
982 : */
983 : std::map<std::string, std::shared_ptr<FilterBase>> _filters;
984 :
985 : /// A vector of the tally objects created by the [Problem/Tallies] block.
986 : std::vector<std::shared_ptr<TallyBase>> _local_tallies;
987 :
988 : /// A list of all of the scores contained by the local tallies added in the [Tallies] block.
989 : std::vector<std::string> _all_tally_scores;
990 :
991 : /**
992 : * The [Tallies] block allows tallies with different scores, and so we can't assume they have the
993 : * same indices in each tally's arrays. This variable map between the name of each score and it's
994 : * index in each local tally.
995 : */
996 : std::vector<std::map<std::string, int>> _local_tally_score_map;
997 :
998 : /// A vector of auxvariable ids added by the [Tallies] block.
999 : std::vector<std::vector<unsigned int>> _tally_var_ids;
1000 :
1001 : /// A vector of external (output-based) auxvariable ids added by the [Tallies] block.
1002 : std::vector<std::vector<std::vector<unsigned int>>> _tally_ext_var_ids;
1003 :
1004 : /// Blocks in MOOSE mesh that provide density feedback
1005 : std::vector<SubdomainID> _density_blocks;
1006 :
1007 : /// Blocks in MOOSE mesh that provide temperature feedback
1008 : std::vector<SubdomainID> _temp_blocks;
1009 :
1010 : /// Blocks for which the cell fills are identical
1011 : std::unordered_set<SubdomainID> _identical_cell_fill_blocks;
1012 :
1013 : /// Mapping of MOOSE elements to the OpenMC cell they map to (if any)
1014 : std::vector<cellInfo> _elem_to_cell{};
1015 :
1016 : /// Phase of each cell
1017 : std::map<cellInfo, coupling::CouplingFields> _cell_phase;
1018 :
1019 : /// Number of elements in the MOOSE mesh that exclusively provide density feedback
1020 : int _n_moose_density_elems;
1021 :
1022 : /// Number of elements in the MOOSE mesh that exclusively provide temperature feedback
1023 : int _n_moose_temp_elems;
1024 :
1025 : /// Number of elements in the MOOSE mesh which provide temperature+density feedback
1026 : int _n_moose_temp_density_elems;
1027 :
1028 : /// Number of no-coupling elements in the MOOSE mesh
1029 : int _n_moose_none_elems;
1030 :
1031 : /**
1032 : * Number of MOOSE elements that exclusively provide temperature feedback,
1033 : * and which successfully mapped to OpenMC cells
1034 : */
1035 : int _n_mapped_temp_elems;
1036 :
1037 : /**
1038 : * Number of MOOSE elements that exclusively provide density feedback,
1039 : * and which successfully mapped to OpenMC cells
1040 : */
1041 : int _n_mapped_density_elems;
1042 :
1043 : /**
1044 : * Number of MOOSE elements that provide temperature+density feedback,
1045 : * and which successfully mapped to OpenMC cells
1046 : */
1047 : int _n_mapped_temp_density_elems;
1048 :
1049 : /// Number of no-coupling elements mapped to OpenMC cells
1050 : int _n_mapped_none_elems;
1051 :
1052 : /// Total volume of uncoupled MOOSE mesh elements
1053 : Real _uncoupled_volume;
1054 :
1055 : /// Whether non-material cells are mapped
1056 : bool _material_cells_only{true};
1057 :
1058 : /// Mapping of OpenMC cell indices to a vector of MOOSE element IDs
1059 : std::map<cellInfo, std::vector<unsigned int>> _cell_to_elem;
1060 :
1061 : /// Mapping of OpenMC cell indices to a vector of MOOSE element IDs, on each local rank
1062 : std::map<cellInfo, std::vector<unsigned int>> _local_cell_to_elem;
1063 :
1064 : /// Mapping of OpenMC cell indices to the subdomain IDs each maps to
1065 : std::map<cellInfo, std::unordered_set<SubdomainID>> _cell_to_elem_subdomain;
1066 :
1067 : /// Mapping of elem subdomains to materials
1068 : std::map<SubdomainID, std::set<int32_t>> _subdomain_to_material;
1069 :
1070 : /**
1071 : * A point inside the cell, taken simply as the centroid of the first global
1072 : * element inside the cell. This is stored to accelerate the particle search.
1073 : */
1074 : std::map<cellInfo, Point> _cell_to_point;
1075 :
1076 : /**
1077 : * Volume associated with the mapped element space for each OpenMC cell; the unit
1078 : * for this volume is whatever is used in the [Mesh] block
1079 : */
1080 : std::map<cellInfo, Real> _cell_to_elem_volume;
1081 :
1082 : /**
1083 : * Volume associated with the actual OpenMC cell, computed by an optional
1084 : * OpenMCVolumeCalculation user object
1085 : */
1086 : std::map<cellInfo, Real> _cell_volume;
1087 :
1088 : /**
1089 : * Material filling each cell to receive density feedback. We enforce that these
1090 : * cells are filled with a material (cannot be filled with a lattice or universe).
1091 : */
1092 : std::map<cellInfo, int32_t> _cell_to_material;
1093 :
1094 : /**
1095 : * Material-type cells contained within a cell; this is only populated if a cell
1096 : * is NOT indicated as having an identical fill
1097 : */
1098 : std::map<cellInfo, containedCells> _cell_to_contained_material_cells;
1099 :
1100 : /// Number of material-type cells contained within a cell
1101 : std::map<cellInfo, int32_t> _cell_to_n_contained;
1102 :
1103 : /**
1104 : * Global tallies. We add one per tally added in the [Tallies] block to
1105 : * enable global noramlization.
1106 : */
1107 : std::vector<openmc::Tally *> _global_tallies;
1108 :
1109 : /// Global tally scores corresponding to '_global_tallies'.
1110 : std::vector<std::vector<std::string>> _global_tally_scores;
1111 :
1112 : /// Global tally estimators corresponding to '_global_tallies'.
1113 : std::vector<openmc::TallyEstimator> _global_tally_estimators;
1114 :
1115 : /// Sum value of the global tally(s), across all bins
1116 : std::vector<Real> _global_sum_tally;
1117 :
1118 : /// Sum value of the local tally(s), across all bins
1119 : std::vector<Real> _local_sum_tally;
1120 :
1121 : /// Mean value of the local tally(s), across all bins; only used for fixed source mode
1122 : std::vector<Real> _local_mean_tally;
1123 :
1124 : /// Whether the present transfer is the first transfer
1125 : static bool _first_transfer;
1126 :
1127 : /// Whether the diagnostic tables on initialization have already been printed
1128 : static bool _printed_initial;
1129 :
1130 : /// Whether a warning has already been printed about very long setup times (for TRISOs)
1131 : static bool _printed_triso_warning;
1132 :
1133 : /// Dummy particle to reduce number of allocations of particles for cell lookup routines
1134 : openmc::Particle _particle;
1135 :
1136 : /// Number of particles simulated in the first iteration
1137 : unsigned int _n_particles_1;
1138 :
1139 : /// Mapping from temperature variable name to the subdomains on which to read it from
1140 : std::map<std::string, std::vector<SubdomainName>> _temp_vars_to_blocks;
1141 :
1142 : /// Mapping from density variable name to the subdomains on which to read it from
1143 : std::map<std::string, std::vector<SubdomainName>> _density_vars_to_blocks;
1144 :
1145 : /// Optional volume calculation for cells which map to MOOSE
1146 : OpenMCVolumeCalculation * _volume_calc;
1147 :
1148 : /// Userobject that maps from a partial-symmetry OpenMC model to a whole-domain [Mesh]
1149 : const SymmetryPointGenerator * _symmetry;
1150 :
1151 : /// Number of temperature-only feedback elements in each mapped OpenMC cell (global)
1152 : std::map<cellInfo, int> _n_temp;
1153 :
1154 : /// Number of density-only feedback elements in each mapped OpenMC cell (global)
1155 : std::map<cellInfo, int> _n_rho;
1156 :
1157 : /// Number of temperature+density feedback elements in each mapped OpenMC cell (global)
1158 : std::map<cellInfo, int> _n_temp_rho;
1159 :
1160 : /// Number of none elements in each mapped OpenMC cell (global)
1161 : std::map<cellInfo, int> _n_none;
1162 :
1163 : /// Index in OpenMC tallies corresponding to the first global tally added by Cardinal
1164 : unsigned int _global_tally_index = 0;
1165 :
1166 : /// Index in tally_score pointing to the score used for normalizing flux tallies in eigenvalue mode
1167 : unsigned int _source_rate_index;
1168 :
1169 : #ifdef ENABLE_DAGMC
1170 : /// Optional skinner to re-generate the OpenMC geometry on-the-fly for DAGMC models
1171 : MoabSkinner * _skinner = nullptr;
1172 :
1173 : /// Pointer to DAGMC
1174 : std::shared_ptr<moab::DagMC> _dagmc = nullptr;
1175 : #endif
1176 :
1177 : /// Total number of unique OpenMC cell IDs + instances combinations
1178 : long unsigned int _n_openmc_cells;
1179 :
1180 : /// ID of the OpenMC universe corresponding to the DAGMC universe
1181 : int32_t _dagmc_universe_id;
1182 :
1183 : /// Whether the DAGMC universe is the root universe or not.
1184 : bool _dagmc_root_universe = true;
1185 :
1186 : /// ID of the OpenMC cell corresponding to the cell which uses the DAGMC universe as a fill.
1187 : int32_t _cell_using_dagmc_universe_id;
1188 :
1189 : /// The number of OpenMC surfaces before skinning occurs. This is required to properly reinitialize
1190 : /// the CSG geometry contained in the OpenMC model.
1191 : const int32_t _initial_num_openmc_surfaces;
1192 :
1193 : /// Conversion rate from eV to Joule
1194 : static constexpr Real EV_TO_JOULE = 1.6022e-19;
1195 :
1196 : /// Tolerance for setting zero tally
1197 : static constexpr Real ZERO_TALLY_THRESHOLD = 1e-12;
1198 :
1199 : private:
1200 : /**
1201 : * Update the number of particles according to the Dufek-Gudowski relaxation scheme
1202 : */
1203 : void dufekGudowskiParticleUpdate();
1204 :
1205 : /// Flattened cell IDs collected after parallel communication
1206 : std::vector<int32_t> _flattened_ids;
1207 :
1208 : /// Flattened cell instancess collected after parallel communication
1209 : std::vector<int32_t> _flattened_instances;
1210 :
1211 : /// Offsets for each cell instance in an identically-repeated universe
1212 : containedCells _instance_offsets;
1213 :
1214 : /// Offset for each cell relative to the first identical-fill cell
1215 : std::map<cellInfo, int32_t> _n_offset;
1216 :
1217 : /// First identical-fill cell
1218 : cellInfo _first_identical_cell;
1219 :
1220 : /// Materials in the first identical-fill cell
1221 : std::vector<int32_t> _first_identical_cell_materials;
1222 :
1223 : /// Whether OpenMCCellAverageProblem should use the displaced mesh
1224 : bool _use_displaced;
1225 :
1226 : /// Mapping from subdomain IDs to which aux variable to read temperature (K) from
1227 : std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_temp_vars;
1228 :
1229 : /// Mapping from subdomain IDs to which aux variable to read density (kg/m3) from
1230 : std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_density_vars;
1231 : };
|