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"
24 
26 #include "TallyBase.h"
27 
28 #include "openmc/tallies/filter_mesh.h"
29 
30 #ifdef ENABLE_DAGMC
31 #include "MoabSkinner.h"
32 #include "DagMC.hpp"
33 #endif
34 
66 {
67 public:
68  OpenMCCellAverageProblem(const InputParameters & params);
69  static InputParameters validParams();
70 
71  virtual void initialSetup() override;
72  virtual void externalSolve() override;
73  virtual void syncSolutions(ExternalProblem::Direction direction) override;
74  virtual bool converged(unsigned int) override { return true; }
75 
83  void read2DBlockParameters(const std::string name,
84  std::vector<std::vector<SubdomainName>> & names,
85  std::vector<SubdomainID> & flattened_ids);
86 
93  void checkBlocksInMesh(const std::string name,
94  const std::vector<SubdomainID> & ids,
95  const std::vector<SubdomainName> & names) const;
96 
98  void setupProblem();
99 
104  virtual void addExternalVariables() override;
105 
111  virtual Real cellVolume(const cellInfo & cell_info) const;
112 
117  virtual const OpenMCVolumeCalculation * volumeCalculation() const { return _volume_calc; }
118 
123  virtual const std::map<cellInfo, std::vector<unsigned int>> & cellToElem() const
124  {
125  return _cell_to_elem;
126  }
127 
133  virtual std::unordered_set<SubdomainID> getCellToElementSub(const cellInfo & info)
134  {
135  return _cell_to_elem_subdomain.at(info);
136  }
137 
142  virtual bool hasPointTransformations() const { return _symmetry != nullptr; }
143 
148  virtual const std::vector<std::string> & getTallyScores() const { return _all_tally_scores; }
149 
155  virtual Point transformPoint(const Point & pt) const
156  {
157  return this->hasPointTransformations() ? _symmetry->transformPoint(pt) : pt;
158  }
159 
181  virtual void createQRules(QuadratureType type,
182  Order order,
183  Order volume_order,
184  Order face_order,
185  SubdomainID block,
186  bool allow_negative_weights = true) override;
187 
192  typedef std::unordered_map<int32_t, std::vector<int32_t>> containedCells;
193 
199  int32_t elemToCellIndex(const int & elem_id) const { return elemToCellInfo(elem_id).first; }
200 
207  int32_t elemToCellID(const int & elem_id) const { return cellID(elemToCellIndex(elem_id)); }
208 
214  int32_t elemToCellInstance(const int & elem_id) const { return elemToCellInfo(elem_id).second; }
215 
222  cellInfo elemToCellInfo(const int & elem_id) const { return _elem_to_cell[elem_id]; }
223 
231  int32_t cellToMaterialIndex(const cellInfo & cell_info) { return _cell_to_material[cell_info]; }
232 
241  coupling::CouplingFields cellFeedback(const cellInfo & cell_info) const;
242 
248  bool hasDensityFeedback(const cellInfo & cell_info) const
249  {
250  std::vector<coupling::CouplingFields> phase = {coupling::density,
252  return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
253  }
254 
260  bool hasTemperatureFeedback(const cellInfo & cell_info) const
261  {
262  std::vector<coupling::CouplingFields> phase = {coupling::temperature,
264  return std::find(phase.begin(), phase.end(), cellFeedback(cell_info)) != phase.end();
265  }
266 
271  const std::vector<std::shared_ptr<TallyBase>> & getLocalTally() const { return _local_tallies; }
272 
278  double cellTemperature(const cellInfo & cell_info) const;
279 
284  double cellMappedVolume(const cellInfo & cell_info) const;
285 
287  void reloadDAGMC();
288 
296  void addTallyObject(const std::string & type,
297  const std::string & name,
298  InputParameters & moose_object_pars);
299 
308  Real tallyMultiplier(unsigned int global_score) const;
309 
311 
313  static constexpr int32_t UNMAPPED{-1};
314 
316  static constexpr int DIMENSION{3};
317 
318 protected:
324  unsigned int getCellLevel(const Point & c) const;
325 
334  void
335  readBlockVariables(const std::string & param,
336  const std::string & default_name,
337  std::map<std::string, std::vector<SubdomainName>> & vars_to_specified_blocks,
338  std::vector<SubdomainID> & specified_blocks);
339 
345  bool cellHasIdenticalFill(const cellInfo & cell_info) const;
346 
353  containedCells shiftCellInstances(const cellInfo & cell_info) const;
354 
361  bool cellMapsToSubdomain(const cellInfo & cell_info,
362  const std::unordered_set<SubdomainID> & id) const;
363 
369  cellInfo firstContainedMaterialCell(const cellInfo & cell_info) const;
370 
376  containedCells containedMaterialCells(const cellInfo & cell_info) const;
377 
382  virtual void updateMaterials();
383 
389  std::vector<std::string> getMaterialInEachSubdomain() const;
390 
396  Point transformPointToOpenMC(const Point & pt) const;
397 
404  void checkNormalization(const Real & sum, unsigned int global_score) const;
405 
415  void
416  printTrisoHelp(const std::chrono::time_point<std::chrono::high_resolution_clock> & start) const;
417 
424  void printAuxVariableIO();
425 
431  std::vector<int32_t> materialsInCells(const containedCells & contained_cells) const;
432 
434  void subdomainsToMaterials();
435 
440  std::set<SubdomainID> coupledSubdomains() const;
441 
447  template <typename T>
448  void gatherCellSum(std::vector<T> & local, std::map<cellInfo, T> & global) const;
449 
456  template <typename T>
457  void gatherCellVector(std::vector<T> & local, std::vector<unsigned int> & n_local, std::map<cellInfo, std::vector<T>> & global);
458 
464  coupling::CouplingFields elemFeedback(const Elem * elem) const;
465 
470  void getTallyTriggerParameters(const InputParameters & params);
471 
477  void readBlockParameters(const std::string name, std::unordered_set<SubdomainID> & blocks);
478 
484  void cacheContainedCells();
485 
492  void setContainedCells(const cellInfo & cell_info,
493  const Point & hint,
494  std::map<cellInfo, containedCells> & map);
495 
504  void checkContainedCellsStructure(const cellInfo & cell_info,
505  containedCells & reference,
506  containedCells & compare) const;
507 
513  void setMinimumVolumeQRules(Order & volume_order, const std::string & type);
514 
516  std::string printNewline() const
517  {
518  if (_verbose)
519  return "\n";
520  else
521  return "";
522  }
523 
529  template <typename T>
530  void checkEmptyVector(const std::vector<T> & vector, const std::string & name) const;
531 
533  void storeElementPhase();
534 
556  void relaxAndNormalizeTally(unsigned int global_score,
557  unsigned int local_score,
558  std::shared_ptr<TallyBase> local_tally);
559 
564  void getCellMappedPhase();
565 
567  void checkCellMappedPhase();
568 
571 
577 
580 
582  void mapElemsToCells();
583 
589  void validateLocalTallies();
590 
592  void initializeTallies();
593 
598  void resetTallies();
599 
601  void getMaterialFills();
602 
608  void getPointInCell();
609 
616  std::map<cellInfo, Real> computeVolumeWeightedCellInput(
617  const std::map<SubdomainID, std::pair<unsigned int, std::string>> & var_num,
618  const std::vector<coupling::CouplingFields> * phase) const;
619 
624  void sendTemperatureToOpenMC() const;
625 
630  void sendDensityToOpenMC() const;
631 
637  Real tallyNormalization(unsigned int global_score) const;
638 
643  void checkTallySum(const unsigned int & score) const;
644 
650  void latticeOuterCheck(const Point & c, int level) const;
651 
657  void latticeOuterError(const Point & c, int level) const;
658 
664  bool findCell(const Point & point);
665 
674  void compareContainedCells(std::map<cellInfo, containedCells> & reference,
675  std::map<cellInfo, containedCells> & compare) const;
676 
677  std::unique_ptr<NumericVector<Number>> _serialized_solution;
678 
683  const bool & _output_cell_mapping;
684 
692 
695 
703 
709  unsigned int _cell_level;
710 
715  const bool & _export_properties;
716 
737 
745 
759  const bool _check_tally_sum;
760 
762  const Real & _relaxation_factor;
763 
796 
804 
813 
819 
825 
832 
834  bool _has_cell_tallies = false;
835 
838 
844 
846  std::vector<std::shared_ptr<TallyBase>> _local_tallies;
847 
849  std::vector<std::string> _all_tally_scores;
850 
856  std::vector<std::map<std::string, int>> _local_tally_score_map;
857 
859  std::vector<std::vector<unsigned int>> _tally_var_ids;
860 
862  std::vector<std::vector<std::vector<unsigned int>>> _tally_ext_var_ids;
863 
865  bool _contains_cell_tally = false;
866 
868  std::vector<SubdomainID> _density_blocks;
869 
871  std::vector<SubdomainID> _temp_blocks;
872 
874  std::unordered_set<SubdomainID> _identical_cell_fill_blocks;
875 
877  std::vector<cellInfo> _elem_to_cell{};
878 
880  std::map<cellInfo, coupling::CouplingFields> _cell_phase;
881 
884 
887 
890 
893 
899 
905 
911 
914 
917 
920 
922  std::map<cellInfo, std::vector<unsigned int>> _cell_to_elem;
923 
925  std::map<cellInfo, std::vector<unsigned int>> _local_cell_to_elem;
926 
928  std::map<cellInfo, std::unordered_set<SubdomainID>> _cell_to_elem_subdomain;
929 
931  std::map<SubdomainID, std::set<int32_t>> _subdomain_to_material;
932 
937  std::map<cellInfo, Point> _cell_to_point;
938 
943  std::map<cellInfo, Real> _cell_to_elem_volume;
944 
949  std::map<cellInfo, Real> _cell_volume;
950 
955  std::map<cellInfo, int32_t> _cell_to_material;
956 
961  std::map<cellInfo, containedCells> _cell_to_contained_material_cells;
962 
964  std::map<cellInfo, int32_t> _cell_to_n_contained;
965 
970  std::vector<openmc::Tally *> _global_tallies;
971 
973  std::vector<std::vector<std::string>> _global_tally_scores;
974 
976  std::vector<openmc::TallyEstimator> _global_tally_estimators;
977 
979  std::vector<Real> _global_sum_tally;
980 
982  std::vector<Real> _local_sum_tally;
983 
985  std::vector<Real> _local_mean_tally;
986 
988  static bool _first_transfer;
989 
991  static bool _printed_initial;
992 
995 
997  openmc::Particle _particle;
998 
1000  unsigned int _n_particles_1;
1001 
1003  std::map<std::string, std::vector<SubdomainName>> _temp_vars_to_blocks;
1004 
1006  std::map<std::string, std::vector<SubdomainName>> _density_vars_to_blocks;
1007 
1010 
1013 
1015  std::map<cellInfo, int> _n_temp;
1016 
1018  std::map<cellInfo, int> _n_rho;
1019 
1021  std::map<cellInfo, int> _n_temp_rho;
1022 
1024  std::map<cellInfo, int> _n_none;
1025 
1027  unsigned int _global_tally_index = 0;
1028 
1030  unsigned int _source_rate_index;
1031 
1032 #ifdef ENABLE_DAGMC
1033  MoabSkinner * _skinner = nullptr;
1035 
1037  std::shared_ptr<moab::DagMC> _dagmc = nullptr;
1038 #endif
1039 
1041  long unsigned int _n_openmc_cells;
1042 
1045 
1047  static constexpr Real EV_TO_JOULE = 1.6022e-19;
1048 
1050  static constexpr Real ZERO_TALLY_THRESHOLD = 1e-12;
1051 
1052 private:
1056  void dufekGudowskiParticleUpdate();
1057 
1059  std::vector<int32_t> _flattened_ids;
1060 
1062  std::vector<int32_t> _flattened_instances;
1063 
1065  containedCells _instance_offsets;
1066 
1068  std::map<cellInfo, int32_t> _n_offset;
1069 
1071  cellInfo _first_identical_cell;
1072 
1074  std::vector<int32_t> _first_identical_cell_materials;
1075 
1077  std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_temp_vars;
1078 
1080  std::map<SubdomainID, std::pair<unsigned int, std::string>> _subdomain_to_density_vars;
1081 };
OpenMCCellAverageProblem::validateLocalTallies
void validateLocalTallies()
OpenMCCellAverageProblem::_cell_to_elem
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:922
OpenMCCellAverageProblem::_n_mapped_none_elems
int _n_mapped_none_elems
Number of no-coupling elements mapped to OpenMC cells.
Definition: OpenMCCellAverageProblem.h:913
OpenMCCellAverageProblem::UNMAPPED
static constexpr int32_t UNMAPPED
Constant flag to indicate that a cell/element was unmapped.
Definition: OpenMCCellAverageProblem.h:313
OpenMCCellAverageProblem::hasDensityFeedback
bool hasDensityFeedback(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:248
relaxation::RelaxationEnum
RelaxationEnum
Type of relaxation.
Definition: CardinalEnums.h:206
OpenMCCellAverageProblem::setupProblem
void setupProblem()
Initialize the mapping of OpenMC to the MooseMesh and perform additional setup actions.
OpenMCCellAverageProblem::setContainedCells
void setContainedCells(const cellInfo &cell_info, const Point &hint, std::map< cellInfo, containedCells > &map)
OpenMCCellAverageProblem::_need_to_reinit_coupling
const bool _need_to_reinit_coupling
Definition: OpenMCCellAverageProblem.h:744
OpenMCCellAverageProblem::elemToCellID
int32_t elemToCellID(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:207
OpenMCCellAverageProblem::_printed_triso_warning
static bool _printed_triso_warning
Whether a warning has already been printed about very long setup times (for TRISOs)
Definition: OpenMCCellAverageProblem.h:994
OpenMCCellAverageProblem::gatherCellSum
void gatherCellSum(std::vector< T > &local, std::map< cellInfo, T > &global) const
OpenMCCellAverageProblem::createQRules
virtual void createQRules(QuadratureType type, Order order, Order volume_order, Order face_order, SubdomainID block, bool allow_negative_weights=true) override
OpenMCCellAverageProblem::_cell_level
unsigned int _cell_level
Definition: OpenMCCellAverageProblem.h:709
OpenMCCellAverageProblem::_n_temp
std::map< cellInfo, int > _n_temp
Number of temperature-only feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1015
OpenMCCellAverageProblem::_tally_var_ids
std::vector< std::vector< unsigned int > > _tally_var_ids
A vector of auxvariable ids added by the [Tallies] block.
Definition: OpenMCCellAverageProblem.h:859
OpenMCCellAverageProblem::_has_identical_cell_fills
const bool _has_identical_cell_fills
Definition: OpenMCCellAverageProblem.h:795
OpenMCCellAverageProblem::checkTallySum
void checkTallySum(const unsigned int &score) const
OpenMCCellAverageProblem::coupledSubdomains
std::set< SubdomainID > coupledSubdomains() const
OpenMCCellAverageProblem::DIMENSION
static constexpr int DIMENSION
Spatial dimension of the Monte Carlo problem.
Definition: OpenMCCellAverageProblem.h:316
OpenMCCellAverageProblem::_contains_cell_tally
bool _contains_cell_tally
Whether the problem contains a cell tally or not.
Definition: OpenMCCellAverageProblem.h:865
OpenMCCellAverageProblem::readBlockParameters
void readBlockParameters(const std::string name, std::unordered_set< SubdomainID > &blocks)
OpenMCCellAverageProblem::fixedPointIteration
int fixedPointIteration() const
Definition: OpenMCCellAverageProblem.h:310
OpenMCCellAverageProblem::elemToCellInfo
cellInfo elemToCellInfo(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:222
OpenMCCellAverageProblem::_cell_to_material
std::map< cellInfo, int32_t > _cell_to_material
Definition: OpenMCCellAverageProblem.h:955
OpenMCCellAverageProblem::elemToCellInstance
int32_t elemToCellInstance(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:214
OpenMCProblemBase::cellID
int32_t cellID(const int32_t index) const
trigger::TallyTriggerTypeEnum
TallyTriggerTypeEnum
Type of trigger to apply.
Definition: CardinalEnums.h:162
OpenMCCellAverageProblem::_all_tally_scores
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:849
OpenMCCellAverageProblem::_temp_blocks
std::vector< SubdomainID > _temp_blocks
Blocks in MOOSE mesh that provide temperature feedback.
Definition: OpenMCCellAverageProblem.h:871
OpenMCCellAverageProblem::_assume_separate_tallies
const bool & _assume_separate_tallies
Definition: OpenMCCellAverageProblem.h:812
OpenMCCellAverageProblem::_check_identical_cell_fills
const bool & _check_identical_cell_fills
Definition: OpenMCCellAverageProblem.h:803
OpenMCCellAverageProblem::getCellToElementSub
virtual std::unordered_set< SubdomainID > getCellToElementSub(const cellInfo &info)
Definition: OpenMCCellAverageProblem.h:133
OpenMCCellAverageProblem::_local_sum_tally
std::vector< Real > _local_sum_tally
Sum value of the local tally(s), across all bins.
Definition: OpenMCCellAverageProblem.h:982
OpenMCCellAverageProblem::_normalize_by_global
const bool _normalize_by_global
Definition: OpenMCCellAverageProblem.h:736
coupling::density
@ density
Definition: CardinalEnums.h:177
OpenMCCellAverageProblem::sendDensityToOpenMC
void sendDensityToOpenMC() const
OpenMCCellAverageProblem::checkContainedCellsStructure
void checkContainedCellsStructure(const cellInfo &cell_info, containedCells &reference, containedCells &compare) const
OpenMCCellAverageProblem::checkBlocksInMesh
void checkBlocksInMesh(const std::string name, const std::vector< SubdomainID > &ids, const std::vector< SubdomainName > &names) const
OpenMCCellAverageProblem::latticeOuterCheck
void latticeOuterCheck(const Point &c, int level) const
OpenMCCellAverageProblem::setMinimumVolumeQRules
void setMinimumVolumeQRules(Order &volume_order, const std::string &type)
OpenMCCellAverageProblem::_n_mapped_temp_density_elems
int _n_mapped_temp_density_elems
Definition: OpenMCCellAverageProblem.h:910
OpenMCCellAverageProblem::read2DBlockParameters
void read2DBlockParameters(const std::string name, std::vector< std::vector< SubdomainName >> &names, std::vector< SubdomainID > &flattened_ids)
OpenMCCellAverageProblem::_n_particles_1
unsigned int _n_particles_1
Number of particles simulated in the first iteration.
Definition: OpenMCCellAverageProblem.h:1000
coupling::density_and_temperature
@ density_and_temperature
Definition: CardinalEnums.h:178
OpenMCVolumeCalculation
Definition: OpenMCVolumeCalculation.h:28
MoabSkinner
Skins the [Mesh] according to individual bins for temperature, density, and subdomain ID.
Definition: MoabSkinner.h:19
OpenMCProblemBase::cellInfo
std::pair< int32_t, int32_t > cellInfo
Definition: OpenMCProblemBase.h:184
OpenMCCellAverageProblem::reloadDAGMC
void reloadDAGMC()
Reconstruct the DAGMC geometry after skinning.
OpenMCCellAverageProblem::cacheContainedCells
void cacheContainedCells()
OpenMCCellAverageProblem::latticeOuterError
void latticeOuterError(const Point &c, int level) const
OpenMCProblemBase::_verbose
const bool & _verbose
Whether to print diagnostic information about model setup and the transfers.
Definition: OpenMCProblemBase.h:370
OpenMCCellAverageProblem::cellHasIdenticalFill
bool cellHasIdenticalFill(const cellInfo &cell_info) const
OpenMCCellAverageProblem::_has_cell_tallies
bool _has_cell_tallies
Whether any cell tallies exist.
Definition: OpenMCCellAverageProblem.h:834
OpenMCCellAverageProblem::gatherCellVector
void gatherCellVector(std::vector< T > &local, std::vector< unsigned int > &n_local, std::map< cellInfo, std::vector< T >> &global)
OpenMCCellAverageProblem::updateMaterials
virtual void updateMaterials()
coupling::temperature
@ temperature
Definition: CardinalEnums.h:176
OpenMCCellAverageProblem
Tally includes.
Definition: OpenMCCellAverageProblem.h:65
OpenMCCellAverageProblem::_cell_to_n_contained
std::map< cellInfo, int32_t > _cell_to_n_contained
Number of material-type cells contained within a cell.
Definition: OpenMCCellAverageProblem.h:964
OpenMCCellAverageProblem::_local_tally_score_map
std::vector< std::map< std::string, int > > _local_tally_score_map
Definition: OpenMCCellAverageProblem.h:856
OpenMCCellAverageProblem::_global_tally_index
unsigned int _global_tally_index
Index in OpenMC tallies corresponding to the first global tally added by Cardinal.
Definition: OpenMCCellAverageProblem.h:1027
OpenMCCellAverageProblem::_identical_cell_fill_blocks
std::unordered_set< SubdomainID > _identical_cell_fill_blocks
Blocks for which the cell fills are identical.
Definition: OpenMCCellAverageProblem.h:874
OpenMCCellAverageProblem::_global_tally_scores
std::vector< std::vector< std::string > > _global_tally_scores
Global tally scores corresponding to '_global_tallies'.
Definition: OpenMCCellAverageProblem.h:973
TallyBase.h
OpenMCCellAverageProblem::cellToElem
virtual const std::map< cellInfo, std::vector< unsigned int > > & cellToElem() const
Definition: OpenMCCellAverageProblem.h:123
OpenMCCellAverageProblem::computeVolumeWeightedCellInput
std::map< cellInfo, Real > computeVolumeWeightedCellInput(const std::map< SubdomainID, std::pair< unsigned int, std::string >> &var_num, const std::vector< coupling::CouplingFields > *phase) const
OpenMCCellAverageProblem::cellVolume
virtual Real cellVolume(const cellInfo &cell_info) const
OpenMCCellAverageProblem::containedCells
std::unordered_map< int32_t, std::vector< int32_t > > containedCells
Definition: OpenMCCellAverageProblem.h:192
OpenMCCellAverageProblem::_n_temp_rho
std::map< cellInfo, int > _n_temp_rho
Number of temperature+density feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1021
OpenMCCellAverageProblem::_local_cell_to_elem
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:925
OpenMCCellAverageProblem::relaxAndNormalizeTally
void relaxAndNormalizeTally(unsigned int global_score, unsigned int local_score, std::shared_ptr< TallyBase > local_tally)
OpenMCCellAverageProblem::initializeTallies
void initializeTallies()
Add OpenMC tallies to facilitate the coupling.
OpenMCCellAverageProblem::ZERO_TALLY_THRESHOLD
static constexpr Real ZERO_TALLY_THRESHOLD
Tolerance for setting zero tally.
Definition: OpenMCCellAverageProblem.h:1050
OpenMCCellAverageProblem::containedMaterialCells
containedCells containedMaterialCells(const cellInfo &cell_info) const
OpenMCCellAverageProblem::cellMappedVolume
double cellMappedVolume(const cellInfo &cell_info) const
OpenMCCellAverageProblem::_n_moose_density_elems
int _n_moose_density_elems
Number of elements in the MOOSE mesh that exclusively provide density feedback.
Definition: OpenMCCellAverageProblem.h:883
OpenMCCellAverageProblem::sendTemperatureToOpenMC
void sendTemperatureToOpenMC() const
OpenMCCellAverageProblem::_export_properties
const bool & _export_properties
Definition: OpenMCCellAverageProblem.h:715
OpenMCCellAverageProblem::mapElemsToCells
void mapElemsToCells()
Populate maps of MOOSE elements to OpenMC cells.
SymmetryPointGenerator::transformPoint
Point transformPoint(const Point &p) const
Definition: SymmetryPointGenerator.C:99
OpenMCCellAverageProblem::tallyNormalization
Real tallyNormalization(unsigned int global_score) const
OpenMCCellAverageProblem::readBlockVariables
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)
OpenMCProblemBase
Definition: OpenMCProblemBase.h:49
OpenMCCellAverageProblem::_cell_to_elem_volume
std::map< cellInfo, Real > _cell_to_elem_volume
Definition: OpenMCCellAverageProblem.h:943
OpenMCCellAverageProblem::addTallyObject
void addTallyObject(const std::string &type, const std::string &name, InputParameters &moose_object_pars)
OpenMCCellAverageProblem::materialsInCells
std::vector< int32_t > materialsInCells(const containedCells &contained_cells) const
OpenMCCellAverageProblem::_n_moose_temp_density_elems
int _n_moose_temp_density_elems
Number of elements in the MOOSE mesh which provide temperature+density feedback.
Definition: OpenMCCellAverageProblem.h:889
OpenMCCellAverageProblem::_n_moose_temp_elems
int _n_moose_temp_elems
Number of elements in the MOOSE mesh that exclusively provide temperature feedback.
Definition: OpenMCCellAverageProblem.h:886
OpenMCCellAverageProblem::transformPoint
virtual Point transformPoint(const Point &pt) const
Definition: OpenMCCellAverageProblem.h:155
OpenMCCellAverageProblem::firstContainedMaterialCell
cellInfo firstContainedMaterialCell(const cellInfo &cell_info) const
OpenMCCellAverageProblem::_initial_condition
const coupling::OpenMCInitialCondition _initial_condition
Definition: OpenMCCellAverageProblem.h:691
OpenMCCellAverageProblem::elemFeedback
coupling::CouplingFields elemFeedback(const Elem *elem) const
OpenMCCellAverageProblem::tallyMultiplier
Real tallyMultiplier(unsigned int global_score) const
OpenMCCellAverageProblem::getLocalTally
const std::vector< std::shared_ptr< TallyBase > > & getLocalTally() const
Definition: OpenMCCellAverageProblem.h:271
OpenMCCellAverageProblem::_specified_temperature_feedback
const bool _specified_temperature_feedback
Definition: OpenMCCellAverageProblem.h:831
OpenMCCellAverageProblem::converged
virtual bool converged(unsigned int) override
Definition: OpenMCCellAverageProblem.h:74
OpenMCCellAverageProblem::_particle
openmc::Particle _particle
Dummy particle to reduce number of allocations of particles for cell lookup routines.
Definition: OpenMCCellAverageProblem.h:997
OpenMCCellAverageProblem::checkEmptyVector
void checkEmptyVector(const std::vector< T > &vector, const std::string &name) const
OpenMCProblemBase::_fixed_point_iteration
int _fixed_point_iteration
Definition: OpenMCProblemBase.h:420
OpenMCCellAverageProblem::checkNormalization
void checkNormalization(const Real &sum, unsigned int global_score) const
OpenMCCellAverageProblem::_n_mapped_density_elems
int _n_mapped_density_elems
Definition: OpenMCCellAverageProblem.h:904
OpenMCCellAverageProblem::_subdomain_to_material
std::map< SubdomainID, std::set< int32_t > > _subdomain_to_material
Mapping of elem subdomains to materials.
Definition: OpenMCCellAverageProblem.h:931
OpenMCCellAverageProblem::_uncoupled_volume
Real _uncoupled_volume
Total volume of uncoupled MOOSE mesh elements.
Definition: OpenMCCellAverageProblem.h:916
OpenMCCellAverageProblem::_density_blocks
std::vector< SubdomainID > _density_blocks
Blocks in MOOSE mesh that provide density feedback.
Definition: OpenMCCellAverageProblem.h:868
OpenMCCellAverageProblem::EV_TO_JOULE
static constexpr Real EV_TO_JOULE
Conversion rate from eV to Joule.
Definition: OpenMCCellAverageProblem.h:1047
OpenMCCellAverageProblem::findCell
bool findCell(const Point &point)
OpenMCCellAverageProblem::syncSolutions
virtual void syncSolutions(ExternalProblem::Direction direction) override
OpenMCCellAverageProblem::_map_density_by_cell
bool _map_density_by_cell
Definition: OpenMCCellAverageProblem.h:818
OpenMCCellAverageProblem::_elem_to_cell
std::vector< cellInfo > _elem_to_cell
Mapping of MOOSE elements to the OpenMC cell they map to (if any)
Definition: OpenMCCellAverageProblem.h:877
OpenMCCellAverageProblem::getPointInCell
void getPointInCell()
OpenMCCellAverageProblem::_n_none
std::map< cellInfo, int > _n_none
Number of none elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1024
OpenMCCellAverageProblem::_cell_to_contained_material_cells
std::map< cellInfo, containedCells > _cell_to_contained_material_cells
Definition: OpenMCCellAverageProblem.h:961
OpenMCCellAverageProblem::_n_mapped_temp_elems
int _n_mapped_temp_elems
Definition: OpenMCCellAverageProblem.h:898
coupling::CouplingFields
CouplingFields
Type of feedback in Monte Carlo simulation.
Definition: CardinalEnums.h:174
OpenMCCellAverageProblem::cellToMaterialIndex
int32_t cellToMaterialIndex(const cellInfo &cell_info)
Definition: OpenMCCellAverageProblem.h:231
OpenMCCellAverageProblem::getMaterialFills
void getMaterialFills()
Find the material filling each cell which receives density feedback.
OpenMCCellAverageProblem::_relaxation
const relaxation::RelaxationEnum _relaxation
Type of relaxation to apply to the OpenMC tallies.
Definition: OpenMCCellAverageProblem.h:694
OpenMCCellAverageProblem::_cell_phase
std::map< cellInfo, coupling::CouplingFields > _cell_phase
Phase of each cell.
Definition: OpenMCCellAverageProblem.h:880
OpenMCCellAverageProblem::_n_rho
std::map< cellInfo, int > _n_rho
Number of density-only feedback elements in each mapped OpenMC cell (global)
Definition: OpenMCCellAverageProblem.h:1018
OpenMCCellAverageProblem::_needs_to_map_cells
bool _needs_to_map_cells
Whether any spatial mapping from OpenMC's cells to the mesh is needed.
Definition: OpenMCCellAverageProblem.h:837
coupling::OpenMCInitialCondition
OpenMCInitialCondition
Where to get the initial temperature and density settings for OpenMC.
Definition: CardinalEnums.h:183
OpenMCCellAverageProblem::_material_cells_only
bool _material_cells_only
Whether non-material cells are mapped.
Definition: OpenMCCellAverageProblem.h:919
OpenMCCellAverageProblem::initialSetup
virtual void initialSetup() override
OpenMCCellAverageProblem::_serialized_solution
std::unique_ptr< NumericVector< Number > > _serialized_solution
Definition: OpenMCCellAverageProblem.h:677
OpenMCCellAverageProblem::_dagmc_universe_index
int32_t _dagmc_universe_index
Index in the OpenMC universes corresponding to the DAGMC universe.
Definition: OpenMCCellAverageProblem.h:1044
OpenMCCellAverageProblem::_tally_ext_var_ids
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:862
OpenMCCellAverageProblem::printNewline
std::string printNewline() const
For keeping the output neat when using verbose.
Definition: OpenMCCellAverageProblem.h:516
OpenMCCellAverageProblem::_printed_initial
static bool _printed_initial
Whether the diagnostic tables on initialization have already been printed.
Definition: OpenMCCellAverageProblem.h:991
OpenMCCellAverageProblem::_k_trigger
const trigger::TallyTriggerTypeEnum _k_trigger
Definition: OpenMCCellAverageProblem.h:702
OpenMCCellAverageProblem::shiftCellInstances
containedCells shiftCellInstances(const cellInfo &cell_info) const
OpenMCCellAverageProblem::_output_cell_mapping
const bool & _output_cell_mapping
Definition: OpenMCCellAverageProblem.h:683
OpenMCProblemBase.h
OpenMCCellAverageProblem::_volume_calc
OpenMCVolumeCalculation * _volume_calc
Optional volume calculation for cells which map to MOOSE.
Definition: OpenMCCellAverageProblem.h:1009
OpenMCVolumeCalculation.h
MoabSkinner.h
OpenMCCellAverageProblem::_first_transfer
static bool _first_transfer
Whether the present transfer is the first transfer.
Definition: OpenMCCellAverageProblem.h:988
OpenMCCellAverageProblem::_local_tallies
std::vector< std::shared_ptr< TallyBase > > _local_tallies
A vector of the tally objects created by the [Tallies] block.
Definition: OpenMCCellAverageProblem.h:846
OpenMCCellAverageProblem::computeCellMappedVolumes
void computeCellMappedVolumes()
order
Definition: CardinalEnums.h:65
OpenMCCellAverageProblem::getTallyScores
virtual const std::vector< std::string > & getTallyScores() const
Definition: OpenMCCellAverageProblem.h:148
OpenMCCellAverageProblem::validParams
static InputParameters validParams()
OpenMCCellAverageProblem::storeElementPhase
void storeElementPhase()
Loop over the elements in the MOOSE mesh and store the type of feedback applied by each.
OpenMCCellAverageProblem::_global_tallies
std::vector< openmc::Tally * > _global_tallies
Definition: OpenMCCellAverageProblem.h:970
OpenMCCellAverageProblem::cellFeedback
coupling::CouplingFields cellFeedback(const cellInfo &cell_info) const
OpenMCCellAverageProblem::_source_rate_index
unsigned int _source_rate_index
Index in tally_score pointing to the score used for normalizing flux tallies in eigenvalue mode.
Definition: OpenMCCellAverageProblem.h:1030
OpenMCCellAverageProblem::elemToCellIndex
int32_t elemToCellIndex(const int &elem_id) const
Definition: OpenMCCellAverageProblem.h:199
OpenMCCellAverageProblem::compareContainedCells
void compareContainedCells(std::map< cellInfo, containedCells > &reference, std::map< cellInfo, containedCells > &compare) const
OpenMCCellAverageProblem::getMaterialInEachSubdomain
std::vector< std::string > getMaterialInEachSubdomain() const
OpenMCCellAverageProblem::_global_tally_estimators
std::vector< openmc::TallyEstimator > _global_tally_estimators
Global tally estimators corresponding to '_global_tallies'.
Definition: OpenMCCellAverageProblem.h:976
OpenMCCellAverageProblem::_temp_vars_to_blocks
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:1003
OpenMCCellAverageProblem::checkCellMappedPhase
void checkCellMappedPhase()
This function is used to ensure that each OpenMC cell only maps to a single phase.
OpenMCCellAverageProblem::printTrisoHelp
void printTrisoHelp(const std::chrono::time_point< std::chrono::high_resolution_clock > &start) const
OpenMCCellAverageProblem::cellTemperature
double cellTemperature(const cellInfo &cell_info) const
OpenMCCellAverageProblem::volumeCalculation
virtual const OpenMCVolumeCalculation * volumeCalculation() const
Definition: OpenMCCellAverageProblem.h:117
OpenMCCellAverageProblem::OpenMCCellAverageProblem
OpenMCCellAverageProblem(const InputParameters &params)
OpenMCCellAverageProblem::_global_sum_tally
std::vector< Real > _global_sum_tally
Sum value of the global tally(s), across all bins.
Definition: OpenMCCellAverageProblem.h:979
OpenMCCellAverageProblem::cellMapsToSubdomain
bool cellMapsToSubdomain(const cellInfo &cell_info, const std::unordered_set< SubdomainID > &id) const
OpenMCCellAverageProblem::_check_tally_sum
const bool _check_tally_sum
Definition: OpenMCCellAverageProblem.h:759
OpenMCCellAverageProblem::resetTallies
void resetTallies()
OpenMCCellAverageProblem::_cell_volume
std::map< cellInfo, Real > _cell_volume
Definition: OpenMCCellAverageProblem.h:949
OpenMCCellAverageProblem::externalSolve
virtual void externalSolve() override
OpenMCCellAverageProblem::getCellMappedPhase
void getCellMappedPhase()
OpenMCCellAverageProblem::_local_mean_tally
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:985
OpenMCCellAverageProblem::_symmetry
const SymmetryPointGenerator * _symmetry
Userobject that maps from a partial-symmetry OpenMC model to a whole-domain [Mesh].
Definition: OpenMCCellAverageProblem.h:1012
OpenMCCellAverageProblem::printAuxVariableIO
void printAuxVariableIO()
OpenMCCellAverageProblem::_n_openmc_cells
long unsigned int _n_openmc_cells
Total number of unique OpenMC cell IDs + instances combinations.
Definition: OpenMCCellAverageProblem.h:1041
OpenMCCellAverageProblem::_needs_global_tally
const bool _needs_global_tally
Definition: OpenMCCellAverageProblem.h:843
OpenMCCellAverageProblem::hasTemperatureFeedback
bool hasTemperatureFeedback(const cellInfo &cell_info) const
Definition: OpenMCCellAverageProblem.h:260
OpenMCCellAverageProblem::initializeElementToCellMapping
void initializeElementToCellMapping()
Set up the mapping from MOOSE elements to OpenMC cells.
SymmetryPointGenerator.h
SymmetryPointGenerator
Definition: SymmetryPointGenerator.h:28
OpenMCCellAverageProblem::_cell_to_point
std::map< cellInfo, Point > _cell_to_point
Definition: OpenMCCellAverageProblem.h:937
OpenMCCellAverageProblem::getCellMappedSubdomains
void getCellMappedSubdomains()
Loop over all the OpenMC cells and get the element subdomain IDs that map to each cell.
OpenMCCellAverageProblem::_n_moose_none_elems
int _n_moose_none_elems
Number of no-coupling elements in the MOOSE mesh.
Definition: OpenMCCellAverageProblem.h:892
OpenMCCellAverageProblem::_density_vars_to_blocks
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:1006
OpenMCCellAverageProblem::transformPointToOpenMC
Point transformPointToOpenMC(const Point &pt) const
OpenMCCellAverageProblem::_relaxation_factor
const Real & _relaxation_factor
Constant relaxation factor.
Definition: OpenMCCellAverageProblem.h:762
OpenMCCellAverageProblem::hasPointTransformations
virtual bool hasPointTransformations() const
Definition: OpenMCCellAverageProblem.h:142
OpenMCCellAverageProblem::_specified_density_feedback
const bool _specified_density_feedback
Definition: OpenMCCellAverageProblem.h:824
OpenMCCellAverageProblem::_cell_to_elem_subdomain
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:928
OpenMCCellAverageProblem::addExternalVariables
virtual void addExternalVariables() override
OpenMCCellAverageProblem::subdomainsToMaterials
void subdomainsToMaterials()
Loop over the mapped cells, and build a map between subdomains to OpenMC materials.
OpenMCCellAverageProblem::getTallyTriggerParameters
void getTallyTriggerParameters(const InputParameters &params)
OpenMCCellAverageProblem::getCellLevel
unsigned int getCellLevel(const Point &c) const