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  bool hasScore(const std::string & score)
157  {
158  return std::find(_all_tally_scores.begin(), _all_tally_scores.end(), score) !=
159  _all_tally_scores.end();
160  }
161 
167  std::vector<const TallyBase *> getTalliesByScore(const std::string & score);
168 
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 
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 
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 
219  bool hasOutput(const std::string & score, const std::string & output) const;
220 
226  virtual Point transformPoint(const Point & pt) const
227  {
228  return this->hasPointTransformations() ? _symmetry->transformPoint(pt) : pt;
229  }
230 
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 
263  typedef std::unordered_map<int32_t, std::vector<int32_t>> containedCells;
264 
270  int32_t elemToCellIndex(const int & elem_id) const { return elemToCellInfo(elem_id).first; }
271 
278  int32_t elemToCellID(const int & elem_id) const { return cellID(elemToCellIndex(elem_id)); }
279 
285  int32_t elemToCellInstance(const int & elem_id) const { return elemToCellInfo(elem_id).second; }
286 
293  cellInfo elemToCellInfo(const int & elem_id) const { return _elem_to_cell[elem_id]; }
294 
302  int32_t cellToMaterialIndex(const cellInfo & cell_info) const
303  {
304  return _cell_to_material.at(cell_info);
305  }
306 
315  coupling::CouplingFields cellFeedback(const cellInfo & cell_info) const;
316 
322  bool hasDensityFeedback(const cellInfo & cell_info) const
323  {
324  std::vector<coupling::CouplingFields> phase = {coupling::density,
326  return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
327  }
328 
334  bool hasTemperatureFeedback(const cellInfo & cell_info) const
335  {
336  std::vector<coupling::CouplingFields> phase = {coupling::temperature,
338  return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
339  }
340 
346  bool hasFilter(const std::string & filter_name) const { return _filters.count(filter_name) > 0; }
347 
353  std::shared_ptr<FilterBase> & getFilter(const std::string & filter_name)
354  {
355  return _filters.at(filter_name);
356  }
357 
362  const std::vector<std::shared_ptr<TallyBase>> & getLocalTally() const { return _local_tallies; }
363 
369  double cellTemperature(const cellInfo & cell_info) const;
370 
375  double cellMappedVolume(const cellInfo & cell_info) const;
376 
378  void reloadDAGMC();
379 
386  void addFilter(const std::string & type,
387  const std::string & name,
388  InputParameters & moose_object_pars);
389 
396  void
397  addTally(const std::string & type, const std::string & name, InputParameters & moose_object_pars);
398 
407  Real tallyMultiplier(unsigned int global_score) const;
408 
414  template <typename T>
415  void checkEmptyVector(const std::vector<T> & vector, const std::string & name) const
416  {
417  if (vector.empty())
418  mooseError(name + " cannot be empty!");
419  }
420 
422 
427  bool hasAdaptivity() const { return _has_adaptivity; }
428 
430  static constexpr int32_t UNMAPPED{-1};
431 
433  static constexpr int DIMENSION{3};
434 
436  virtual MooseMesh & getMooseMesh();
437 
439  virtual const MooseMesh & getMooseMesh() const;
440 
445  const bool & useDisplaced() const { return _use_displaced; }
446 
447 protected:
453  unsigned int getCellLevel(const Point & c) const;
454 
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 
474  bool cellHasIdenticalFill(const cellInfo & cell_info) const;
475 
482  containedCells shiftCellInstances(const cellInfo & cell_info) const;
483 
490  bool cellMapsToSubdomain(const cellInfo & cell_info,
491  const std::unordered_set<SubdomainID> & id) const;
492 
498  cellInfo firstContainedMaterialCell(const cellInfo & cell_info) const;
499 
505  containedCells containedMaterialCells(const cellInfo & cell_info) const;
506 
510  void updateOpenMCGeometry();
511 
516  virtual void updateMaterials();
517 
523  std::vector<std::string> getMaterialInEachSubdomain() const;
524 
530  Point transformPointToOpenMC(const Point & pt) const;
531 
538  void checkNormalization(const Real & sum, unsigned int global_score) const;
539 
549  void
550  printTrisoHelp(const std::chrono::time_point<std::chrono::high_resolution_clock> & start) const;
551 
558  void printAuxVariableIO();
559 
565  std::vector<int32_t> materialsInCells(const containedCells & contained_cells) const;
566 
568  void subdomainsToMaterials();
569 
574  std::set<SubdomainID> coupledSubdomains() const;
575 
581  template <typename T>
582  void gatherCellSum(std::vector<T> & local, std::map<cellInfo, T> & global) const;
583 
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 
598  coupling::CouplingFields elemFeedback(const Elem * elem) const;
599 
604  void getTallyTriggerParameters(const InputParameters & params);
605 
611  void readBlockParameters(const std::string name, std::unordered_set<SubdomainID> & blocks);
612 
618  void cacheContainedCells();
619 
626  void setContainedCells(const cellInfo & cell_info,
627  const Point & hint,
628  std::map<cellInfo, containedCells> & map);
629 
638  void checkContainedCellsStructure(const cellInfo & cell_info,
639  containedCells & reference,
640  containedCells & compare) const;
641 
647  void setMinimumVolumeQRules(Order & volume_order, const std::string & type);
648 
650  std::string printNewline() const
651  {
652  if (_verbose)
653  return "\n";
654  else
655  return "";
656  }
657 
659  void storeElementPhase();
660 
682  void relaxAndNormalizeTally(unsigned int global_score,
683  unsigned int local_score,
684  std::shared_ptr<TallyBase> local_tally);
685 
690  void getCellMappedPhase();
691 
693  void checkCellMappedPhase();
694 
697 
703 
706 
708  void mapElemsToCells();
709 
715  void validateLocalTallies();
716 
718  void initializeTallies();
719 
724  void resetTallies();
725 
727  void getMaterialFills();
728 
734  void getPointInCell();
735 
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 
750  void sendTemperatureToOpenMC() const;
751 
756  void sendDensityToOpenMC() const;
757 
763  Real tallyNormalization(unsigned int global_score) const;
764 
769  void checkTallySum(const unsigned int & score) const;
770 
776  void latticeOuterCheck(const Point & c, int level) const;
777 
783  void latticeOuterError(const Point & c, int level) const;
784 
790  bool findCell(const Point & point);
791 
800  void compareContainedCells(std::map<cellInfo, containedCells> & reference,
801  std::map<cellInfo, containedCells> & compare) const;
802 
803  NumericVector<Number> & _serialized_solution;
804 
809  virtual std::vector<int32_t> getMappedTallyIDs() const override;
810 
815  const bool & _output_cell_mapping;
816 
824 
827 
835 
841  unsigned int _cell_level;
842 
847  const bool & _export_properties;
848 
869 
871  const bool _using_skinner;
872 
879 
893  const bool _check_tally_sum;
894 
896  const Real & _relaxation_factor;
897 
930 
938 
947 
953 
959 
966 
968  bool _has_cell_tallies = false;
969 
972 
978 
983  std::map<std::string, std::shared_ptr<FilterBase>> _filters;
984 
986  std::vector<std::shared_ptr<TallyBase>> _local_tallies;
987 
989  std::vector<std::string> _all_tally_scores;
990 
996  std::vector<std::map<std::string, int>> _local_tally_score_map;
997 
999  std::vector<std::vector<unsigned int>> _tally_var_ids;
1000 
1002  std::vector<std::vector<std::vector<unsigned int>>> _tally_ext_var_ids;
1003 
1005  std::vector<SubdomainID> _density_blocks;
1006 
1008  std::vector<SubdomainID> _temp_blocks;
1009 
1011  std::unordered_set<SubdomainID> _identical_cell_fill_blocks;
1012 
1014  std::vector<cellInfo> _elem_to_cell{};
1015 
1017  std::map<cellInfo, coupling::CouplingFields> _cell_phase;
1018 
1021 
1024 
1027 
1030 
1036 
1042 
1048 
1051 
1054 
1057 
1059  std::map<cellInfo, std::vector<unsigned int>> _cell_to_elem;
1060 
1062  std::map<cellInfo, std::vector<unsigned int>> _local_cell_to_elem;
1063 
1065  std::map<cellInfo, std::unordered_set<SubdomainID>> _cell_to_elem_subdomain;
1066 
1068  std::map<SubdomainID, std::set<int32_t>> _subdomain_to_material;
1069 
1074  std::map<cellInfo, Point> _cell_to_point;
1075 
1080  std::map<cellInfo, Real> _cell_to_elem_volume;
1081 
1086  std::map<cellInfo, Real> _cell_volume;
1087 
1092  std::map<cellInfo, int32_t> _cell_to_material;
1093 
1098  std::map<cellInfo, containedCells> _cell_to_contained_material_cells;
1099 
1101  std::map<cellInfo, int32_t> _cell_to_n_contained;
1102 
1107  std::vector<openmc::Tally *> _global_tallies;
1108 
1110  std::vector<std::vector<std::string>> _global_tally_scores;
1111 
1113  std::vector<openmc::TallyEstimator> _global_tally_estimators;
1114 
1116  std::vector<Real> _global_sum_tally;
1117 
1119  std::vector<Real> _local_sum_tally;
1120 
1122  std::vector<Real> _local_mean_tally;
1123 
1125  static bool _first_transfer;
1126 
1128  static bool _printed_initial;
1129 
1132 
1134  openmc::Particle _particle;
1135 
1137  unsigned int _n_particles_1;
1138 
1140  std::map<std::string, std::vector<SubdomainName>> _temp_vars_to_blocks;
1141 
1143  std::map<std::string, std::vector<SubdomainName>> _density_vars_to_blocks;
1144 
1147 
1150 
1152  std::map<cellInfo, int> _n_temp;
1153 
1155  std::map<cellInfo, int> _n_rho;
1156 
1158  std::map<cellInfo, int> _n_temp_rho;
1159 
1161  std::map<cellInfo, int> _n_none;
1162 
1164  unsigned int _global_tally_index = 0;
1165 
1167  unsigned int _source_rate_index;
1168 
1169 #ifdef ENABLE_DAGMC
1170  MoabSkinner * _skinner = nullptr;
1172 
1174  std::shared_ptr<moab::DagMC> _dagmc = nullptr;
1175 #endif
1176 
1178  long unsigned int _n_openmc_cells;
1179 
1182 
1185 
1188 
1192 
1194  static constexpr Real EV_TO_JOULE = 1.6022e-19;
1195 
1197  static constexpr Real ZERO_TALLY_THRESHOLD = 1e-12;
1198 
1199 private:
1203  void dufekGudowskiParticleUpdate();
1204 
1206  std::vector<int32_t> _flattened_ids;
1207 
1209  std::vector<int32_t> _flattened_instances;
1210 
1212  containedCells _instance_offsets;
1213 
1215  std::map<cellInfo, int32_t> _n_offset;
1216 
1218  cellInfo _first_identical_cell;
1219 
1221  std::vector<int32_t> _first_identical_cell_materials;
1222 
1224  bool _use_displaced;
1225 
1227  std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_temp_vars;
1228 
1230  std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_density_vars;
1231 };
TallyTriggerTypeEnum
Type of trigger to apply.
Definition: CardinalEnums.h:189
void checkNormalization(const Real &sum, unsigned int global_score) const
std::vector< const MooseVariableFE< Real > * > getTallyScoreVariables(const std::string &score, THREAD_ID tid, const std::string &output="", bool skip_func_exp=false)
const bool _has_identical_cell_fills
Definition: OpenMCCellAverageProblem.h:929
const trigger::TallyTriggerTypeEnum _k_trigger
Definition: OpenMCCellAverageProblem.h:834
int32_t elemToCellIndex(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:270
bool findCell(const Point &point)
int _n_moose_none_elems
Number of no-coupling elements in the MOOSE mesh.
Definition: OpenMCCellAverageProblem.h:1029
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:1143
static bool _first_transfer
Whether the present transfer is the first transfer.
Definition: OpenMCCellAverageProblem.h:1125
Definition: CardinalEnums.h:233
int _fixed_point_iteration
Definition: OpenMCProblemBase.h:466
const bool _specified_temperature_feedback
Definition: OpenMCCellAverageProblem.h:965
void checkEmptyVector(const std::vector< T > &vector, const std::string &name) const
Definition: OpenMCCellAverageProblem.h:415
bool hasScore(const std::string &score)
Definition: OpenMCCellAverageProblem.h:156
std::unordered_map< int32_t, std::vector< int32_t > > containedCells
Definition: OpenMCCellAverageProblem.h:263
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:1140
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:1002
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:1062
const relaxation::RelaxationEnum _relaxation
Type of relaxation to apply to the OpenMC tallies.
Definition: OpenMCCellAverageProblem.h:826
void relaxAndNormalizeTally(unsigned int global_score, unsigned int local_score, std::shared_ptr< TallyBase > local_tally)
std::vector< std::vector< unsigned int > > _tally_var_ids
A vector of auxvariable ids added by the [Tallies] block.
Definition: OpenMCCellAverageProblem.h:999
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:1020
static InputParameters validParams()
int _n_mapped_temp_density_elems
Definition: OpenMCCellAverageProblem.h:1047
Real _uncoupled_volume
Total volume of uncoupled MOOSE mesh elements.
Definition: OpenMCCellAverageProblem.h:1053
std::map< cellInfo, int > _n_rho
Number of density-only feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1155
NumericVector< Number > & _serialized_solution
Definition: OpenMCCellAverageProblem.h:803
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:1026
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:650
std::map< cellInfo, int > _n_none
Number of none elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1161
const bool _needs_global_tally
Definition: OpenMCCellAverageProblem.h:977
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
std::vector< Real > _local_mean_tally
Mean value of the local tally(s), across all bins; only used for fixed source mode.
Definition: OpenMCCellAverageProblem.h:1122
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:1008
bool hasTemperatureFeedback(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:334
virtual std::unordered_set< SubdomainID > getCellToElementSub(const cellInfo &info)
Definition: OpenMCCellAverageProblem.h:134
bool hasFilter(const std::string &filter_name) const
Definition: OpenMCCellAverageProblem.h:346
bool hasOutput(const std::string &score, const std::string &output) const
std::vector< Real > _local_sum_tally
Sum value of the local tally(s), across all bins.
Definition: OpenMCCellAverageProblem.h:1119
bool _needs_to_map_cells
Whether any spatial mapping from OpenMC's cells to the mesh is needed.
Definition: OpenMCCellAverageProblem.h:971
Definition: CardinalEnums.h:231
Definition: OpenMCCellAverageProblem.h:66
int _n_mapped_density_elems
Definition: OpenMCCellAverageProblem.h:1041
bool hasDensityFeedback(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:322
bool _map_density_by_cell
Definition: OpenMCCellAverageProblem.h:952
std::vector< Real > _global_sum_tally
Sum value of the global tally(s), across all bins.
Definition: OpenMCCellAverageProblem.h:1116
const bool _has_adaptivity
Whether or not the problem contains mesh adaptivity.
Definition: OpenMCProblemBase.h:498
std::map< cellInfo, int > _n_temp
Number of temperature-only feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1152
void sendDensityToOpenMC() const
virtual const std::vector< std::string > & getTallyScores() const
Definition: OpenMCCellAverageProblem.h:149
const bool & useDisplaced() const
Definition: OpenMCCellAverageProblem.h:445
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:1050
int _n_mapped_temp_elems
Definition: OpenMCCellAverageProblem.h:1035
int32_t cellToMaterialIndex(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:302
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:1056
Point transformPointToOpenMC(const Point &pt) const
unsigned int _cell_level
Definition: OpenMCCellAverageProblem.h:841
containedCells shiftCellInstances(const cellInfo &cell_info) const
int32_t elemToCellID(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:278
std::map< cellInfo, containedCells > _cell_to_contained_material_cells
Definition: OpenMCCellAverageProblem.h:1098
long unsigned int _n_openmc_cells
Total number of unique OpenMC cell IDs + instances combinations.
Definition: OpenMCCellAverageProblem.h:1178
void gatherCellSum(std::vector< T > &local, std::map< cellInfo, T > &global) const
bool hasAdaptivity() const
Definition: OpenMCCellAverageProblem.h:427
std::map< cellInfo, int32_t > _cell_to_n_contained
Number of material-type cells contained within a cell.
Definition: OpenMCCellAverageProblem.h:1101
std::vector< const TallyBase * > getTalliesByScore(const std::string &score)
OpenMCInitialCondition
Where to get the initial temperature and density settings for OpenMC.
Definition: CardinalEnums.h:238
void subdomainsToMaterials()
Loop over the mapped cells, and build a map between subdomains to OpenMC materials.
const bool & _assume_separate_tallies
Definition: OpenMCCellAverageProblem.h:946
std::vector< openmc::TallyEstimator > _global_tally_estimators
Global tally estimators corresponding to '_global_tallies'.
Definition: OpenMCCellAverageProblem.h:1113
const bool & _check_identical_cell_fills
Definition: OpenMCCellAverageProblem.h:937
void readBlockParameters(const std::string name, std::unordered_set< SubdomainID > &blocks)
int fixedPointIteration() const
Definition: OpenMCCellAverageProblem.h:421
const SymmetryPointGenerator * _symmetry
Userobject that maps from a partial-symmetry OpenMC model to a whole-domain [Mesh].
Definition: OpenMCCellAverageProblem.h:1149
unsigned int _source_rate_index
Index in tally_score pointing to the score used for normalizing flux tallies in eigenvalue mode.
Definition: OpenMCCellAverageProblem.h:1167
const bool & _output_cell_mapping
Definition: OpenMCCellAverageProblem.h:815
Definition: CardinalEnums.h:232
void checkTallySum(const unsigned int &score) const
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:986
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:1059
virtual std::vector< int32_t > getMappedTallyIDs() const override
double cellMappedVolume(const cellInfo &cell_info) const
std::vector< openmc::Tally * > _global_tallies
Definition: OpenMCCellAverageProblem.h:1107
coupling::CouplingFields cellFeedback(const cellInfo &cell_info) const
std::vector< SubdomainID > _density_blocks
Blocks in MOOSE mesh that provide density feedback.
Definition: OpenMCCellAverageProblem.h:1005
int32_t cellID(const int32_t index) const
static constexpr int DIMENSION
Spatial dimension of the Monte Carlo problem.
Definition: OpenMCCellAverageProblem.h:433
RelaxationEnum
Type of relaxation.
Definition: CardinalEnums.h:272
std::map< cellInfo, int32_t > _cell_to_material
Definition: OpenMCCellAverageProblem.h:1092
Real tallyMultiplier(unsigned int global_score) const
std::vector< cellInfo > _elem_to_cell
Mapping of MOOSE elements to the OpenMC cell they map to (if any)
Definition: OpenMCCellAverageProblem.h:1014
void addTally(const std::string &type, const std::string &name, InputParameters &moose_object_pars)
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:1187
openmc::Particle _particle
Dummy particle to reduce number of allocations of particles for cell lookup routines.
Definition: OpenMCCellAverageProblem.h:1134
const bool _check_tally_sum
Definition: OpenMCCellAverageProblem.h:893
Definition: OpenMCProblemBase.h:52
bool cellHasIdenticalFill(const cellInfo &cell_info) const
virtual const OpenMCVolumeCalculation * volumeCalculation() const
Definition: OpenMCCellAverageProblem.h:118
const bool _normalize_by_global
Definition: OpenMCCellAverageProblem.h:868
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:71
static bool _printed_initial
Whether the diagnostic tables on initialization have already been printed.
Definition: OpenMCCellAverageProblem.h:1128
OpenMCVolumeCalculation * _volume_calc
Optional volume calculation for cells which map to MOOSE.
Definition: OpenMCCellAverageProblem.h:1146
void getTallyTriggerParameters(const InputParameters &params)
bool _has_cell_tallies
Whether any cell tallies exist.
Definition: OpenMCCellAverageProblem.h:968
std::map< cellInfo, Real > _cell_volume
Definition: OpenMCCellAverageProblem.h:1086
std::vector< const VariableValue * > getTallyScoreNeighborVariableValues(const std::string &score, THREAD_ID tid, const std::string &output="", bool skip_func_exp=false)
Definition: OpenMCVolumeCalculation.h:30
const bool _specified_density_feedback
Definition: OpenMCCellAverageProblem.h:958
virtual bool converged(unsigned int) override
Definition: OpenMCCellAverageProblem.h:75
std::vector< std::map< std::string, int > > _local_tally_score_map
Definition: OpenMCCellAverageProblem.h:996
std::map< cellInfo, int > _n_temp_rho
Number of temperature+density feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1158
void setupProblem()
Initialize the mapping of OpenMC to the MooseMesh and perform additional setup actions.
bool _need_to_reinit_coupling
Definition: OpenMCCellAverageProblem.h:878
virtual bool hasPointTransformations() const
Definition: OpenMCCellAverageProblem.h:143
std::map< cellInfo, Point > _cell_to_point
Definition: OpenMCCellAverageProblem.h:1074
static constexpr Real EV_TO_JOULE
Conversion rate from eV to Joule.
Definition: OpenMCCellAverageProblem.h:1194
virtual void createQRules(QuadratureType type, Order order, Order volume_order, Order face_order, SubdomainID block, bool allow_negative_weights=true) override
std::shared_ptr< FilterBase > & getFilter(const std::string &filter_name)
Definition: OpenMCCellAverageProblem.h:353
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:229
const bool _using_skinner
Whether or not the problem uses a skinner to regenerate the OpenMC geometry.
Definition: OpenMCCellAverageProblem.h:871
unsigned int _n_particles_1
Number of particles simulated in the first iteration.
Definition: OpenMCCellAverageProblem.h:1137
static bool _printed_triso_warning
Whether a warning has already been printed about very long setup times (for TRISOs)
Definition: OpenMCCellAverageProblem.h:1131
void compareContainedCells(std::map< cellInfo, containedCells > &reference, std::map< cellInfo, containedCells > &compare) const
int32_t elemToCellInstance(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:285
unsigned int _global_tally_index
Index in OpenMC tallies corresponding to the first global tally added by Cardinal.
Definition: OpenMCCellAverageProblem.h:1164
std::unordered_set< SubdomainID > _identical_cell_fill_blocks
Blocks for which the cell fills are identical.
Definition: OpenMCCellAverageProblem.h:1011
static constexpr Real ZERO_TALLY_THRESHOLD
Tolerance for setting zero tally.
Definition: OpenMCCellAverageProblem.h:1197
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:983
cellInfo elemToCellInfo(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:293
const bool & _verbose
Whether to print diagnostic information about model setup and the transfers.
Definition: OpenMCProblemBase.h:416
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:989
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:1017
Point transformPoint(const Point &p) const
Definition: SymmetryPointGenerator.C:100
Definition: SymmetryPointGenerator.h:28
virtual Point transformPoint(const Point &pt) const
Definition: OpenMCCellAverageProblem.h:226
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:1191
const bool & _export_properties
Definition: OpenMCCellAverageProblem.h:847
const std::vector< std::shared_ptr< TallyBase > > & getLocalTally() const
Definition: OpenMCCellAverageProblem.h:362
virtual void updateMaterials()
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:1068
std::vector< const VariableValue * > getTallyScoreVariableValues(const std::string &score, THREAD_ID tid, const std::string &output="", bool skip_func_exp=false)
Real tallyNormalization(unsigned int global_score) const
Definition: CardinalEnums.h:129
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:1065
int32_t _dagmc_universe_id
ID of the OpenMC universe corresponding to the DAGMC universe.
Definition: OpenMCCellAverageProblem.h:1181
virtual void addExternalVariables() override
const coupling::OpenMCInitialCondition _initial_condition
Definition: OpenMCCellAverageProblem.h:823
bool _dagmc_root_universe
Whether the DAGMC universe is the root universe or not.
Definition: OpenMCCellAverageProblem.h:1184
void mapElemsToCells()
Populate maps of MOOSE elements to OpenMC cells.
coupling::CouplingFields elemFeedback(const Elem *elem) const
void setMinimumVolumeQRules(Order &volume_order, const std::string &type)
const Real & _relaxation_factor
Constant relaxation factor.
Definition: OpenMCCellAverageProblem.h:896
void latticeOuterError(const Point &c, int level) const
std::map< cellInfo, Real > _cell_to_elem_volume
Definition: OpenMCCellAverageProblem.h:1080
void getMaterialFills()
Find the material filling each cell which receives density feedback.
std::pair< int32_t, int32_t > cellInfo
Definition: OpenMCProblemBase.h:197
void gatherCellVector(std::vector< T > &local, std::vector< unsigned int > &n_local, std::map< cellInfo, std::vector< T >> &global)
int _n_moose_temp_elems
Number of elements in the MOOSE mesh that exclusively provide temperature feedback.
Definition: OpenMCCellAverageProblem.h:1023
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:430
std::vector< std::vector< std::string > > _global_tally_scores
Global tally scores corresponding to '_global_tallies'.
Definition: OpenMCCellAverageProblem.h:1110
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