Cardinal
OpenMCCellAverageProblem.h
Go to the documentation of this file.
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 
25 #include "TallyBase.h"
26 #include "FilterBase.h"
27 
28 #ifdef ENABLE_DAGMC
29 #include "MoabSkinner.h"
30 #include "DagMC.hpp"
31 #endif
32 
35 
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  virtual bool converged(unsigned int) override { return true; }
76 
84  void read2DBlockParameters(const std::string name,
85  std::vector<std::vector<SubdomainName>> & names,
86  std::vector<SubdomainID> & flattened_ids);
87 
94  void checkBlocksInMesh(const std::string name,
95  const std::vector<SubdomainID> & ids,
96  const std::vector<SubdomainName> & names) const;
97 
99  void setupProblem();
100 
105  virtual void addExternalVariables() override;
106 
112  virtual Real cellVolume(const cellInfo & cell_info) const;
113 
118  virtual const OpenMCVolumeCalculation * volumeCalculation() const { return _volume_calc; }
119 
124  virtual const std::map<cellInfo, std::vector<unsigned int>> & cellToElem() const
125  {
126  return _cell_to_elem;
127  }
128 
134  virtual std::unordered_set<SubdomainID> getCellToElementSub(const cellInfo & info)
135  {
136  return _cell_to_elem_subdomain.at(info);
137  }
138 
143  virtual bool hasPointTransformations() const { return _symmetry != nullptr; }
144 
149  virtual const std::vector<std::string> & getTallyScores() const { return _all_tally_scores; }
150 
156  unsigned int getNumScoringTallies(const std::string & score) const
157  {
158  return _score_count.count(score) > 0 ? _score_count.at(score) : 0;
159  }
160 
166  bool hasScore(const std::string & score)
167  {
168  return std::find(_all_tally_scores.begin(), _all_tally_scores.end(), score) !=
169  _all_tally_scores.end();
170  }
171 
177  const TallyBase * getTally(const std::string & name);
178 
189  std::vector<const MooseVariableFE<Real> *> getTallyScoreVariables(const std::string & score,
190  const std::string & tally_name,
191  THREAD_ID tid,
192  const std::string & output = "",
193  bool skip_func_exp = false);
194 
205  std::vector<const VariableValue *> getTallyScoreVariableValues(const std::string & score,
206  const std::string & tally_name,
207  THREAD_ID tid,
208  const std::string & output = "",
209  bool skip_func_exp = false);
210 
221  std::vector<const VariableValue *>
222  getTallyScoreNeighborVariableValues(const std::string & score,
223  const std::string & tally_name,
224  THREAD_ID tid,
225  const std::string & output = "",
226  bool skip_func_exp = false);
227 
235  bool hasOutput(const std::string & score, const std::string & output) const;
236 
242  virtual Point transformPoint(const Point & pt) const
243  {
244  return this->hasPointTransformations() ? _symmetry->transformPoint(pt) : pt;
245  }
246 
268  virtual void createQRules(QuadratureType type,
269  Order order,
270  Order volume_order,
271  Order face_order,
272  SubdomainID block,
273  bool allow_negative_weights = true) override;
274 
279  typedef std::unordered_map<int32_t, std::vector<int32_t>> containedCells;
280 
286  int32_t elemToCellIndex(const int & elem_id) const { return elemToCellInfo(elem_id).first; }
287 
294  int32_t elemToCellID(const int & elem_id) const { return cellID(elemToCellIndex(elem_id)); }
295 
301  int32_t elemToCellInstance(const int & elem_id) const { return elemToCellInfo(elem_id).second; }
302 
309  cellInfo elemToCellInfo(const int & elem_id) const { return _elem_to_cell[elem_id]; }
310 
318  int32_t cellToMaterialIndex(const cellInfo & cell_info) const
319  {
320  return _cell_to_material.at(cell_info);
321  }
322 
331  coupling::CouplingFields cellFeedback(const cellInfo & cell_info) const;
332 
338  bool hasDensityFeedback(const cellInfo & cell_info) const
339  {
340  std::vector<coupling::CouplingFields> phase = {coupling::density,
342  return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
343  }
344 
350  bool hasTemperatureFeedback(const cellInfo & cell_info) const
351  {
352  std::vector<coupling::CouplingFields> phase = {coupling::temperature,
354  return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
355  }
356 
362  bool hasFilter(const std::string & filter_name) const { return _filters.count(filter_name) > 0; }
363 
369  std::shared_ptr<FilterBase> & getFilter(const std::string & filter_name)
370  {
371  return _filters.at(filter_name);
372  }
373 
378  const std::vector<std::shared_ptr<TallyBase>> & getLocalTallies() const { return _local_tallies; }
379 
385  double cellTemperature(const cellInfo & cell_info) const;
386 
391  double cellMappedVolume(const cellInfo & cell_info) const;
392 
394  void reloadDAGMC();
395 
402  void addFilter(const std::string & type,
403  const std::string & name,
404  InputParameters & moose_object_pars);
405 
413  std::shared_ptr<TallyBase>
414  addTally(const std::string & type, const std::string & name, InputParameters & moose_object_pars);
415 
425  Real tallyMultiplier(const std::string & score_name, const Real & local_mean_tally) const;
426 
432  template <typename T>
433  void checkEmptyVector(const std::vector<T> & vector, const std::string & name) const
434  {
435  if (vector.empty())
436  mooseError(name + " cannot be empty!");
437  }
438 
440 
445  bool hasAdaptivity() const { return _has_adaptivity; }
446 
448  static constexpr int32_t UNMAPPED{-1};
449 
451  static constexpr int DIMENSION{3};
452 
454  virtual MooseMesh & getMooseMesh();
455 
457  virtual const MooseMesh & getMooseMesh() const;
458 
463  const bool & useDisplaced() const { return _use_displaced; }
464 
465 protected:
471  unsigned int getCellLevel(const Point & c) const;
472 
481  void
482  readBlockVariables(const std::string & param,
483  const std::string & default_name,
484  std::map<std::string, std::vector<SubdomainName>> & vars_to_specified_blocks,
485  std::vector<SubdomainID> & specified_blocks);
486 
492  bool cellHasIdenticalFill(const cellInfo & cell_info) const;
493 
500  containedCells shiftCellInstances(const cellInfo & cell_info) const;
501 
508  bool cellMapsToSubdomain(const cellInfo & cell_info,
509  const std::unordered_set<SubdomainID> & id) const;
510 
516  cellInfo firstContainedMaterialCell(const cellInfo & cell_info) const;
517 
523  containedCells containedMaterialCells(const cellInfo & cell_info) const;
524 
528  void updateOpenMCGeometry();
529 
535  std::vector<std::string> getMaterialInEachSubdomain() const;
536 
542  Point transformPointToOpenMC(const Point & pt) const;
543 
553  void
554  printTrisoHelp(const std::chrono::time_point<std::chrono::high_resolution_clock> & start) const;
555 
562  void printAuxVariableIO();
563 
569  std::vector<int32_t> materialsInCells(const containedCells & contained_cells) const;
570 
572  void subdomainsToMaterials();
573 
578  std::set<SubdomainID> coupledSubdomains() const;
579 
585  template <typename T>
586  void gatherCellSum(std::vector<T> & local, std::map<cellInfo, T> & global) const;
587 
594  template <typename T>
595  void gatherCellVector(std::vector<T> & local, std::vector<unsigned int> & n_local, std::map<cellInfo, std::vector<T>> & global);
596 
602  coupling::CouplingFields elemFeedback(const Elem * elem) const;
603 
608  void getTallyTriggerParameters(const InputParameters & params);
609 
615  void readBlockParameters(const std::string name, std::unordered_set<SubdomainID> & blocks);
616 
622  void cacheContainedCells();
623 
630  void setContainedCells(const cellInfo & cell_info,
631  const Point & hint,
632  std::map<cellInfo, containedCells> & map);
633 
642  void checkContainedCellsStructure(const cellInfo & cell_info,
643  containedCells & reference,
644  containedCells & compare) const;
645 
651  void setMinimumVolumeQRules(Order & volume_order, const std::string & type);
652 
654  std::string printNewline() const
655  {
656  if (_verbose)
657  return "\n";
658  else
659  return "";
660  }
661 
663  void storeElementPhase();
664 
669  void getCellMappedPhase();
670 
672  void checkCellMappedPhase();
673 
676 
682 
685 
687  void mapElemsToCells();
688 
694  void validateLocalTallies();
695 
697  void initializeTallies();
698 
703  void resetTallies();
704 
706  void getMaterialFills();
707 
713  void getPointInCell();
714 
721  std::map<cellInfo, Real> computeVolumeWeightedCellInput(
722  const std::map<SubdomainID, std::pair<unsigned int, std::string>> & var_num,
723  const std::vector<coupling::CouplingFields> * phase) const;
724 
729  void sendTemperatureToOpenMC() const;
730 
735  void sendDensityToOpenMC() const;
736 
742  void latticeOuterCheck(const Point & c, int level) const;
743 
749  void latticeOuterError(const Point & c, int level) const;
750 
756  bool findCell(const Point & point);
757 
766  void compareContainedCells(std::map<cellInfo, containedCells> & reference,
767  std::map<cellInfo, containedCells> & compare) const;
768 
769  NumericVector<Number> & _serialized_solution;
770 
775  virtual std::vector<int32_t> getMappedTallyIDs() const override;
776 
781  const bool & _output_cell_mapping;
782 
790 
793 
801 
807  unsigned int _cell_level;
808 
813  const bool & _export_properties;
814 
816  const bool _using_skinner;
817 
824 
857 
865 
874 
880 
886 
893 
895  bool _has_cell_tallies = false;
896 
899 
904  std::map<std::string, std::shared_ptr<FilterBase>> _filters;
905 
907  std::vector<std::shared_ptr<TallyBase>> _local_tallies;
908 
910  std::vector<std::string> _all_tally_scores;
911 
913  std::map<std::string, unsigned int> _score_count;
914 
916  std::vector<std::vector<unsigned int>> _tally_var_ids;
917 
919  std::vector<std::vector<std::vector<unsigned int>>> _tally_ext_var_ids;
920 
922  std::vector<SubdomainID> _density_blocks;
923 
925  std::vector<SubdomainID> _temp_blocks;
926 
928  std::unordered_set<SubdomainID> _identical_cell_fill_blocks;
929 
931  std::vector<cellInfo> _elem_to_cell{};
932 
934  std::map<cellInfo, coupling::CouplingFields> _cell_phase;
935 
938 
941 
944 
947 
953 
959 
965 
968 
971 
974 
976  std::map<cellInfo, std::vector<unsigned int>> _cell_to_elem;
977 
979  std::map<cellInfo, std::vector<unsigned int>> _local_cell_to_elem;
980 
982  std::map<cellInfo, std::unordered_set<SubdomainID>> _cell_to_elem_subdomain;
983 
985  std::map<SubdomainID, std::set<int32_t>> _subdomain_to_material;
986 
991  std::map<cellInfo, Point> _cell_to_point;
992 
997  std::map<cellInfo, Real> _cell_to_elem_volume;
998 
1003  std::map<cellInfo, Real> _cell_volume;
1004 
1009  std::map<cellInfo, int32_t> _cell_to_material;
1010 
1015  std::map<cellInfo, containedCells> _cell_to_contained_material_cells;
1016 
1018  std::map<cellInfo, int32_t> _cell_to_n_contained;
1019 
1021  static bool _first_transfer;
1022 
1024  static bool _printed_initial;
1025 
1028 
1030  openmc::Particle _particle;
1031 
1033  unsigned int _n_particles_1;
1034 
1036  std::map<std::string, std::vector<SubdomainName>> _temp_vars_to_blocks;
1037 
1039  std::map<std::string, std::vector<SubdomainName>> _density_vars_to_blocks;
1040 
1043 
1046 
1048  std::map<cellInfo, int> _n_temp;
1049 
1051  std::map<cellInfo, int> _n_rho;
1052 
1054  std::map<cellInfo, int> _n_temp_rho;
1055 
1057  std::map<cellInfo, int> _n_none;
1058 
1060  std::shared_ptr<TallyBase> _source_rate_norm_tally;
1061 
1064 
1065 #ifdef ENABLE_DAGMC
1066  MoabSkinner * _skinner = nullptr;
1068 
1070  std::shared_ptr<moab::DagMC> _dagmc = nullptr;
1071 #endif
1072 
1074  long unsigned int _n_openmc_cells;
1075 
1078 
1081 
1084 
1088 
1090  static constexpr Real EV_TO_JOULE = 1.6022e-19;
1091 
1093  static constexpr Real ZERO_TALLY_THRESHOLD = 1e-12;
1094 
1095 private:
1099  void dufekGudowskiParticleUpdate();
1100 
1102  std::vector<int32_t> _flattened_ids;
1103 
1105  std::vector<int32_t> _flattened_instances;
1106 
1108  containedCells _instance_offsets;
1109 
1111  std::map<cellInfo, int32_t> _n_offset;
1112 
1114  cellInfo _first_identical_cell;
1115 
1117  std::vector<int32_t> _first_identical_cell_materials;
1118 
1120  bool _use_displaced;
1121 
1123  std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_temp_vars;
1124 
1126  std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_density_vars;
1127 };
std::vector< const VariableValue * > getTallyScoreNeighborVariableValues(const std::string &score, const std::string &tally_name, THREAD_ID tid, const std::string &output="", bool skip_func_exp=false)
TallyTriggerTypeEnum
Type of trigger to apply.
Definition: CardinalEnums.h:190
int _source_rate_score
The score index into "_source_rate_norm_tally".
Definition: OpenMCCellAverageProblem.h:1063
const bool _has_identical_cell_fills
Definition: OpenMCCellAverageProblem.h:856
const trigger::TallyTriggerTypeEnum _k_trigger
Definition: OpenMCCellAverageProblem.h:800
int32_t elemToCellIndex(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:286
bool findCell(const Point &point)
int _n_moose_none_elems
Number of no-coupling elements in the MOOSE mesh.
Definition: OpenMCCellAverageProblem.h:946
void sendTemperatureToOpenMC() const
std::map< std::string, std::vector< SubdomainName > > _density_vars_to_blocks
Mapping from density variable name to the subdomains on which to read it from.
Definition: OpenMCCellAverageProblem.h:1039
static bool _first_transfer
Whether the present transfer is the first transfer.
Definition: OpenMCCellAverageProblem.h:1021
Definition: CardinalEnums.h:234
int _fixed_point_iteration
Definition: OpenMCProblemBase.h:499
const bool _specified_temperature_feedback
Definition: OpenMCCellAverageProblem.h:892
void checkEmptyVector(const std::vector< T > &vector, const std::string &name) const
Definition: OpenMCCellAverageProblem.h:433
bool hasScore(const std::string &score)
Definition: OpenMCCellAverageProblem.h:166
std::unordered_map< int32_t, std::vector< int32_t > > containedCells
Definition: OpenMCCellAverageProblem.h:279
std::map< std::string, std::vector< SubdomainName > > _temp_vars_to_blocks
Mapping from temperature variable name to the subdomains on which to read it from.
Definition: OpenMCCellAverageProblem.h:1036
std::vector< std::vector< std::vector< unsigned int > > > _tally_ext_var_ids
A vector of external (output-based) auxvariable ids added by the [Tallies] block.
Definition: OpenMCCellAverageProblem.h:919
void storeElementPhase()
Loop over the elements in the MOOSE mesh and store the type of feedback applied by each.
std::map< cellInfo, std::vector< unsigned int > > _local_cell_to_elem
Mapping of OpenMC cell indices to a vector of MOOSE element IDs, on each local rank.
Definition: OpenMCCellAverageProblem.h:979
const relaxation::RelaxationEnum _relaxation
Type of relaxation to apply to the OpenMC tallies.
Definition: OpenMCCellAverageProblem.h:792
std::shared_ptr< TallyBase > _source_rate_norm_tally
The tally to be used for normalizing all other tallies when running an eigenvalue calculation.
Definition: OpenMCCellAverageProblem.h:1060
std::vector< std::vector< unsigned int > > _tally_var_ids
A vector of auxvariable ids added by the [Tallies] block.
Definition: OpenMCCellAverageProblem.h:916
void checkContainedCellsStructure(const cellInfo &cell_info, containedCells &reference, containedCells &compare) const
int _n_moose_density_elems
Number of elements in the MOOSE mesh that exclusively provide density feedback.
Definition: OpenMCCellAverageProblem.h:937
static InputParameters validParams()
int _n_mapped_temp_density_elems
Definition: OpenMCCellAverageProblem.h:964
Real _uncoupled_volume
Total volume of uncoupled MOOSE mesh elements.
Definition: OpenMCCellAverageProblem.h:970
std::map< cellInfo, int > _n_rho
Number of density-only feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1051
NumericVector< Number > & _serialized_solution
Definition: OpenMCCellAverageProblem.h:769
OpenMCCellAverageProblem(const InputParameters &params)
int _n_moose_temp_density_elems
Number of elements in the MOOSE mesh which provide temperature+density feedback.
Definition: OpenMCCellAverageProblem.h:943
void initializeTallies()
Add OpenMC tallies to facilitate the coupling.
void readBlockVariables(const std::string &param, const std::string &default_name, std::map< std::string, std::vector< SubdomainName >> &vars_to_specified_blocks, std::vector< SubdomainID > &specified_blocks)
std::string printNewline() const
For keeping the output neat when using verbose.
Definition: OpenMCCellAverageProblem.h:654
std::map< cellInfo, int > _n_none
Number of none elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1057
containedCells containedMaterialCells(const cellInfo &cell_info) const
std::vector< int32_t > materialsInCells(const containedCells &contained_cells) const
void printTrisoHelp(const std::chrono::time_point< std::chrono::high_resolution_clock > &start) const
virtual void externalSolve() override
double cellTemperature(const cellInfo &cell_info) const
virtual void syncSolutions(ExternalProblem::Direction direction) override
Set the 'mesh changed' adaptivity flag.
std::vector< SubdomainID > _temp_blocks
Blocks in MOOSE mesh that provide temperature feedback.
Definition: OpenMCCellAverageProblem.h:925
bool hasTemperatureFeedback(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:350
virtual std::unordered_set< SubdomainID > getCellToElementSub(const cellInfo &info)
Definition: OpenMCCellAverageProblem.h:134
bool hasFilter(const std::string &filter_name) const
Definition: OpenMCCellAverageProblem.h:362
bool hasOutput(const std::string &score, const std::string &output) const
std::vector< const MooseVariableFE< Real > * > getTallyScoreVariables(const std::string &score, const std::string &tally_name, THREAD_ID tid, const std::string &output="", bool skip_func_exp=false)
bool _needs_to_map_cells
Whether any spatial mapping from OpenMC's cells to the mesh is needed.
Definition: OpenMCCellAverageProblem.h:898
Definition: CardinalEnums.h:232
Definition: OpenMCCellAverageProblem.h:66
int _n_mapped_density_elems
Definition: OpenMCCellAverageProblem.h:958
bool hasDensityFeedback(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:338
bool _map_density_by_cell
Definition: OpenMCCellAverageProblem.h:879
const bool _has_adaptivity
Whether or not the problem contains mesh adaptivity.
Definition: OpenMCProblemBase.h:534
std::map< cellInfo, int > _n_temp
Number of temperature-only feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1048
void sendDensityToOpenMC() const
virtual const std::vector< std::string > & getTallyScores() const
Definition: OpenMCCellAverageProblem.h:149
const bool & useDisplaced() const
Definition: OpenMCCellAverageProblem.h:463
void initializeElementToCellMapping()
Set up the mapping from MOOSE elements to OpenMC cells.
int _n_mapped_none_elems
Number of no-coupling elements mapped to OpenMC cells.
Definition: OpenMCCellAverageProblem.h:967
int _n_mapped_temp_elems
Definition: OpenMCCellAverageProblem.h:952
unsigned int getNumScoringTallies(const std::string &score) const
Definition: OpenMCCellAverageProblem.h:156
int32_t cellToMaterialIndex(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:318
void checkCellMappedPhase()
This function is used to ensure that each OpenMC cell only maps to a single phase.
bool _material_cells_only
Whether non-material cells are mapped.
Definition: OpenMCCellAverageProblem.h:973
Point transformPointToOpenMC(const Point &pt) const
unsigned int _cell_level
Definition: OpenMCCellAverageProblem.h:807
containedCells shiftCellInstances(const cellInfo &cell_info) const
int32_t elemToCellID(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:294
std::map< cellInfo, containedCells > _cell_to_contained_material_cells
Definition: OpenMCCellAverageProblem.h:1015
long unsigned int _n_openmc_cells
Total number of unique OpenMC cell IDs + instances combinations.
Definition: OpenMCCellAverageProblem.h:1074
void gatherCellSum(std::vector< T > &local, std::map< cellInfo, T > &global) const
bool hasAdaptivity() const
Definition: OpenMCCellAverageProblem.h:445
std::map< cellInfo, int32_t > _cell_to_n_contained
Number of material-type cells contained within a cell.
Definition: OpenMCCellAverageProblem.h:1018
Real tallyMultiplier(const std::string &score_name, const Real &local_mean_tally) const
OpenMCInitialCondition
Where to get the initial temperature and density settings for OpenMC.
Definition: CardinalEnums.h:239
void subdomainsToMaterials()
Loop over the mapped cells, and build a map between subdomains to OpenMC materials.
const bool & _assume_separate_tallies
Definition: OpenMCCellAverageProblem.h:873
const bool & _check_identical_cell_fills
Definition: OpenMCCellAverageProblem.h:864
void readBlockParameters(const std::string name, std::unordered_set< SubdomainID > &blocks)
int fixedPointIteration() const
Definition: OpenMCCellAverageProblem.h:439
const SymmetryPointGenerator * _symmetry
Userobject that maps from a partial-symmetry OpenMC model to a whole-domain [Mesh].
Definition: OpenMCCellAverageProblem.h:1045
const bool & _output_cell_mapping
Definition: OpenMCCellAverageProblem.h:781
Definition: CardinalEnums.h:233
std::set< SubdomainID > coupledSubdomains() const
std::vector< std::string > getMaterialInEachSubdomain() const
std::vector< std::shared_ptr< TallyBase > > _local_tallies
A vector of the tally objects created by the [Problem/Tallies] block.
Definition: OpenMCCellAverageProblem.h:907
std::map< cellInfo, std::vector< unsigned int > > _cell_to_elem
Mapping of OpenMC cell indices to a vector of MOOSE element IDs.
Definition: OpenMCCellAverageProblem.h:976
virtual std::vector< int32_t > getMappedTallyIDs() const override
double cellMappedVolume(const cellInfo &cell_info) const
coupling::CouplingFields cellFeedback(const cellInfo &cell_info) const
std::vector< SubdomainID > _density_blocks
Blocks in MOOSE mesh that provide density feedback.
Definition: OpenMCCellAverageProblem.h:922
int32_t cellID(const int32_t index) const
static constexpr int DIMENSION
Spatial dimension of the Monte Carlo problem.
Definition: OpenMCCellAverageProblem.h:451
RelaxationEnum
Type of relaxation.
Definition: CardinalEnums.h:273
std::map< cellInfo, int32_t > _cell_to_material
Definition: OpenMCCellAverageProblem.h:1009
std::vector< cellInfo > _elem_to_cell
Mapping of MOOSE elements to the OpenMC cell they map to (if any)
Definition: OpenMCCellAverageProblem.h:931
const std::vector< std::shared_ptr< TallyBase > > & getLocalTallies() const
Definition: OpenMCCellAverageProblem.h:378
int32_t _cell_using_dagmc_universe_id
ID of the OpenMC cell corresponding to the cell which uses the DAGMC universe as a fill.
Definition: OpenMCCellAverageProblem.h:1083
openmc::Particle _particle
Dummy particle to reduce number of allocations of particles for cell lookup routines.
Definition: OpenMCCellAverageProblem.h:1030
Definition: OpenMCProblemBase.h:55
bool cellHasIdenticalFill(const cellInfo &cell_info) const
Definition: TallyBase.h:35
virtual const OpenMCVolumeCalculation * volumeCalculation() const
Definition: OpenMCCellAverageProblem.h:118
void setContainedCells(const cellInfo &cell_info, const Point &hint, std::map< cellInfo, containedCells > &map)
bool cellMapsToSubdomain(const cellInfo &cell_info, const std::unordered_set< SubdomainID > &id) const
Definition: CardinalEnums.h:72
static bool _printed_initial
Whether the diagnostic tables on initialization have already been printed.
Definition: OpenMCCellAverageProblem.h:1024
std::map< std::string, unsigned int > _score_count
Number of tallies scoring a particular score.
Definition: OpenMCCellAverageProblem.h:913
OpenMCVolumeCalculation * _volume_calc
Optional volume calculation for cells which map to MOOSE.
Definition: OpenMCCellAverageProblem.h:1042
void getTallyTriggerParameters(const InputParameters &params)
bool _has_cell_tallies
Whether any cell tallies exist.
Definition: OpenMCCellAverageProblem.h:895
std::map< cellInfo, Real > _cell_volume
Definition: OpenMCCellAverageProblem.h:1003
Definition: OpenMCVolumeCalculation.h:30
const bool _specified_density_feedback
Definition: OpenMCCellAverageProblem.h:885
virtual bool converged(unsigned int) override
Definition: OpenMCCellAverageProblem.h:75
std::map< cellInfo, int > _n_temp_rho
Number of temperature+density feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1054
void setupProblem()
Initialize the mapping of OpenMC to the MooseMesh and perform additional setup actions.
bool _need_to_reinit_coupling
Definition: OpenMCCellAverageProblem.h:823
virtual bool hasPointTransformations() const
Definition: OpenMCCellAverageProblem.h:143
std::map< cellInfo, Point > _cell_to_point
Definition: OpenMCCellAverageProblem.h:991
static constexpr Real EV_TO_JOULE
Conversion rate from eV to Joule.
Definition: OpenMCCellAverageProblem.h:1090
virtual void createQRules(QuadratureType type, Order order, Order volume_order, Order face_order, SubdomainID block, bool allow_negative_weights=true) override
const TallyBase * getTally(const std::string &name)
std::shared_ptr< FilterBase > & getFilter(const std::string &filter_name)
Definition: OpenMCCellAverageProblem.h:369
virtual Real cellVolume(const cellInfo &cell_info) const
void reloadDAGMC()
Reconstruct the DAGMC geometry after skinning.
CouplingFields
Type of feedback in Monte Carlo simulation.
Definition: CardinalEnums.h:230
const bool _using_skinner
Whether or not the problem uses a skinner to regenerate the OpenMC geometry.
Definition: OpenMCCellAverageProblem.h:816
unsigned int _n_particles_1
Number of particles simulated in the first iteration in Dufek-Gudowski relaxation.
Definition: OpenMCCellAverageProblem.h:1033
static bool _printed_triso_warning
Whether a warning has already been printed about very long setup times (for TRISOs)
Definition: OpenMCCellAverageProblem.h:1027
void compareContainedCells(std::map< cellInfo, containedCells > &reference, std::map< cellInfo, containedCells > &compare) const
int32_t elemToCellInstance(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:301
std::unordered_set< SubdomainID > _identical_cell_fill_blocks
Blocks for which the cell fills are identical.
Definition: OpenMCCellAverageProblem.h:928
static constexpr Real ZERO_TALLY_THRESHOLD
Tolerance for setting zero tally.
Definition: OpenMCCellAverageProblem.h:1093
virtual const std::map< cellInfo, std::vector< unsigned int > > & cellToElem() const
Definition: OpenMCCellAverageProblem.h:124
unsigned int getCellLevel(const Point &c) const
std::map< std::string, std::shared_ptr< FilterBase > > _filters
Definition: OpenMCCellAverageProblem.h:904
cellInfo elemToCellInfo(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:309
const bool & _verbose
Whether to print diagnostic information about model setup and the transfers.
Definition: OpenMCProblemBase.h:449
std::vector< std::string > _all_tally_scores
A list of all of the scores contained by the local tallies added in the [Tallies] block.
Definition: OpenMCCellAverageProblem.h:910
void read2DBlockParameters(const std::string name, std::vector< std::vector< SubdomainName >> &names, std::vector< SubdomainID > &flattened_ids)
std::map< cellInfo, coupling::CouplingFields > _cell_phase
Phase of each cell.
Definition: OpenMCCellAverageProblem.h:934
Point transformPoint(const Point &p) const
Definition: SymmetryPointGenerator.C:100
Definition: SymmetryPointGenerator.h:28
virtual Point transformPoint(const Point &pt) const
Definition: OpenMCCellAverageProblem.h:242
void addFilter(const std::string &type, const std::string &name, InputParameters &moose_object_pars)
void checkBlocksInMesh(const std::string name, const std::vector< SubdomainID > &ids, const std::vector< SubdomainName > &names) const
void latticeOuterCheck(const Point &c, int level) const
void getCellMappedSubdomains()
Loop over all the OpenMC cells and get the element subdomain IDs that map to each cell.
Skins the [Mesh] according to individual bins for temperature, density, and subdomain ID.
Definition: MoabSkinner.h:19
const int32_t _initial_num_openmc_surfaces
Definition: OpenMCCellAverageProblem.h:1087
const bool & _export_properties
Definition: OpenMCCellAverageProblem.h:813
virtual MooseMesh & getMooseMesh()
Get a modifyable non-const reference to the Moose mesh.
std::map< SubdomainID, std::set< int32_t > > _subdomain_to_material
Mapping of elem subdomains to materials.
Definition: OpenMCCellAverageProblem.h:985
Definition: CardinalEnums.h:130
std::map< cellInfo, std::unordered_set< SubdomainID > > _cell_to_elem_subdomain
Mapping of OpenMC cell indices to the subdomain IDs each maps to.
Definition: OpenMCCellAverageProblem.h:982
int32_t _dagmc_universe_id
ID of the OpenMC universe corresponding to the DAGMC universe.
Definition: OpenMCCellAverageProblem.h:1077
virtual void addExternalVariables() override
const coupling::OpenMCInitialCondition _initial_condition
Definition: OpenMCCellAverageProblem.h:789
bool _dagmc_root_universe
Whether the DAGMC universe is the root universe or not.
Definition: OpenMCCellAverageProblem.h:1080
void mapElemsToCells()
Populate maps of MOOSE elements to OpenMC cells.
coupling::CouplingFields elemFeedback(const Elem *elem) const
std::shared_ptr< TallyBase > addTally(const std::string &type, const std::string &name, InputParameters &moose_object_pars)
void setMinimumVolumeQRules(Order &volume_order, const std::string &type)
void latticeOuterError(const Point &c, int level) const
std::map< cellInfo, Real > _cell_to_elem_volume
Definition: OpenMCCellAverageProblem.h:997
void getMaterialFills()
Find the material filling each cell which receives density feedback.
std::pair< int32_t, int32_t > cellInfo
Definition: OpenMCProblemBase.h:201
void gatherCellVector(std::vector< T > &local, std::vector< unsigned int > &n_local, std::map< cellInfo, std::vector< T >> &global)
std::vector< const VariableValue * > getTallyScoreVariableValues(const std::string &score, const std::string &tally_name, THREAD_ID tid, const std::string &output="", bool skip_func_exp=false)
int _n_moose_temp_elems
Number of elements in the MOOSE mesh that exclusively provide temperature feedback.
Definition: OpenMCCellAverageProblem.h:940
cellInfo firstContainedMaterialCell(const cellInfo &cell_info) const
static constexpr int32_t UNMAPPED
Constant flag to indicate that a cell/element was unmapped.
Definition: OpenMCCellAverageProblem.h:448
std::map< cellInfo, Real > computeVolumeWeightedCellInput(const std::map< SubdomainID, std::pair< unsigned int, std::string >> &var_num, const std::vector< coupling::CouplingFields > *phase) const
virtual void initialSetup() override
Set up an OpenMC simulation.